2
0

min_quad_with_fixed.cpp 44 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2016 Alec Jacobson <[email protected]>
  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 "EPS.h"
  15. #include "cat.h"
  16. //#include <Eigen/SparseExtra>
  17. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  18. #include <iostream>
  19. #include <unsupported/Eigen/SparseExtra>
  20. #include <cassert>
  21. #include <cstdio>
  22. #include <igl/matlab_format.h>
  23. #include <type_traits>
  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(known.cols() == 1 && "known should be a vector");
  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. // FIXME: This is *NOT* generic and introduces a copy.
  60. data.known = known.template cast<int>();
  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, 0)] = 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. // Would like to do:
  87. //data.unknown_lagrange << data.unknown, data.lagrange;
  88. // but Eigen can't handle empty vectors in comma initialization
  89. // https://forum.kde.org/viewtopic.php?f=74&t=107974&p=364947#p364947
  90. if(data.unknown.size() > 0)
  91. {
  92. data.unknown_lagrange.head(data.unknown.size()) = data.unknown;
  93. }
  94. if(data.lagrange.size() > 0)
  95. {
  96. data.unknown_lagrange.tail(data.lagrange.size()) = data.lagrange;
  97. }
  98. SparseMatrix<T> Auu;
  99. slice(A,data.unknown,data.unknown,Auu);
  100. assert(Auu.size() != 0 && Auu.rows() > 0 && "There should be at least one unknown.");
  101. // Positive definiteness is *not* determined, rather it is given as a
  102. // parameter
  103. data.Auu_pd = pd;
  104. if(data.Auu_pd)
  105. {
  106. // PD implies symmetric
  107. data.Auu_sym = true;
  108. // This is an annoying assertion unless EPS can be chosen in a nicer way.
  109. //assert(is_symmetric(Auu,EPS<double>()));
  110. assert(is_symmetric(Auu,1.0) &&
  111. "Auu should be symmetric if positive definite");
  112. }else
  113. {
  114. // determine if A(unknown,unknown) is symmetric and/or positive definite
  115. VectorXi AuuI,AuuJ;
  116. MatrixXd AuuV;
  117. find(Auu,AuuI,AuuJ,AuuV);
  118. data.Auu_sym = is_symmetric(Auu,EPS<double>()*AuuV.maxCoeff());
  119. }
  120. // Determine number of linearly independent constraints
  121. int nc = 0;
  122. if(neq>0)
  123. {
  124. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  125. cout<<" qr"<<endl;
  126. #endif
  127. // QR decomposition to determine row rank in Aequ
  128. slice(Aeq,data.unknown,2,data.Aequ);
  129. assert(data.Aequ.rows() == neq &&
  130. "#Rows in Aequ should match #constraints");
  131. assert(data.Aequ.cols() == data.unknown.size() &&
  132. "#cols in Aequ should match #unknowns");
  133. data.AeqTQR.compute(data.Aequ.transpose().eval());
  134. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  135. //cout<<endl<<matlab_format(SparseMatrix<T>(data.Aequ.transpose().eval()),"AeqT")<<endl<<endl;
  136. #endif
  137. switch(data.AeqTQR.info())
  138. {
  139. case Eigen::Success:
  140. break;
  141. case Eigen::NumericalIssue:
  142. cerr<<"Error: Numerical issue."<<endl;
  143. return false;
  144. case Eigen::InvalidInput:
  145. cerr<<"Error: Invalid input."<<endl;
  146. return false;
  147. default:
  148. cerr<<"Error: Other."<<endl;
  149. return false;
  150. }
  151. nc = data.AeqTQR.rank();
  152. assert(nc<=neq &&
  153. "Rank of reduced constraints should be <= #original constraints");
  154. data.Aeq_li = nc == neq;
  155. //cout<<"data.Aeq_li: "<<data.Aeq_li<<endl;
  156. }else
  157. {
  158. data.Aeq_li = true;
  159. }
  160. if(data.Aeq_li)
  161. {
  162. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  163. cout<<" Aeq_li=true"<<endl;
  164. #endif
  165. // Append lagrange multiplier quadratic terms
  166. SparseMatrix<T> new_A;
  167. SparseMatrix<T> AeqT = Aeq.transpose();
  168. SparseMatrix<T> Z(neq,neq);
  169. // This is a bit slower. But why isn't cat fast?
  170. new_A = cat(1, cat(2, A, AeqT ),
  171. cat(2, Aeq, Z ));
  172. // precompute RHS builders
  173. if(kr > 0)
  174. {
  175. SparseMatrix<T> Aulk,Akul;
  176. // Slow
  177. slice(new_A,data.unknown_lagrange,data.known,Aulk);
  178. //// This doesn't work!!!
  179. //data.preY = Aulk + Akul.transpose();
  180. // Slow
  181. if(data.Auu_sym)
  182. {
  183. data.preY = Aulk*2;
  184. }else
  185. {
  186. slice(new_A,data.known,data.unknown_lagrange,Akul);
  187. SparseMatrix<T> AkulT = Akul.transpose();
  188. data.preY = Aulk + AkulT;
  189. }
  190. }else
  191. {
  192. data.preY.resize(data.unknown_lagrange.size(),0);
  193. }
  194. // Positive definite and no equality constraints (Positive definiteness
  195. // implies symmetric)
  196. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  197. cout<<" factorize"<<endl;
  198. #endif
  199. if(data.Auu_pd && neq == 0)
  200. {
  201. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  202. cout<<" llt"<<endl;
  203. #endif
  204. data.llt.compute(Auu);
  205. switch(data.llt.info())
  206. {
  207. case Eigen::Success:
  208. break;
  209. case Eigen::NumericalIssue:
  210. cerr<<"Error: Numerical issue."<<endl;
  211. return false;
  212. default:
  213. cerr<<"Error: Other."<<endl;
  214. return false;
  215. }
  216. data.solver_type = min_quad_with_fixed_data<T>::LLT;
  217. }else
  218. {
  219. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  220. cout<<" ldlt"<<endl;
  221. #endif
  222. // Either not PD or there are equality constraints
  223. SparseMatrix<T> NA;
  224. slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  225. data.NA = NA;
  226. // Ideally we'd use LDLT but Eigen doesn't support positive semi-definite
  227. // matrices:
  228. // http://forum.kde.org/viewtopic.php?f=74&t=106962&p=291990#p291990
  229. if(data.Auu_sym && false)
  230. {
  231. data.ldlt.compute(NA);
  232. switch(data.ldlt.info())
  233. {
  234. case Eigen::Success:
  235. break;
  236. case Eigen::NumericalIssue:
  237. cerr<<"Error: Numerical issue."<<endl;
  238. return false;
  239. default:
  240. cerr<<"Error: Other."<<endl;
  241. return false;
  242. }
  243. data.solver_type = min_quad_with_fixed_data<T>::LDLT;
  244. }else
  245. {
  246. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  247. cout<<" lu"<<endl;
  248. #endif
  249. // Resort to LU
  250. // Bottleneck >1/2
  251. data.lu.compute(NA);
  252. //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
  253. switch(data.lu.info())
  254. {
  255. case Eigen::Success:
  256. break;
  257. case Eigen::NumericalIssue:
  258. cerr<<"Error: Numerical issue."<<endl;
  259. return false;
  260. case Eigen::InvalidInput:
  261. cerr<<"Error: Invalid Input."<<endl;
  262. return false;
  263. default:
  264. cerr<<"Error: Other."<<endl;
  265. return false;
  266. }
  267. data.solver_type = min_quad_with_fixed_data<T>::LU;
  268. }
  269. }
  270. }else
  271. {
  272. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  273. cout<<" Aeq_li=false"<<endl;
  274. #endif
  275. data.neq = neq;
  276. const int nu = data.unknown.size();
  277. //cout<<"nu: "<<nu<<endl;
  278. //cout<<"neq: "<<neq<<endl;
  279. //cout<<"nc: "<<nc<<endl;
  280. //cout<<" matrixR"<<endl;
  281. SparseMatrix<T> AeqTR,AeqTQ;
  282. AeqTR = data.AeqTQR.matrixR();
  283. // This shouldn't be necessary
  284. AeqTR.prune(static_cast<T>(0.0));
  285. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  286. cout<<" matrixQ"<<endl;
  287. #endif
  288. // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
  289. // http://forum.kde.org/viewtopic.php?f=74&t=117500
  290. AeqTQ = data.AeqTQR.matrixQ();
  291. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  292. cout<<" prune"<<endl;
  293. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  294. #endif
  295. // This shouldn't be necessary
  296. AeqTQ.prune(static_cast<T>(0.0));
  297. //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
  298. //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
  299. //cout<<" perms"<<endl;
  300. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  301. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  302. cout<<" perm"<<endl;
  303. #endif
  304. SparseMatrix<T> I(neq,neq);
  305. I.setIdentity();
  306. data.AeqTE = data.AeqTQR.colsPermutation() * I;
  307. data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
  308. assert(AeqTR.rows() == nu && "#rows in AeqTR should match #unknowns");
  309. assert(AeqTR.cols() == neq && "#cols in AeqTR should match #constraints");
  310. assert(AeqTQ.rows() == nu && "#rows in AeqTQ should match #unknowns");
  311. assert(AeqTQ.cols() == nu && "#cols in AeqTQ should match #unknowns");
  312. //cout<<" slice"<<endl;
  313. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  314. cout<<" slice"<<endl;
  315. #endif
  316. data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
  317. data.AeqTQ1T = data.AeqTQ1.transpose().eval();
  318. // ALREADY TRIM (Not 100% sure about this)
  319. data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
  320. data.AeqTR1T = data.AeqTR1.transpose().eval();
  321. //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
  322. // Null space
  323. data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
  324. data.AeqTQ2T = data.AeqTQ2.transpose().eval();
  325. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  326. cout<<" proj"<<endl;
  327. #endif
  328. // Projected hessian
  329. SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
  330. {
  331. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  332. cout<<" factorize"<<endl;
  333. #endif
  334. // QRAuu should always be PD
  335. data.llt.compute(QRAuu);
  336. switch(data.llt.info())
  337. {
  338. case Eigen::Success:
  339. break;
  340. case Eigen::NumericalIssue:
  341. cerr<<"Error: Numerical issue."<<endl;
  342. return false;
  343. default:
  344. cerr<<"Error: Other."<<endl;
  345. return false;
  346. }
  347. data.solver_type = min_quad_with_fixed_data<T>::QR_LLT;
  348. }
  349. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  350. cout<<" smash"<<endl;
  351. #endif
  352. // Known value multiplier
  353. SparseMatrix<T> Auk;
  354. slice(A,data.unknown,data.known,Auk);
  355. SparseMatrix<T> Aku;
  356. slice(A,data.known,data.unknown,Aku);
  357. SparseMatrix<T> AkuT = Aku.transpose();
  358. data.preY = Auk + AkuT;
  359. // Needed during solve
  360. data.Auu = Auu;
  361. slice(Aeq,data.known,2,data.Aeqk);
  362. assert(data.Aeqk.rows() == neq);
  363. assert(data.Aeqk.cols() == data.known.size());
  364. }
  365. return true;
  366. }
  367. template <
  368. typename T,
  369. typename DerivedB,
  370. typename DerivedY,
  371. typename DerivedBeq,
  372. typename DerivedZ,
  373. typename Derivedsol>
  374. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  375. const min_quad_with_fixed_data<T> & data,
  376. const Eigen::MatrixBase<DerivedB> & B,
  377. const Eigen::MatrixBase<DerivedY> & Y,
  378. const Eigen::MatrixBase<DerivedBeq> & Beq,
  379. Eigen::PlainObjectBase<DerivedZ> & Z,
  380. Eigen::PlainObjectBase<Derivedsol> & sol)
  381. {
  382. using namespace std;
  383. using namespace Eigen;
  384. typedef Matrix<T,Dynamic,1> VectorXT;
  385. typedef Matrix<T,Dynamic,Dynamic> MatrixXT;
  386. // number of known rows
  387. int kr = data.known.size();
  388. if(kr!=0)
  389. {
  390. assert(kr == Y.rows());
  391. }
  392. // number of columns to solve
  393. int cols = Y.cols();
  394. assert(B.cols() == 1 || B.cols() == cols);
  395. assert(Beq.size() == 0 || Beq.cols() == 1 || Beq.cols() == cols);
  396. // resize output
  397. Z.resize(data.n,cols);
  398. // Set known values
  399. for(int i = 0;i < kr;i++)
  400. {
  401. for(int j = 0;j < cols;j++)
  402. {
  403. Z(data.known(i),j) = Y(i,j);
  404. }
  405. }
  406. if(data.Aeq_li)
  407. {
  408. // number of lagrange multipliers aka linear equality constraints
  409. int neq = data.lagrange.size();
  410. // append lagrange multiplier rhs's
  411. MatrixXT BBeq(B.rows() + Beq.rows(),cols);
  412. if(B.size() > 0)
  413. {
  414. BBeq.topLeftCorner(B.rows(),cols) = B.replicate(1,B.cols()==cols?1:cols);
  415. }
  416. if(Beq.size() > 0)
  417. {
  418. BBeq.bottomLeftCorner(Beq.rows(),cols) = -2.0*Beq.replicate(1,Beq.cols()==cols?1:cols);
  419. }
  420. // Build right hand side
  421. MatrixXT BBequlcols;
  422. igl::slice(BBeq,data.unknown_lagrange,1,BBequlcols);
  423. MatrixXT NB;
  424. if(kr == 0)
  425. {
  426. NB = BBequlcols;
  427. }else
  428. {
  429. NB = data.preY * Y + BBequlcols;
  430. }
  431. //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
  432. //cout<<matlab_format(NB,"NB")<<endl;
  433. switch(data.solver_type)
  434. {
  435. case igl::min_quad_with_fixed_data<T>::LLT:
  436. sol = data.llt.solve(NB);
  437. break;
  438. case igl::min_quad_with_fixed_data<T>::LDLT:
  439. sol = data.ldlt.solve(NB);
  440. break;
  441. case igl::min_quad_with_fixed_data<T>::LU:
  442. // Not a bottleneck
  443. sol = data.lu.solve(NB);
  444. break;
  445. default:
  446. cerr<<"Error: invalid solver type"<<endl;
  447. return false;
  448. }
  449. //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
  450. // Now sol contains sol/-0.5
  451. sol *= -0.5;
  452. // Now sol contains solution
  453. // Place solution in Z
  454. for(int i = 0;i<(sol.rows()-neq);i++)
  455. {
  456. for(int j = 0;j<sol.cols();j++)
  457. {
  458. Z(data.unknown_lagrange(i),j) = sol(i,j);
  459. }
  460. }
  461. }else
  462. {
  463. assert(data.solver_type == min_quad_with_fixed_data<T>::QR_LLT);
  464. MatrixXT eff_Beq;
  465. // Adjust Aeq rhs to include known parts
  466. eff_Beq =
  467. //data.AeqTQR.colsPermutation().transpose() * (-data.Aeqk * Y + Beq);
  468. data.AeqTET * (-data.Aeqk * Y + Beq.replicate(1,Beq.cols()==cols?1:cols));
  469. // Where did this -0.5 come from? Probably the same place as above.
  470. MatrixXT Bu;
  471. slice(B,data.unknown,1,Bu);
  472. MatrixXT NB;
  473. NB = -0.5*(Bu.replicate(1,B.cols()==cols?1:cols) + data.preY * Y);
  474. // Trim eff_Beq
  475. const int nc = data.AeqTQR.rank();
  476. const int neq = Beq.rows();
  477. eff_Beq = eff_Beq.topLeftCorner(nc,cols).eval();
  478. data.AeqTR1T.template triangularView<Lower>().solveInPlace(eff_Beq);
  479. // Now eff_Beq = (data.AeqTR1T \ (data.AeqTET * (-data.Aeqk * Y + Beq)))
  480. MatrixXT lambda_0;
  481. lambda_0 = data.AeqTQ1 * eff_Beq;
  482. //cout<<matlab_format(lambda_0,"lambda_0")<<endl;
  483. MatrixXT QRB;
  484. QRB = -data.AeqTQ2T * (data.Auu * lambda_0) + data.AeqTQ2T * NB;
  485. Derivedsol lambda;
  486. lambda = data.llt.solve(QRB);
  487. // prepare output
  488. Derivedsol solu;
  489. solu = data.AeqTQ2 * lambda + lambda_0;
  490. // http://www.math.uh.edu/~rohop/fall_06/Chapter3.pdf
  491. Derivedsol solLambda;
  492. {
  493. Derivedsol temp1,temp2;
  494. temp1 = (data.AeqTQ1T * NB - data.AeqTQ1T * data.Auu * solu);
  495. data.AeqTR1.template triangularView<Upper>().solveInPlace(temp1);
  496. //cout<<matlab_format(temp1,"temp1")<<endl;
  497. temp2 = Derivedsol::Zero(neq,cols);
  498. temp2.topLeftCorner(nc,cols) = temp1;
  499. //solLambda = data.AeqTQR.colsPermutation() * temp2;
  500. solLambda = data.AeqTE * temp2;
  501. }
  502. // sol is [Z(unknown);Lambda]
  503. assert(data.unknown.size() == solu.rows());
  504. assert(cols == solu.cols());
  505. assert(data.neq == neq);
  506. assert(data.neq == solLambda.rows());
  507. assert(cols == solLambda.cols());
  508. sol.resize(data.unknown.size()+data.neq,cols);
  509. sol.block(0,0,solu.rows(),solu.cols()) = solu;
  510. sol.block(solu.rows(),0,solLambda.rows(),solLambda.cols()) = solLambda;
  511. for(int u = 0;u<data.unknown.size();u++)
  512. {
  513. for(int j = 0;j<Z.cols();j++)
  514. {
  515. Z(data.unknown(u),j) = solu(u,j);
  516. }
  517. }
  518. }
  519. return true;
  520. }
  521. template <
  522. typename T,
  523. typename DerivedB,
  524. typename DerivedY,
  525. typename DerivedBeq,
  526. typename DerivedZ>
  527. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  528. const min_quad_with_fixed_data<T> & data,
  529. const Eigen::MatrixBase<DerivedB> & B,
  530. const Eigen::MatrixBase<DerivedY> & Y,
  531. const Eigen::MatrixBase<DerivedBeq> & Beq,
  532. Eigen::PlainObjectBase<DerivedZ> & Z)
  533. {
  534. Eigen::Matrix<typename DerivedZ::Scalar, Eigen::Dynamic, Eigen::Dynamic> sol;
  535. return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
  536. }
  537. template <
  538. typename T,
  539. typename Derivedknown,
  540. typename DerivedB,
  541. typename DerivedY,
  542. typename DerivedBeq,
  543. typename DerivedZ>
  544. IGL_INLINE bool igl::min_quad_with_fixed(
  545. const Eigen::SparseMatrix<T>& A,
  546. const Eigen::MatrixBase<DerivedB> & B,
  547. const Eigen::MatrixBase<Derivedknown> & known,
  548. const Eigen::MatrixBase<DerivedY> & Y,
  549. const Eigen::SparseMatrix<T>& Aeq,
  550. const Eigen::MatrixBase<DerivedBeq> & Beq,
  551. const bool pd,
  552. Eigen::PlainObjectBase<DerivedZ> & Z)
  553. {
  554. min_quad_with_fixed_data<T> data;
  555. if(!min_quad_with_fixed_precompute(A,known,Aeq,pd,data))
  556. {
  557. return false;
  558. }
  559. return min_quad_with_fixed_solve(data,B,Y,Beq,Z);
  560. }
  561. template <typename Scalar, int n, int m, bool Hpd>
  562. IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
  563. const Eigen::Matrix<Scalar,n,n> & H,
  564. const Eigen::Matrix<Scalar,n,1> & f,
  565. const Eigen::Array<bool,n,1> & k,
  566. const Eigen::Matrix<Scalar,n,1> & bc,
  567. const Eigen::Matrix<Scalar,m,n> & A,
  568. const Eigen::Matrix<Scalar,m,1> & b)
  569. {
  570. //printf("igl::min_quad_with_fixed<S,n,m,Hpd>\n");
  571. //std::cout<<igl::matlab_format(H,"H")<<std::endl;
  572. //std::cout<<igl::matlab_format(f,"f")<<std::endl;
  573. //std::cout<<igl::matlab_format(k,"k")<<std::endl;
  574. //std::cout<<igl::matlab_format(bc,"bc")<<std::endl;
  575. //std::cout<<igl::matlab_format(A,"A")<<std::endl;
  576. //std::cout<<igl::matlab_format(b,"b")<<std::endl;
  577. const auto dyn_n = n == Eigen::Dynamic ? H.rows() : n;
  578. const auto dyn_m = m == Eigen::Dynamic ? A.rows() : m;
  579. const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  580. const auto dyn_nn = nn == Eigen::Dynamic ? dyn_n+dyn_m : nn;
  581. if(dyn_m == 0)
  582. {
  583. return igl::min_quad_with_fixed<Scalar,n,Hpd>(H,f,k,bc);
  584. }
  585. // min_x ½ xᵀ H x + xᵀ f subject to A x = b and x(k) = bc(k)
  586. // let zᵀ = [xᵀ λᵀ]
  587. // min_z ½ zᵀ [H Aᵀ;A 0] z + zᵀ [f;-b] z(k) = bc(k)
  588. const auto make_HH = [&]()
  589. {
  590. Eigen::Matrix<Scalar,nn,nn> HH =
  591. Eigen::Matrix<Scalar,nn,nn>::Zero(dyn_nn,dyn_nn);
  592. HH.topLeftCorner(dyn_n,dyn_n) = H;
  593. HH.bottomLeftCorner(dyn_m,dyn_n) = A;
  594. HH.topRightCorner(dyn_n,dyn_m) = A.transpose();
  595. return HH;
  596. };
  597. const Eigen::Matrix<Scalar,nn,nn> HH = make_HH();
  598. const auto make_ff = [&]()
  599. {
  600. Eigen::Matrix<Scalar,nn,1> ff(dyn_nn);
  601. ff.head(dyn_n) = f;
  602. ff.tail(dyn_m) = -b;
  603. return ff;
  604. };
  605. const Eigen::Matrix<Scalar,nn,1> ff = make_ff();
  606. const auto make_kk = [&]()
  607. {
  608. Eigen::Array<bool,nn,1> kk =
  609. Eigen::Array<bool,nn,1>::Constant(dyn_nn,1,false);
  610. kk.head(dyn_n) = k;
  611. return kk;
  612. };
  613. const Eigen::Array<bool,nn,1> kk = make_kk();
  614. const auto make_bcbc= [&]()
  615. {
  616. Eigen::Matrix<Scalar,nn,1> bcbc(dyn_nn);
  617. bcbc.head(dyn_n) = bc;
  618. return bcbc;
  619. };
  620. const Eigen::Matrix<double,nn,1> bcbc = make_bcbc();
  621. const Eigen::Matrix<Scalar,nn,1> xx =
  622. min_quad_with_fixed<Scalar,nn,false>(HH,ff,kk,bcbc);
  623. //std::cout<<igl::matlab_format(xx.head(dyn_n).eval(),"x")<<std::endl;
  624. return xx.head(dyn_n);
  625. }
  626. template <typename Scalar, int n, bool Hpd>
  627. IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
  628. const Eigen::Matrix<Scalar,n,n> & H,
  629. const Eigen::Matrix<Scalar,n,1> & f,
  630. const Eigen::Array<bool,n,1> & k,
  631. const Eigen::Matrix<Scalar,n,1> & bc)
  632. {
  633. //printf("igl::min_quad_with_fixed<S,n,Hpd>\n");
  634. //std::cout<<igl::matlab_format(H,"H")<<std::endl;
  635. //std::cout<<igl::matlab_format(f,"f")<<std::endl;
  636. //std::cout<<igl::matlab_format(k,"k")<<std::endl;
  637. //std::cout<<igl::matlab_format(bc,"bc")<<std::endl;
  638. assert(H.isApprox(H.transpose(),1e-7));
  639. assert(H.rows() == H.cols());
  640. assert(H.rows() == f.size());
  641. assert(H.rows() == k.size());
  642. assert(H.rows() == bc.size());
  643. const auto kcount = k.count();
  644. // Everything fixed
  645. if(kcount == (Eigen::Dynamic?H.rows():n))
  646. {
  647. //std::cout<<igl::matlab_format(bc,"x")<<std::endl;
  648. return bc;
  649. }
  650. // Nothing fixed
  651. if(kcount == 0)
  652. {
  653. // avoid function call
  654. typedef Eigen::Matrix<Scalar,n,n> MatrixSn;
  655. typedef typename
  656. std::conditional<Hpd,Eigen::LLT<MatrixSn>,Eigen::CompleteOrthogonalDecomposition<MatrixSn>>::type
  657. Solver;
  658. //std::cout<<igl::matlab_format( Solver(H).solve(-f).eval() ,"x")<<std::endl;
  659. return Solver(H).solve(-f);
  660. }
  661. // All-but-one fixed
  662. if( (Eigen::Dynamic?H.rows():n)-kcount == 1)
  663. {
  664. // which one is not fixed?
  665. int u = -1;
  666. for(int i=0;i<k.size();i++){ if(!k(i)){ u=i; break; } }
  667. assert(u>=0);
  668. // min ½ x(u) Huu x(u) + x(u)(fu + H(u,k)bc(k))
  669. // Huu x(u) = -(fu + H(u,k) bc(k))
  670. // x(u) = (-fu + ∑ -Huj bcj)/Huu
  671. Eigen::Matrix<Scalar,n,1> x = bc;
  672. x(u) = -f(u);
  673. for(int i=0;i<k.size();i++){ if(i!=u){ x(u)-=bc(i)*H(i,u); } }
  674. x(u) /= H(u,u);
  675. //std::cout<<igl::matlab_format(x,"x")<<std::endl;
  676. return x;
  677. }
  678. // Is there a smart template way to do this?
  679. switch(kcount)
  680. {
  681. case 0: assert(false && "Handled above."); return Eigen::Matrix<Scalar,n,1>();
  682. // maxi=16;for i=1:maxi;fprintf(' case %d:\n {\n const bool D = (n-%d<=0)||(%d>=n)||(n>%d);\n return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:%d,Hpd>(H,f,k,bc);\n }\n',[i i i maxi i]);end
  683. case 1:
  684. {
  685. const bool D = (n-1<=0)||(1>=n)||(n>16);
  686. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:1,Hpd>(H,f,k,bc);
  687. }
  688. case 2:
  689. {
  690. const bool D = (n-2<=0)||(2>=n)||(n>16);
  691. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:2,Hpd>(H,f,k,bc);
  692. }
  693. case 3:
  694. {
  695. const bool D = (n-3<=0)||(3>=n)||(n>16);
  696. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:3,Hpd>(H,f,k,bc);
  697. }
  698. case 4:
  699. {
  700. const bool D = (n-4<=0)||(4>=n)||(n>16);
  701. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:4,Hpd>(H,f,k,bc);
  702. }
  703. case 5:
  704. {
  705. const bool D = (n-5<=0)||(5>=n)||(n>16);
  706. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:5,Hpd>(H,f,k,bc);
  707. }
  708. case 6:
  709. {
  710. const bool D = (n-6<=0)||(6>=n)||(n>16);
  711. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:6,Hpd>(H,f,k,bc);
  712. }
  713. case 7:
  714. {
  715. const bool D = (n-7<=0)||(7>=n)||(n>16);
  716. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:7,Hpd>(H,f,k,bc);
  717. }
  718. case 8:
  719. {
  720. const bool D = (n-8<=0)||(8>=n)||(n>16);
  721. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:8,Hpd>(H,f,k,bc);
  722. }
  723. case 9:
  724. {
  725. const bool D = (n-9<=0)||(9>=n)||(n>16);
  726. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:9,Hpd>(H,f,k,bc);
  727. }
  728. case 10:
  729. {
  730. const bool D = (n-10<=0)||(10>=n)||(n>16);
  731. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:10,Hpd>(H,f,k,bc);
  732. }
  733. case 11:
  734. {
  735. const bool D = (n-11<=0)||(11>=n)||(n>16);
  736. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:11,Hpd>(H,f,k,bc);
  737. }
  738. case 12:
  739. {
  740. const bool D = (n-12<=0)||(12>=n)||(n>16);
  741. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:12,Hpd>(H,f,k,bc);
  742. }
  743. case 13:
  744. {
  745. const bool D = (n-13<=0)||(13>=n)||(n>16);
  746. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:13,Hpd>(H,f,k,bc);
  747. }
  748. case 14:
  749. {
  750. const bool D = (n-14<=0)||(14>=n)||(n>16);
  751. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:14,Hpd>(H,f,k,bc);
  752. }
  753. case 15:
  754. {
  755. const bool D = (n-15<=0)||(15>=n)||(n>16);
  756. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:15,Hpd>(H,f,k,bc);
  757. }
  758. case 16:
  759. {
  760. const bool D = (n-16<=0)||(16>=n)||(n>16);
  761. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:16,Hpd>(H,f,k,bc);
  762. }
  763. default:
  764. return min_quad_with_fixed<Scalar,Eigen::Dynamic,Eigen::Dynamic,Hpd>(H,f,k,bc);
  765. }
  766. }
  767. IGL_INLINE template <typename Scalar, int n, int kcount, bool Hpd>
  768. Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
  769. const Eigen::Matrix<Scalar,n,n> & H,
  770. const Eigen::Matrix<Scalar,n,1> & f,
  771. const Eigen::Array<bool,n,1> & k,
  772. const Eigen::Matrix<Scalar,n,1> & bc)
  773. {
  774. //printf("igl::min_quad_with_fixed<S,n,kcount,Hpd>\n");
  775. // 0 and n should be handle outside this function
  776. static_assert(kcount==Eigen::Dynamic || kcount>0 ,"");
  777. static_assert(kcount==Eigen::Dynamic || kcount<n ,"");
  778. const int ucount = n==Eigen::Dynamic ? Eigen::Dynamic : n-kcount;
  779. static_assert(kcount==Eigen::Dynamic || ucount+kcount == n ,"");
  780. static_assert((n==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
  781. static_assert((kcount==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
  782. assert((n==Eigen::Dynamic) || n == H.rows());
  783. assert((kcount==Eigen::Dynamic) || kcount == k.count());
  784. typedef Eigen::Matrix<Scalar,ucount,ucount> MatrixSuu;
  785. typedef Eigen::Matrix<Scalar,ucount,kcount> MatrixSuk;
  786. typedef Eigen::Matrix<Scalar,n,1> VectorSn;
  787. typedef Eigen::Matrix<Scalar,ucount,1> VectorSu;
  788. typedef Eigen::Matrix<Scalar,kcount,1> VectorSk;
  789. const auto dyn_n = n==Eigen::Dynamic ? H.rows() : n;
  790. const auto dyn_kcount = kcount==Eigen::Dynamic ? k.count() : kcount;
  791. const auto dyn_ucount = ucount==Eigen::Dynamic ? dyn_n- dyn_kcount : ucount;
  792. // For ucount==2 or kcount==2 this calls the coefficient initiliazer rather
  793. // than the size initilizer, but I guess that's ok.
  794. MatrixSuu Huu(dyn_ucount,dyn_ucount);
  795. MatrixSuk Huk(dyn_ucount,dyn_kcount);
  796. VectorSu mrhs(dyn_ucount);
  797. VectorSk bck(dyn_kcount);
  798. {
  799. int ui = 0;
  800. int ki = 0;
  801. for(int i = 0;i<dyn_n;i++)
  802. {
  803. if(k(i))
  804. {
  805. bck(ki) = bc(i);
  806. ki++;
  807. }else
  808. {
  809. mrhs(ui) = f(i);
  810. int uj = 0;
  811. int kj = 0;
  812. for(int j = 0;j<dyn_n;j++)
  813. {
  814. if(k(j))
  815. {
  816. Huk(ui,kj) = H(i,j);
  817. kj++;
  818. }else
  819. {
  820. Huu(ui,uj) = H(i,j);
  821. uj++;
  822. }
  823. }
  824. ui++;
  825. }
  826. }
  827. }
  828. mrhs += Huk * bck;
  829. typedef typename
  830. std::conditional<Hpd,
  831. Eigen::LLT<MatrixSuu>,
  832. // LDLT should be faster for indefinite problems but already found some
  833. // cases where it was too inaccurate when called via quadprog_primal.
  834. // Ideally this function takes LLT,LDLT, or
  835. // CompleteOrthogonalDecomposition as a template parameter. "template
  836. // template" parameters did work because LLT,LDLT have different number of
  837. // template parameters from CompleteOrthogonalDecomposition. Perhaps
  838. // there's a way to take advantage of LLT and LDLT's default template
  839. // parameters (I couldn't figure out how).
  840. Eigen::CompleteOrthogonalDecomposition<MatrixSuu>>::type
  841. Solver;
  842. VectorSu xu = Solver(Huu).solve(-mrhs);
  843. VectorSn x(dyn_n);
  844. {
  845. int ui = 0;
  846. int ki = 0;
  847. for(int i = 0;i<dyn_n;i++)
  848. {
  849. if(k(i))
  850. {
  851. x(i) = bck(ki);
  852. ki++;
  853. }else
  854. {
  855. x(i) = xu(ui);
  856. ui++;
  857. }
  858. }
  859. }
  860. //std::cout<<igl::matlab_format(x,"x")<<std::endl;
  861. return x;
  862. }
  863. #ifdef IGL_STATIC_LIBRARY
  864. // Explicit template instantiation
  865. // generated by autoexplicit.sh
  866. template Eigen::Matrix<double, -1, 1, ((Eigen::StorageOptions)0) | ((((-1) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), -1, 1> igl::min_quad_with_fixed<double, -1, -1, true>(Eigen::Matrix<double, -1, -1, ((Eigen::StorageOptions)0) | ((((-1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)1) : ((((-1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), -1, -1> const&, Eigen::Matrix<double, -1, 1, ((Eigen::StorageOptions)0) | ((((-1) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), -1, 1> const&, Eigen::Array<bool, -1, 1, ((Eigen::StorageOptions)0) | ((((-1) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), -1, 1> const&, Eigen::Matrix<double, -1, 1, ((Eigen::StorageOptions)0) | ((((-1) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), -1, 1> const&, Eigen::Matrix<double, -1, -1, ((Eigen::StorageOptions)0) | ((((-1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)1) : ((((-1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), -1, -1> const&, Eigen::Matrix<double, -1, 1, ((Eigen::StorageOptions)0) | ((((-1) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((-1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), -1, 1> const&);
  867. // generated by autoexplicit.sh
  868. template Eigen::Matrix<double, 8, 1, ((Eigen::StorageOptions)0) | ((((8) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((8) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 8, 1> igl::min_quad_with_fixed<double, 8, 2, true>(Eigen::Matrix<double, 8, 8, ((Eigen::StorageOptions)0) | ((((8) == (1)) && ((8) != (1))) ? ((Eigen::StorageOptions)1) : ((((8) == (1)) && ((8) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 8, 8> const&, Eigen::Matrix<double, 8, 1, ((Eigen::StorageOptions)0) | ((((8) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((8) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 8, 1> const&, Eigen::Array<bool, 8, 1, ((Eigen::StorageOptions)0) | ((((8) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((8) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 8, 1> const&, Eigen::Matrix<double, 8, 1, ((Eigen::StorageOptions)0) | ((((8) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((8) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 8, 1> const&, Eigen::Matrix<double, 2, 8, ((Eigen::StorageOptions)0) | ((((2) == (1)) && ((8) != (1))) ? ((Eigen::StorageOptions)1) : ((((8) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 2, 8> const&, Eigen::Matrix<double, 2, 1, ((Eigen::StorageOptions)0) | ((((2) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 2, 1> const&);
  869. // generated by autoexplicit.sh
  870. template Eigen::Matrix<double, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> igl::min_quad_with_fixed<double, 3, 1, true>(Eigen::Matrix<double, 3, 3, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)1) : ((((3) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 3> const&, Eigen::Matrix<double, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> const&, Eigen::Array<bool, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> const&, Eigen::Matrix<double, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> const&, Eigen::Matrix<double, 1, 3, ((Eigen::StorageOptions)0) | ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)1) : ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 1, 3> const&, Eigen::Matrix<double, 1, 1, ((Eigen::StorageOptions)0) | ((((1) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 1, 1> const&);
  871. template Eigen::Matrix<double, 2, 1, ((Eigen::StorageOptions)0) | ((((2) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 2, 1> igl::min_quad_with_fixed<double, 2, 0, true>(Eigen::Matrix<double, 2, 2, ((Eigen::StorageOptions)0) | ((((2) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)1) : ((((2) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 2, 2> const&, Eigen::Matrix<double, 2, 1, ((Eigen::StorageOptions)0) | ((((2) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 2, 1> const&, Eigen::Array<bool, 2, 1, ((Eigen::StorageOptions)0) | ((((2) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 2, 1> const&, Eigen::Matrix<double, 2, 1, ((Eigen::StorageOptions)0) | ((((2) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 2, 1> const&, Eigen::Matrix<double, 0, 2, ((Eigen::StorageOptions)0) | ((((0) == (1)) && ((2) != (1))) ? ((Eigen::StorageOptions)1) : ((((2) == (1)) && ((0) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 0, 2> const&, Eigen::Matrix<double, 0, 1, ((Eigen::StorageOptions)0) | ((((0) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((0) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 0, 1> const&);
  872. template Eigen::Matrix<double, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> igl::min_quad_with_fixed<double, 3, 0, true>(Eigen::Matrix<double, 3, 3, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)1) : ((((3) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 3> const&, Eigen::Matrix<double, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> const&, Eigen::Array<bool, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> const&, Eigen::Matrix<double, 3, 1, ((Eigen::StorageOptions)0) | ((((3) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 3, 1> const&, Eigen::Matrix<double, 0, 3, ((Eigen::StorageOptions)0) | ((((0) == (1)) && ((3) != (1))) ? ((Eigen::StorageOptions)1) : ((((3) == (1)) && ((0) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 0, 3> const&, Eigen::Matrix<double, 0, 1, ((Eigen::StorageOptions)0) | ((((0) == (1)) && ((1) != (1))) ? ((Eigen::StorageOptions)1) : ((((1) == (1)) && ((0) != (1))) ? ((Eigen::StorageOptions)0) : ((Eigen::StorageOptions)0))), 0, 1> const&);
  873. 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> >&);
  874. #if EIGEN_VERSION_AT_LEAST(3,3,0)
  875. #else
  876. 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> >&);
  877. 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> >&);
  878. #endif
  879. 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> >&);
  880. 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> >&);
  881. 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>&);
  882. 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> >&);
  883. 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> >&);
  884. 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>&);
  885. 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> >&);
  886. 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> >&);
  887. 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> >&);
  888. 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> >&);
  889. 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> >&);
  890. 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> >&);
  891. #endif