min_quad_with_fixed.cpp 16 KB

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