min_quad_with_fixed.cpp 17 KB

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