min_quad_with_fixed.cpp 17 KB

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