min_quad_with_fixed.cpp 17 KB

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