min_quad_with_fixed.cpp 23 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 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 "matlab_format.h"
  15. #include "EPS.h"
  16. #include "cat.h"
  17. //#include <Eigen/SparseExtra>
  18. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  19. #include <iostream>
  20. #include <unsupported/Eigen/SparseExtra>
  21. #include <cassert>
  22. #include <cstdio>
  23. #include <iostream>
  24. template <typename T, typename Derivedknown>
  25. IGL_INLINE bool igl::min_quad_with_fixed_precompute(
  26. const Eigen::SparseMatrix<T>& A2,
  27. const Eigen::PlainObjectBase<Derivedknown> & known,
  28. const Eigen::SparseMatrix<T>& Aeq,
  29. const bool pd,
  30. min_quad_with_fixed_data<T> & data
  31. )
  32. {
  33. //#define MIN_QUAD_WITH_FIXED_CPP_DEBUG
  34. using namespace Eigen;
  35. using namespace std;
  36. const Eigen::SparseMatrix<T> A = 0.5*A2;
  37. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  38. cout<<" pre"<<endl;
  39. #endif
  40. // number of rows
  41. int n = A.rows();
  42. // cache problem size
  43. data.n = n;
  44. int neq = Aeq.rows();
  45. // default is to have 0 linear equality constraints
  46. if(Aeq.size() != 0)
  47. {
  48. assert(n == Aeq.cols() && "#Aeq.cols() should match A.rows()");
  49. }
  50. assert(A.rows() == n && "A should be square");
  51. assert(A.cols() == n && "A should be square");
  52. // number of known rows
  53. int kr = known.size();
  54. assert((kr == 0 || known.minCoeff() >= 0)&& "known indices should be in [0,n)");
  55. assert((kr == 0 || known.maxCoeff() < n) && "known indices should be in [0,n)");
  56. assert(neq <= n && "Number of equality constraints should be less than DOFs");
  57. // cache known
  58. data.known = known;
  59. // get list of unknown indices
  60. data.unknown.resize(n-kr);
  61. std::vector<bool> unknown_mask;
  62. unknown_mask.resize(n,true);
  63. for(int i = 0;i<kr;i++)
  64. {
  65. unknown_mask[known(i)] = false;
  66. }
  67. int u = 0;
  68. for(int i = 0;i<n;i++)
  69. {
  70. if(unknown_mask[i])
  71. {
  72. data.unknown(u) = i;
  73. u++;
  74. }
  75. }
  76. // get list of lagrange multiplier indices
  77. data.lagrange.resize(neq);
  78. for(int i = 0;i<neq;i++)
  79. {
  80. data.lagrange(i) = n + i;
  81. }
  82. // cache unknown followed by lagrange indices
  83. data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
  84. data.unknown_lagrange << data.unknown, data.lagrange;
  85. SparseMatrix<T> Auu;
  86. slice(A,data.unknown,data.unknown,Auu);
  87. assert(Auu.size() > 0 && "There should be at least one unknown.");
  88. // Positive definiteness is *not* determined, rather it is given as a
  89. // parameter
  90. data.Auu_pd = pd;
  91. if(data.Auu_pd)
  92. {
  93. // PD implies symmetric
  94. data.Auu_sym = true;
  95. // This is an annoying assertion unless EPS can be chosen in a nicer way.
  96. //assert(is_symmetric(Auu,EPS<double>()));
  97. assert(is_symmetric(Auu,1.0) &&
  98. "Auu should be symmetric if positive definite");
  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. "#Rows in Aequ should match #constraints");
  118. assert(data.Aequ.cols() == data.unknown.size() &&
  119. "#cols in Aequ should match #unknowns");
  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. "Rank of reduced constraints should be <= #original constraints");
  141. data.Aeq_li = nc == neq;
  142. //cout<<"data.Aeq_li: "<<data.Aeq_li<<endl;
  143. }else
  144. {
  145. data.Aeq_li = true;
  146. }
  147. if(data.Aeq_li)
  148. {
  149. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  150. cout<<" Aeq_li=true"<<endl;
  151. #endif
  152. // Append lagrange multiplier quadratic terms
  153. SparseMatrix<T> new_A;
  154. SparseMatrix<T> AeqT = Aeq.transpose();
  155. SparseMatrix<T> Z(neq,neq);
  156. // This is a bit slower. But why isn't cat fast?
  157. new_A = cat(1, cat(2, A, AeqT ),
  158. cat(2, Aeq, Z ));
  159. // precompute RHS builders
  160. if(kr > 0)
  161. {
  162. SparseMatrix<T> Aulk,Akul;
  163. // Slow
  164. slice(new_A,data.unknown_lagrange,data.known,Aulk);
  165. //// This doesn't work!!!
  166. //data.preY = Aulk + Akul.transpose();
  167. // Slow
  168. if(data.Auu_sym)
  169. {
  170. data.preY = Aulk*2;
  171. }else
  172. {
  173. slice(new_A,data.known,data.unknown_lagrange,Akul);
  174. SparseMatrix<T> AkulT = Akul.transpose();
  175. data.preY = Aulk + AkulT;
  176. }
  177. }else
  178. {
  179. data.preY.resize(data.unknown_lagrange.size(),0);
  180. }
  181. // Positive definite and no equality constraints (Postive definiteness
  182. // implies symmetric)
  183. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  184. cout<<" factorize"<<endl;
  185. #endif
  186. if(data.Auu_pd && neq == 0)
  187. {
  188. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  189. cout<<" llt"<<endl;
  190. #endif
  191. data.llt.compute(Auu);
  192. switch(data.llt.info())
  193. {
  194. case Eigen::Success:
  195. break;
  196. case Eigen::NumericalIssue:
  197. cerr<<"Error: Numerical issue."<<endl;
  198. return false;
  199. default:
  200. cerr<<"Error: Other."<<endl;
  201. return false;
  202. }
  203. data.solver_type = min_quad_with_fixed_data<T>::LLT;
  204. }else
  205. {
  206. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  207. cout<<" ldlt"<<endl;
  208. #endif
  209. // Either not PD or there are equality constraints
  210. SparseMatrix<T> NA;
  211. slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  212. data.NA = NA;
  213. // Ideally we'd use LDLT but Eigen doesn't support positive semi-definite
  214. // matrices:
  215. // http://forum.kde.org/viewtopic.php?f=74&t=106962&p=291990#p291990
  216. if(data.Auu_sym && false)
  217. {
  218. data.ldlt.compute(NA);
  219. switch(data.ldlt.info())
  220. {
  221. case Eigen::Success:
  222. break;
  223. case Eigen::NumericalIssue:
  224. cerr<<"Error: Numerical issue."<<endl;
  225. return false;
  226. default:
  227. cerr<<"Error: Other."<<endl;
  228. return false;
  229. }
  230. data.solver_type = min_quad_with_fixed_data<T>::LDLT;
  231. }else
  232. {
  233. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  234. cout<<" lu"<<endl;
  235. #endif
  236. // Resort to LU
  237. // Bottleneck >1/2
  238. data.lu.compute(NA);
  239. //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
  240. switch(data.lu.info())
  241. {
  242. case Eigen::Success:
  243. break;
  244. case Eigen::NumericalIssue:
  245. cerr<<"Error: Numerical issue."<<endl;
  246. return false;
  247. case Eigen::InvalidInput:
  248. cerr<<"Error: Invalid Input."<<endl;
  249. return false;
  250. default:
  251. cerr<<"Error: Other."<<endl;
  252. return false;
  253. }
  254. data.solver_type = min_quad_with_fixed_data<T>::LU;
  255. }
  256. }
  257. }else
  258. {
  259. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  260. cout<<" Aeq_li=false"<<endl;
  261. #endif
  262. data.neq = neq;
  263. const int nu = data.unknown.size();
  264. //cout<<"nu: "<<nu<<endl;
  265. //cout<<"neq: "<<neq<<endl;
  266. //cout<<"nc: "<<nc<<endl;
  267. //cout<<" matrixR"<<endl;
  268. SparseMatrix<T> AeqTR,AeqTQ;
  269. AeqTR = data.AeqTQR.matrixR();
  270. // This shouldn't be necessary
  271. AeqTR.prune(0.0);
  272. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  273. cout<<" matrixQ"<<endl;
  274. #endif
  275. // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
  276. // http://forum.kde.org/viewtopic.php?f=74&t=117500
  277. AeqTQ = data.AeqTQR.matrixQ();
  278. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  279. cout<<" prune"<<endl;
  280. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  281. #endif
  282. // This shouldn't be necessary
  283. AeqTQ.prune(0.0);
  284. //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
  285. //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
  286. //cout<<" perms"<<endl;
  287. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  288. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  289. cout<<" perm"<<endl;
  290. #endif
  291. SparseMatrix<double> I(neq,neq);
  292. I.setIdentity();
  293. data.AeqTE = data.AeqTQR.colsPermutation() * I;
  294. data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
  295. assert(AeqTR.rows() == nu && "#rows in AeqTR should match #unknowns");
  296. assert(AeqTR.cols() == neq && "#cols in AeqTR should match #constraints");
  297. assert(AeqTQ.rows() == nu && "#rows in AeqTQ should match #unknowns");
  298. assert(AeqTQ.cols() == nu && "#cols in AeqTQ should match #unknowns");
  299. //cout<<" slice"<<endl;
  300. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  301. cout<<" slice"<<endl;
  302. #endif
  303. data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
  304. data.AeqTQ1T = data.AeqTQ1.transpose().eval();
  305. // ALREADY TRIM (Not 100% sure about this)
  306. data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
  307. data.AeqTR1T = data.AeqTR1.transpose().eval();
  308. //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
  309. // Null space
  310. data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
  311. data.AeqTQ2T = data.AeqTQ2.transpose().eval();
  312. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  313. cout<<" proj"<<endl;
  314. #endif
  315. // Projected hessian
  316. SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
  317. {
  318. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  319. cout<<" factorize"<<endl;
  320. #endif
  321. // QRAuu should always be PD
  322. data.llt.compute(QRAuu);
  323. switch(data.llt.info())
  324. {
  325. case Eigen::Success:
  326. break;
  327. case Eigen::NumericalIssue:
  328. cerr<<"Error: Numerical issue."<<endl;
  329. return false;
  330. default:
  331. cerr<<"Error: Other."<<endl;
  332. return false;
  333. }
  334. data.solver_type = min_quad_with_fixed_data<T>::QR_LLT;
  335. }
  336. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  337. cout<<" smash"<<endl;
  338. #endif
  339. // Known value multiplier
  340. SparseMatrix<T> Auk;
  341. slice(A,data.unknown,data.known,Auk);
  342. SparseMatrix<T> Aku;
  343. slice(A,data.known,data.unknown,Aku);
  344. SparseMatrix<T> AkuT = Aku.transpose();
  345. data.preY = Auk + AkuT;
  346. // Needed during solve
  347. data.Auu = Auu;
  348. slice(Aeq,data.known,2,data.Aeqk);
  349. assert(data.Aeqk.rows() == neq);
  350. assert(data.Aeqk.cols() == data.known.size());
  351. }
  352. return true;
  353. }
  354. template <
  355. typename T,
  356. typename DerivedB,
  357. typename DerivedY,
  358. typename DerivedBeq,
  359. typename DerivedZ,
  360. typename Derivedsol>
  361. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  362. const min_quad_with_fixed_data<T> & data,
  363. const Eigen::PlainObjectBase<DerivedB> & B,
  364. const Eigen::PlainObjectBase<DerivedY> & Y,
  365. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  366. Eigen::PlainObjectBase<DerivedZ> & Z,
  367. Eigen::PlainObjectBase<Derivedsol> & sol)
  368. {
  369. using namespace std;
  370. using namespace Eigen;
  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::Matrix<typename DerivedZ::Scalar, Eigen::Dynamic, Eigen::Dynamic> sol;
  517. return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
  518. }
  519. template <
  520. typename T,
  521. typename Derivedknown,
  522. typename DerivedB,
  523. typename DerivedY,
  524. typename DerivedBeq,
  525. typename DerivedZ>
  526. IGL_INLINE bool igl::min_quad_with_fixed(
  527. const Eigen::SparseMatrix<T>& A,
  528. const Eigen::PlainObjectBase<DerivedB> & B,
  529. const Eigen::PlainObjectBase<Derivedknown> & known,
  530. const Eigen::PlainObjectBase<DerivedY> & Y,
  531. const Eigen::SparseMatrix<T>& Aeq,
  532. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  533. const bool pd,
  534. Eigen::PlainObjectBase<DerivedZ> & Z)
  535. {
  536. min_quad_with_fixed_data<T> data;
  537. if(!min_quad_with_fixed_precompute(A,known,Aeq,pd,data))
  538. {
  539. return false;
  540. }
  541. return min_quad_with_fixed_solve(data,B,Y,Beq,Z);
  542. }
  543. #ifdef IGL_STATIC_LIBRARY
  544. // Explicit template specialization
  545. 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> >&);
  546. 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>&);
  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_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> >&);
  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>, 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> >&);
  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>, 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> >&);
  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<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> >&);
  553. 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> >&);
  554. 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> >&);
  555. 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> >&);
  556. #endif