min_quad_with_fixed.cpp 21 KB

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