min_quad_with_fixed.cpp 24 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Alec Jacobson <alecjacobson@gmail.com>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #include "min_quad_with_fixed.h"
  9. #include "slice.h"
  10. #include "is_symmetric.h"
  11. #include "find.h"
  12. #include "sparse.h"
  13. #include "repmat.h"
  14. #include "EPS.h"
  15. #include "cat.h"
  16. //#include <Eigen/SparseExtra>
  17. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  18. #include <iostream>
  19. #include <unsupported/Eigen/SparseExtra>
  20. #include <cassert>
  21. #include <cstdio>
  22. #include <iostream>
  23. template <typename T, typename Derivedknown>
  24. IGL_INLINE bool igl::min_quad_with_fixed_precompute(
  25. const Eigen::SparseMatrix<T>& A2,
  26. const Eigen::MatrixBase<Derivedknown> & known,
  27. const Eigen::SparseMatrix<T>& Aeq,
  28. const bool pd,
  29. min_quad_with_fixed_data<T> & data
  30. )
  31. {
  32. //#define MIN_QUAD_WITH_FIXED_CPP_DEBUG
  33. using namespace Eigen;
  34. using namespace std;
  35. const Eigen::SparseMatrix<T> A = 0.5*A2;
  36. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  37. cout<<" pre"<<endl;
  38. #endif
  39. // number of rows
  40. int n = A.rows();
  41. // cache problem size
  42. data.n = n;
  43. int neq = Aeq.rows();
  44. // default is to have 0 linear equality constraints
  45. if(Aeq.size() != 0)
  46. {
  47. assert(n == Aeq.cols() && "#Aeq.cols() should match A.rows()");
  48. }
  49. assert(A.rows() == n && "A should be square");
  50. assert(A.cols() == n && "A should be square");
  51. // number of known rows
  52. int kr = known.size();
  53. assert((kr == 0 || known.minCoeff() >= 0)&& "known indices should be in [0,n)");
  54. assert((kr == 0 || known.maxCoeff() < n) && "known indices should be in [0,n)");
  55. assert(neq <= n && "Number of equality constraints should be less than DOFs");
  56. // cache known
  57. data.known = known;
  58. // get list of unknown indices
  59. data.unknown.resize(n-kr);
  60. std::vector<bool> unknown_mask;
  61. unknown_mask.resize(n,true);
  62. for(int i = 0;i<kr;i++)
  63. {
  64. unknown_mask[known(i)] = false;
  65. }
  66. int u = 0;
  67. for(int i = 0;i<n;i++)
  68. {
  69. if(unknown_mask[i])
  70. {
  71. data.unknown(u) = i;
  72. u++;
  73. }
  74. }
  75. // get list of lagrange multiplier indices
  76. data.lagrange.resize(neq);
  77. for(int i = 0;i<neq;i++)
  78. {
  79. data.lagrange(i) = n + i;
  80. }
  81. // cache unknown followed by lagrange indices
  82. data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
  83. // Would like to do:
  84. //data.unknown_lagrange << data.unknown, data.lagrange;
  85. // but Eigen can't handle empty vectors in comma initialization
  86. // https://forum.kde.org/viewtopic.php?f=74&t=107974&p=364947#p364947
  87. if(data.unknown.size() > 0)
  88. {
  89. data.unknown_lagrange.head(data.unknown.size()) = data.unknown;
  90. }
  91. if(data.lagrange.size() > 0)
  92. {
  93. data.unknown_lagrange.tail(data.lagrange.size()) = data.lagrange;
  94. }
  95. SparseMatrix<T> Auu;
  96. slice(A,data.unknown,data.unknown,Auu);
  97. assert(Auu.size() != 0 && Auu.rows() > 0 && "There should be at least one unknown.");
  98. // Positive definiteness is *not* determined, rather it is given as a
  99. // parameter
  100. data.Auu_pd = pd;
  101. if(data.Auu_pd)
  102. {
  103. // PD implies symmetric
  104. data.Auu_sym = true;
  105. // This is an annoying assertion unless EPS can be chosen in a nicer way.
  106. //assert(is_symmetric(Auu,EPS<double>()));
  107. assert(is_symmetric(Auu,1.0) &&
  108. "Auu should be symmetric if positive definite");
  109. }else
  110. {
  111. // determine if A(unknown,unknown) is symmetric and/or positive definite
  112. VectorXi AuuI,AuuJ;
  113. MatrixXd AuuV;
  114. find(Auu,AuuI,AuuJ,AuuV);
  115. data.Auu_sym = is_symmetric(Auu,EPS<double>()*AuuV.maxCoeff());
  116. }
  117. // Determine number of linearly independent constraints
  118. int nc = 0;
  119. if(neq>0)
  120. {
  121. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  122. cout<<" qr"<<endl;
  123. #endif
  124. // QR decomposition to determine row rank in Aequ
  125. slice(Aeq,data.unknown,2,data.Aequ);
  126. assert(data.Aequ.rows() == neq &&
  127. "#Rows in Aequ should match #constraints");
  128. assert(data.Aequ.cols() == data.unknown.size() &&
  129. "#cols in Aequ should match #unknowns");
  130. data.AeqTQR.compute(data.Aequ.transpose().eval());
  131. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  132. //cout<<endl<<matlab_format(SparseMatrix<T>(data.Aequ.transpose().eval()),"AeqT")<<endl<<endl;
  133. #endif
  134. switch(data.AeqTQR.info())
  135. {
  136. case Eigen::Success:
  137. break;
  138. case Eigen::NumericalIssue:
  139. cerr<<"Error: Numerical issue."<<endl;
  140. return false;
  141. case Eigen::InvalidInput:
  142. cerr<<"Error: Invalid input."<<endl;
  143. return false;
  144. default:
  145. cerr<<"Error: Other."<<endl;
  146. return false;
  147. }
  148. nc = data.AeqTQR.rank();
  149. assert(nc<=neq &&
  150. "Rank of reduced constraints should be <= #original constraints");
  151. data.Aeq_li = nc == neq;
  152. //cout<<"data.Aeq_li: "<<data.Aeq_li<<endl;
  153. }else
  154. {
  155. data.Aeq_li = true;
  156. }
  157. if(data.Aeq_li)
  158. {
  159. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  160. cout<<" Aeq_li=true"<<endl;
  161. #endif
  162. // Append lagrange multiplier quadratic terms
  163. SparseMatrix<T> new_A;
  164. SparseMatrix<T> AeqT = Aeq.transpose();
  165. SparseMatrix<T> Z(neq,neq);
  166. // This is a bit slower. But why isn't cat fast?
  167. new_A = cat(1, cat(2, A, AeqT ),
  168. cat(2, Aeq, Z ));
  169. // precompute RHS builders
  170. if(kr > 0)
  171. {
  172. SparseMatrix<T> Aulk,Akul;
  173. // Slow
  174. slice(new_A,data.unknown_lagrange,data.known,Aulk);
  175. //// This doesn't work!!!
  176. //data.preY = Aulk + Akul.transpose();
  177. // Slow
  178. if(data.Auu_sym)
  179. {
  180. data.preY = Aulk*2;
  181. }else
  182. {
  183. slice(new_A,data.known,data.unknown_lagrange,Akul);
  184. SparseMatrix<T> AkulT = Akul.transpose();
  185. data.preY = Aulk + AkulT;
  186. }
  187. }else
  188. {
  189. data.preY.resize(data.unknown_lagrange.size(),0);
  190. }
  191. // Positive definite and no equality constraints (Positive definiteness
  192. // implies symmetric)
  193. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  194. cout<<" factorize"<<endl;
  195. #endif
  196. if(data.Auu_pd && neq == 0)
  197. {
  198. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  199. cout<<" llt"<<endl;
  200. #endif
  201. data.llt.compute(Auu);
  202. switch(data.llt.info())
  203. {
  204. case Eigen::Success:
  205. break;
  206. case Eigen::NumericalIssue:
  207. cerr<<"Error: Numerical issue."<<endl;
  208. return false;
  209. default:
  210. cerr<<"Error: Other."<<endl;
  211. return false;
  212. }
  213. data.solver_type = min_quad_with_fixed_data<T>::LLT;
  214. }else
  215. {
  216. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  217. cout<<" ldlt"<<endl;
  218. #endif
  219. // Either not PD or there are equality constraints
  220. SparseMatrix<T> NA;
  221. slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  222. data.NA = NA;
  223. // Ideally we'd use LDLT but Eigen doesn't support positive semi-definite
  224. // matrices:
  225. // http://forum.kde.org/viewtopic.php?f=74&t=106962&p=291990#p291990
  226. if(data.Auu_sym && false)
  227. {
  228. data.ldlt.compute(NA);
  229. switch(data.ldlt.info())
  230. {
  231. case Eigen::Success:
  232. break;
  233. case Eigen::NumericalIssue:
  234. cerr<<"Error: Numerical issue."<<endl;
  235. return false;
  236. default:
  237. cerr<<"Error: Other."<<endl;
  238. return false;
  239. }
  240. data.solver_type = min_quad_with_fixed_data<T>::LDLT;
  241. }else
  242. {
  243. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  244. cout<<" lu"<<endl;
  245. #endif
  246. // Resort to LU
  247. // Bottleneck >1/2
  248. data.lu.compute(NA);
  249. //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
  250. switch(data.lu.info())
  251. {
  252. case Eigen::Success:
  253. break;
  254. case Eigen::NumericalIssue:
  255. cerr<<"Error: Numerical issue."<<endl;
  256. return false;
  257. case Eigen::InvalidInput:
  258. cerr<<"Error: Invalid Input."<<endl;
  259. return false;
  260. default:
  261. cerr<<"Error: Other."<<endl;
  262. return false;
  263. }
  264. data.solver_type = min_quad_with_fixed_data<T>::LU;
  265. }
  266. }
  267. }else
  268. {
  269. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  270. cout<<" Aeq_li=false"<<endl;
  271. #endif
  272. data.neq = neq;
  273. const int nu = data.unknown.size();
  274. //cout<<"nu: "<<nu<<endl;
  275. //cout<<"neq: "<<neq<<endl;
  276. //cout<<"nc: "<<nc<<endl;
  277. //cout<<" matrixR"<<endl;
  278. SparseMatrix<T> AeqTR,AeqTQ;
  279. AeqTR = data.AeqTQR.matrixR();
  280. // This shouldn't be necessary
  281. AeqTR.prune(static_cast<T>(0.0));
  282. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  283. cout<<" matrixQ"<<endl;
  284. #endif
  285. // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
  286. // http://forum.kde.org/viewtopic.php?f=74&t=117500
  287. AeqTQ = data.AeqTQR.matrixQ();
  288. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  289. cout<<" prune"<<endl;
  290. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  291. #endif
  292. // This shouldn't be necessary
  293. AeqTQ.prune(static_cast<T>(0.0));
  294. //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
  295. //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
  296. //cout<<" perms"<<endl;
  297. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  298. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  299. cout<<" perm"<<endl;
  300. #endif
  301. SparseMatrix<T> I(neq,neq);
  302. I.setIdentity();
  303. data.AeqTE = data.AeqTQR.colsPermutation() * I;
  304. data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
  305. assert(AeqTR.rows() == nu && "#rows in AeqTR should match #unknowns");
  306. assert(AeqTR.cols() == neq && "#cols in AeqTR should match #constraints");
  307. assert(AeqTQ.rows() == nu && "#rows in AeqTQ should match #unknowns");
  308. assert(AeqTQ.cols() == nu && "#cols in AeqTQ should match #unknowns");
  309. //cout<<" slice"<<endl;
  310. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  311. cout<<" slice"<<endl;
  312. #endif
  313. data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
  314. data.AeqTQ1T = data.AeqTQ1.transpose().eval();
  315. // ALREADY TRIM (Not 100% sure about this)
  316. data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
  317. data.AeqTR1T = data.AeqTR1.transpose().eval();
  318. //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
  319. // Null space
  320. data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
  321. data.AeqTQ2T = data.AeqTQ2.transpose().eval();
  322. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  323. cout<<" proj"<<endl;
  324. #endif
  325. // Projected hessian
  326. SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
  327. {
  328. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  329. cout<<" factorize"<<endl;
  330. #endif
  331. // QRAuu should always be PD
  332. data.llt.compute(QRAuu);
  333. switch(data.llt.info())
  334. {
  335. case Eigen::Success:
  336. break;
  337. case Eigen::NumericalIssue:
  338. cerr<<"Error: Numerical issue."<<endl;
  339. return false;
  340. default:
  341. cerr<<"Error: Other."<<endl;
  342. return false;
  343. }
  344. data.solver_type = min_quad_with_fixed_data<T>::QR_LLT;
  345. }
  346. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  347. cout<<" smash"<<endl;
  348. #endif
  349. // Known value multiplier
  350. SparseMatrix<T> Auk;
  351. slice(A,data.unknown,data.known,Auk);
  352. SparseMatrix<T> Aku;
  353. slice(A,data.known,data.unknown,Aku);
  354. SparseMatrix<T> AkuT = Aku.transpose();
  355. data.preY = Auk + AkuT;
  356. // Needed during solve
  357. data.Auu = Auu;
  358. slice(Aeq,data.known,2,data.Aeqk);
  359. assert(data.Aeqk.rows() == neq);
  360. assert(data.Aeqk.cols() == data.known.size());
  361. }
  362. return true;
  363. }
  364. template <
  365. typename T,
  366. typename DerivedB,
  367. typename DerivedY,
  368. typename DerivedBeq,
  369. typename DerivedZ,
  370. typename Derivedsol>
  371. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  372. const min_quad_with_fixed_data<T> & data,
  373. const Eigen::MatrixBase<DerivedB> & B,
  374. const Eigen::MatrixBase<DerivedY> & Y,
  375. const Eigen::MatrixBase<DerivedBeq> & Beq,
  376. Eigen::PlainObjectBase<DerivedZ> & Z,
  377. Eigen::PlainObjectBase<Derivedsol> & sol)
  378. {
  379. using namespace std;
  380. using namespace Eigen;
  381. typedef Matrix<T,Dynamic,1> VectorXT;
  382. typedef Matrix<T,Dynamic,Dynamic> MatrixXT;
  383. // number of known rows
  384. int kr = data.known.size();
  385. if(kr!=0)
  386. {
  387. assert(kr == Y.rows());
  388. }
  389. // number of columns to solve
  390. int cols = Y.cols();
  391. assert(B.cols() == 1 || B.cols() == cols);
  392. assert(Beq.size() == 0 || Beq.cols() == 1 || Beq.cols() == cols);
  393. // resize output
  394. Z.resize(data.n,cols);
  395. // Set known values
  396. for(int i = 0;i < kr;i++)
  397. {
  398. for(int j = 0;j < cols;j++)
  399. {
  400. Z(data.known(i),j) = Y(i,j);
  401. }
  402. }
  403. if(data.Aeq_li)
  404. {
  405. // number of lagrange multipliers aka linear equality constraints
  406. int neq = data.lagrange.size();
  407. // append lagrange multiplier rhs's
  408. MatrixXT BBeq(B.rows() + Beq.rows(),cols);
  409. if(B.size() > 0)
  410. {
  411. BBeq.topLeftCorner(B.rows(),cols) = B.replicate(1,B.cols()==cols?1:cols);
  412. }
  413. if(Beq.size() > 0)
  414. {
  415. BBeq.bottomLeftCorner(Beq.rows(),cols) = -2.0*Beq.replicate(1,Beq.cols()==cols?1:cols);
  416. }
  417. // Build right hand side
  418. MatrixXT BBequlcols;
  419. igl::slice(BBeq,data.unknown_lagrange,1,BBequlcols);
  420. MatrixXT NB;
  421. if(kr == 0)
  422. {
  423. NB = BBequlcols;
  424. }else
  425. {
  426. NB = data.preY * Y + BBequlcols;
  427. }
  428. //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
  429. //cout<<matlab_format(NB,"NB")<<endl;
  430. switch(data.solver_type)
  431. {
  432. case igl::min_quad_with_fixed_data<T>::LLT:
  433. sol = data.llt.solve(NB);
  434. break;
  435. case igl::min_quad_with_fixed_data<T>::LDLT:
  436. sol = data.ldlt.solve(NB);
  437. break;
  438. case igl::min_quad_with_fixed_data<T>::LU:
  439. // Not a bottleneck
  440. sol = data.lu.solve(NB);
  441. break;
  442. default:
  443. cerr<<"Error: invalid solver type"<<endl;
  444. return false;
  445. }
  446. //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
  447. // Now sol contains sol/-0.5
  448. sol *= -0.5;
  449. // Now sol contains solution
  450. // Place solution in Z
  451. for(int i = 0;i<(sol.rows()-neq);i++)
  452. {
  453. for(int j = 0;j<sol.cols();j++)
  454. {
  455. Z(data.unknown_lagrange(i),j) = sol(i,j);
  456. }
  457. }
  458. }else
  459. {
  460. assert(data.solver_type == min_quad_with_fixed_data<T>::QR_LLT);
  461. MatrixXT eff_Beq;
  462. // Adjust Aeq rhs to include known parts
  463. eff_Beq =
  464. //data.AeqTQR.colsPermutation().transpose() * (-data.Aeqk * Y + Beq);
  465. data.AeqTET * (-data.Aeqk * Y + Beq.replicate(1,Beq.cols()==cols?1:cols));
  466. // Where did this -0.5 come from? Probably the same place as above.
  467. MatrixXT Bu;
  468. slice(B,data.unknown,1,Bu);
  469. MatrixXT NB;
  470. NB = -0.5*(Bu.replicate(1,B.cols()==cols?1:cols) + data.preY * Y);
  471. // Trim eff_Beq
  472. const int nc = data.AeqTQR.rank();
  473. const int neq = Beq.rows();
  474. eff_Beq = eff_Beq.topLeftCorner(nc,cols).eval();
  475. data.AeqTR1T.template triangularView<Lower>().solveInPlace(eff_Beq);
  476. // Now eff_Beq = (data.AeqTR1T \ (data.AeqTET * (-data.Aeqk * Y + Beq)))
  477. MatrixXT lambda_0;
  478. lambda_0 = data.AeqTQ1 * eff_Beq;
  479. //cout<<matlab_format(lambda_0,"lambda_0")<<endl;
  480. MatrixXT QRB;
  481. QRB = -data.AeqTQ2T * (data.Auu * lambda_0) + data.AeqTQ2T * NB;
  482. Derivedsol lambda;
  483. lambda = data.llt.solve(QRB);
  484. // prepare output
  485. Derivedsol solu;
  486. solu = data.AeqTQ2 * lambda + lambda_0;
  487. // http://www.math.uh.edu/~rohop/fall_06/Chapter3.pdf
  488. Derivedsol solLambda;
  489. {
  490. Derivedsol temp1,temp2;
  491. temp1 = (data.AeqTQ1T * NB - data.AeqTQ1T * data.Auu * solu);
  492. data.AeqTR1.template triangularView<Upper>().solveInPlace(temp1);
  493. //cout<<matlab_format(temp1,"temp1")<<endl;
  494. temp2 = Derivedsol::Zero(neq,cols);
  495. temp2.topLeftCorner(nc,cols) = temp1;
  496. //solLambda = data.AeqTQR.colsPermutation() * temp2;
  497. solLambda = data.AeqTE * temp2;
  498. }
  499. // sol is [Z(unknown);Lambda]
  500. assert(data.unknown.size() == solu.rows());
  501. assert(cols == solu.cols());
  502. assert(data.neq == neq);
  503. assert(data.neq == solLambda.rows());
  504. assert(cols == solLambda.cols());
  505. sol.resize(data.unknown.size()+data.neq,cols);
  506. sol.block(0,0,solu.rows(),solu.cols()) = solu;
  507. sol.block(solu.rows(),0,solLambda.rows(),solLambda.cols()) = solLambda;
  508. for(int u = 0;u<data.unknown.size();u++)
  509. {
  510. for(int j = 0;j<Z.cols();j++)
  511. {
  512. Z(data.unknown(u),j) = solu(u,j);
  513. }
  514. }
  515. }
  516. return true;
  517. }
  518. template <
  519. typename T,
  520. typename DerivedB,
  521. typename DerivedY,
  522. typename DerivedBeq,
  523. typename DerivedZ>
  524. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  525. const min_quad_with_fixed_data<T> & data,
  526. const Eigen::MatrixBase<DerivedB> & B,
  527. const Eigen::MatrixBase<DerivedY> & Y,
  528. const Eigen::MatrixBase<DerivedBeq> & Beq,
  529. Eigen::PlainObjectBase<DerivedZ> & Z)
  530. {
  531. Eigen::Matrix<typename DerivedZ::Scalar, Eigen::Dynamic, Eigen::Dynamic> sol;
  532. return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
  533. }
  534. template <
  535. typename T,
  536. typename Derivedknown,
  537. typename DerivedB,
  538. typename DerivedY,
  539. typename DerivedBeq,
  540. typename DerivedZ>
  541. IGL_INLINE bool igl::min_quad_with_fixed(
  542. const Eigen::SparseMatrix<T>& A,
  543. const Eigen::MatrixBase<DerivedB> & B,
  544. const Eigen::MatrixBase<Derivedknown> & known,
  545. const Eigen::MatrixBase<DerivedY> & Y,
  546. const Eigen::SparseMatrix<T>& Aeq,
  547. const Eigen::MatrixBase<DerivedBeq> & Beq,
  548. const bool pd,
  549. Eigen::PlainObjectBase<DerivedZ> & Z)
  550. {
  551. min_quad_with_fixed_data<T> data;
  552. if(!min_quad_with_fixed_precompute(A,known,Aeq,pd,data))
  553. {
  554. return false;
  555. }
  556. return min_quad_with_fixed_solve(data,B,Y,Beq,Z);
  557. }
  558. #ifdef IGL_STATIC_LIBRARY
  559. // Explicit template instantiation
  560. #if EIGEN_VERSION_AT_LEAST(3,3,0)
  561. #else
  562. template bool igl::min_quad_with_fixed_solve<double, Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<double>, Eigen::Matrix<double, -1, 1, 0, -1, 1> const>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::CwiseUnaryOp<Eigen::internal::scalar_multiple_op<double>, Eigen::Matrix<double, -1, 1, 0, -1, 1> const> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  563. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<double>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::CwiseNullaryOp<Eigen::internal::scalar_constant_op<double>, Eigen::Matrix<double, -1, 1, 0, -1, 1> > > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  564. #endif
  565. template bool igl::min_quad_with_fixed<double, Eigen::Matrix<int, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  566. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  567. template bool igl::min_quad_with_fixed_precompute<double, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int> const&, bool, igl::min_quad_with_fixed_data<double>&);
  568. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  569. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  570. template bool igl::min_quad_with_fixed_precompute<double, Eigen::Matrix<int, -1, -1, 0, -1, -1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int> const&, bool, igl::min_quad_with_fixed_data<double>&);
  571. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  572. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  573. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  574. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  575. template bool igl::min_quad_with_fixed<double, Eigen::Matrix<int, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, -1, 0, -1, -1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<int, -1, -1, 0, -1, -1> > const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, Eigen::SparseMatrix<double, 0, int> const&, Eigen::MatrixBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> > const&, bool, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&);
  576. #endif