min_quad_with_fixed.cpp 21 KB

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