min_quad_with_fixed.cpp 24 KB

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