min_quad_with_fixed.impl.h 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886
  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. #pragma once
  9. #include "min_quad_with_fixed.h"
  10. #include "slice.h"
  11. #include "is_symmetric.h"
  12. #include "find.h"
  13. #include "sparse.h"
  14. #include "repmat.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 "matlab_format.h"
  24. #include <type_traits>
  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::MatrixBase<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(known.cols() == 1 && "known should be a vector");
  52. assert(A.rows() == n && "A should be square");
  53. assert(A.cols() == n && "A should be square");
  54. // number of known rows
  55. int kr = known.size();
  56. assert((kr == 0 || known.minCoeff() >= 0)&& "known indices should be in [0,n)");
  57. assert((kr == 0 || known.maxCoeff() < n) && "known indices should be in [0,n)");
  58. assert(neq <= n && "Number of equality constraints should be less than DOFs");
  59. // cache known
  60. // FIXME: This is *NOT* generic and introduces a copy.
  61. data.known = known.template cast<int>();
  62. // get list of unknown indices
  63. data.unknown.resize(n-kr);
  64. std::vector<bool> unknown_mask;
  65. unknown_mask.resize(n,true);
  66. for(int i = 0;i<kr;i++)
  67. {
  68. unknown_mask[known(i, 0)] = false;
  69. }
  70. int u = 0;
  71. for(int i = 0;i<n;i++)
  72. {
  73. if(unknown_mask[i])
  74. {
  75. data.unknown(u) = i;
  76. u++;
  77. }
  78. }
  79. // get list of lagrange multiplier indices
  80. data.lagrange.resize(neq);
  81. for(int i = 0;i<neq;i++)
  82. {
  83. data.lagrange(i) = n + i;
  84. }
  85. // cache unknown followed by lagrange indices
  86. data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
  87. // Would like to do:
  88. //data.unknown_lagrange << data.unknown, data.lagrange;
  89. // but Eigen can't handle empty vectors in comma initialization
  90. // https://forum.kde.org/viewtopic.php?f=74&t=107974&p=364947#p364947
  91. if(data.unknown.size() > 0)
  92. {
  93. data.unknown_lagrange.head(data.unknown.size()) = data.unknown;
  94. }
  95. if(data.lagrange.size() > 0)
  96. {
  97. data.unknown_lagrange.tail(data.lagrange.size()) = data.lagrange;
  98. }
  99. SparseMatrix<T> Auu;
  100. slice(A,data.unknown,data.unknown,Auu);
  101. assert(Auu.size() != 0 && Auu.rows() > 0 && "There should be at least one unknown.");
  102. // Positive definiteness is *not* determined, rather it is given as a
  103. // parameter
  104. data.Auu_pd = pd;
  105. if(data.Auu_pd)
  106. {
  107. // PD implies symmetric
  108. data.Auu_sym = true;
  109. // This is an annoying assertion unless EPS can be chosen in a nicer way.
  110. //assert(is_symmetric(Auu,EPS<T>()));
  111. assert(is_symmetric(Auu,1.0) &&
  112. "Auu should be symmetric if positive definite");
  113. }else
  114. {
  115. // determine if A(unknown,unknown) is symmetric and/or positive definite
  116. VectorXi AuuI,AuuJ;
  117. Matrix<T,Eigen::Dynamic,Eigen::Dynamic> AuuV;
  118. find(Auu,AuuI,AuuJ,AuuV);
  119. data.Auu_sym = is_symmetric(Auu,EPS<T>()*AuuV.maxCoeff());
  120. }
  121. // Determine number of linearly independent constraints
  122. int nc = 0;
  123. if(neq>0)
  124. {
  125. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  126. cout<<" qr"<<endl;
  127. #endif
  128. // QR decomposition to determine row rank in Aequ
  129. slice(Aeq,data.unknown,2,data.Aequ);
  130. assert(data.Aequ.rows() == neq &&
  131. "#Rows in Aequ should match #constraints");
  132. assert(data.Aequ.cols() == data.unknown.size() &&
  133. "#cols in Aequ should match #unknowns");
  134. data.AeqTQR.compute(data.Aequ.transpose().eval());
  135. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  136. //cout<<endl<<matlab_format(SparseMatrix<T>(data.Aequ.transpose().eval()),"AeqT")<<endl<<endl;
  137. #endif
  138. switch(data.AeqTQR.info())
  139. {
  140. case Eigen::Success:
  141. break;
  142. case Eigen::NumericalIssue:
  143. cerr<<"Error: Numerical issue."<<endl;
  144. return false;
  145. case Eigen::InvalidInput:
  146. cerr<<"Error: Invalid input."<<endl;
  147. return false;
  148. default:
  149. cerr<<"Error: Other."<<endl;
  150. return false;
  151. }
  152. nc = data.AeqTQR.rank();
  153. assert(nc<=neq &&
  154. "Rank of reduced constraints should be <= #original constraints");
  155. data.Aeq_li = nc == neq;
  156. //cout<<"data.Aeq_li: "<<data.Aeq_li<<endl;
  157. }else
  158. {
  159. data.Aeq_li = true;
  160. }
  161. if(data.Aeq_li)
  162. {
  163. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  164. cout<<" Aeq_li=true"<<endl;
  165. #endif
  166. // Append lagrange multiplier quadratic terms
  167. SparseMatrix<T> new_A;
  168. SparseMatrix<T> AeqT = Aeq.transpose();
  169. SparseMatrix<T> Z(neq,neq);
  170. // This is a bit slower. But why isn't cat fast?
  171. new_A = cat(1, cat(2, A, AeqT ),
  172. cat(2, Aeq, Z ));
  173. // precompute RHS builders
  174. if(kr > 0)
  175. {
  176. SparseMatrix<T> Aulk,Akul;
  177. // Slow
  178. slice(new_A,data.unknown_lagrange,data.known,Aulk);
  179. //// This doesn't work!!!
  180. //data.preY = Aulk + Akul.transpose();
  181. // Slow
  182. if(data.Auu_sym)
  183. {
  184. data.preY = Aulk*2;
  185. }else
  186. {
  187. slice(new_A,data.known,data.unknown_lagrange,Akul);
  188. SparseMatrix<T> AkulT = Akul.transpose();
  189. data.preY = Aulk + AkulT;
  190. }
  191. }else
  192. {
  193. data.preY.resize(data.unknown_lagrange.size(),0);
  194. }
  195. // Positive definite and no equality constraints (Positive definiteness
  196. // implies symmetric)
  197. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  198. cout<<" factorize"<<endl;
  199. #endif
  200. if(data.Auu_pd && neq == 0)
  201. {
  202. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  203. cout<<" llt"<<endl;
  204. #endif
  205. data.llt.compute(Auu);
  206. switch(data.llt.info())
  207. {
  208. case Eigen::Success:
  209. break;
  210. case Eigen::NumericalIssue:
  211. cerr<<"Error: Numerical issue."<<endl;
  212. return false;
  213. default:
  214. cerr<<"Error: Other."<<endl;
  215. return false;
  216. }
  217. data.solver_type = min_quad_with_fixed_data<T>::LLT;
  218. }else
  219. {
  220. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  221. cout<<" ldlt/lu"<<endl;
  222. #endif
  223. // Either not PD or there are equality constraints
  224. SparseMatrix<T> NA;
  225. slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  226. data.NA = NA;
  227. if(data.Auu_pd)
  228. {
  229. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  230. cout<<" ldlt"<<endl;
  231. #endif
  232. data.ldlt.compute(NA);
  233. switch(data.ldlt.info())
  234. {
  235. case Eigen::Success:
  236. break;
  237. case Eigen::NumericalIssue:
  238. cerr<<"Error: Numerical issue."<<endl;
  239. return false;
  240. default:
  241. cerr<<"Error: Other."<<endl;
  242. return false;
  243. }
  244. data.solver_type = min_quad_with_fixed_data<T>::LDLT;
  245. }else
  246. {
  247. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  248. cout<<" lu"<<endl;
  249. #endif
  250. // Resort to LU
  251. // Bottleneck >1/2
  252. data.lu.compute(NA);
  253. //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
  254. switch(data.lu.info())
  255. {
  256. case Eigen::Success:
  257. break;
  258. case Eigen::NumericalIssue:
  259. cerr<<"Error: Numerical issue."<<endl;
  260. return false;
  261. case Eigen::InvalidInput:
  262. cerr<<"Error: Invalid Input."<<endl;
  263. return false;
  264. default:
  265. cerr<<"Error: Other."<<endl;
  266. return false;
  267. }
  268. data.solver_type = min_quad_with_fixed_data<T>::LU;
  269. }
  270. }
  271. }else
  272. {
  273. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  274. cout<<" Aeq_li=false"<<endl;
  275. #endif
  276. data.neq = neq;
  277. const int nu = data.unknown.size();
  278. //cout<<"nu: "<<nu<<endl;
  279. //cout<<"neq: "<<neq<<endl;
  280. //cout<<"nc: "<<nc<<endl;
  281. //cout<<" matrixR"<<endl;
  282. SparseMatrix<T> AeqTR,AeqTQ;
  283. AeqTR = data.AeqTQR.matrixR();
  284. // This shouldn't be necessary
  285. AeqTR.prune(static_cast<T>(0.0));
  286. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  287. cout<<" matrixQ"<<endl;
  288. #endif
  289. // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
  290. // http://forum.kde.org/viewtopic.php?f=74&t=117500
  291. AeqTQ = data.AeqTQR.matrixQ();
  292. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  293. cout<<" prune"<<endl;
  294. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  295. #endif
  296. // This shouldn't be necessary
  297. AeqTQ.prune(static_cast<T>(0.0));
  298. //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
  299. //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
  300. //cout<<" perms"<<endl;
  301. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  302. cout<<" nnz: "<<AeqTQ.nonZeros()<<endl;
  303. cout<<" perm"<<endl;
  304. #endif
  305. SparseMatrix<T> I(neq,neq);
  306. I.setIdentity();
  307. data.AeqTE = data.AeqTQR.colsPermutation() * I;
  308. data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
  309. assert(AeqTR.rows() == nu && "#rows in AeqTR should match #unknowns");
  310. assert(AeqTR.cols() == neq && "#cols in AeqTR should match #constraints");
  311. assert(AeqTQ.rows() == nu && "#rows in AeqTQ should match #unknowns");
  312. assert(AeqTQ.cols() == nu && "#cols in AeqTQ should match #unknowns");
  313. //cout<<" slice"<<endl;
  314. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  315. cout<<" slice"<<endl;
  316. #endif
  317. data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
  318. data.AeqTQ1T = data.AeqTQ1.transpose().eval();
  319. // ALREADY TRIM (Not 100% sure about this)
  320. data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
  321. data.AeqTR1T = data.AeqTR1.transpose().eval();
  322. //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
  323. // Null space
  324. data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
  325. data.AeqTQ2T = data.AeqTQ2.transpose().eval();
  326. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  327. cout<<" proj"<<endl;
  328. #endif
  329. // Projected hessian
  330. SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
  331. {
  332. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  333. cout<<" factorize"<<endl;
  334. #endif
  335. // QRAuu should always be PD
  336. data.llt.compute(QRAuu);
  337. switch(data.llt.info())
  338. {
  339. case Eigen::Success:
  340. break;
  341. case Eigen::NumericalIssue:
  342. cerr<<"Error: Numerical issue."<<endl;
  343. return false;
  344. default:
  345. cerr<<"Error: Other."<<endl;
  346. return false;
  347. }
  348. data.solver_type = min_quad_with_fixed_data<T>::QR_LLT;
  349. }
  350. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  351. cout<<" smash"<<endl;
  352. #endif
  353. // Known value multiplier
  354. SparseMatrix<T> Auk;
  355. slice(A,data.unknown,data.known,Auk);
  356. SparseMatrix<T> Aku;
  357. slice(A,data.known,data.unknown,Aku);
  358. SparseMatrix<T> AkuT = Aku.transpose();
  359. data.preY = Auk + AkuT;
  360. // Needed during solve
  361. data.Auu = Auu;
  362. slice(Aeq,data.known,2,data.Aeqk);
  363. assert(data.Aeqk.rows() == neq);
  364. assert(data.Aeqk.cols() == data.known.size());
  365. }
  366. return true;
  367. }
  368. template <
  369. typename T,
  370. typename DerivedB,
  371. typename DerivedY,
  372. typename DerivedBeq,
  373. typename DerivedZ,
  374. typename Derivedsol>
  375. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  376. const min_quad_with_fixed_data<T> & data,
  377. const Eigen::MatrixBase<DerivedB> & B,
  378. const Eigen::MatrixBase<DerivedY> & Y,
  379. const Eigen::MatrixBase<DerivedBeq> & Beq,
  380. Eigen::PlainObjectBase<DerivedZ> & Z,
  381. Eigen::PlainObjectBase<Derivedsol> & sol)
  382. {
  383. using namespace std;
  384. using namespace Eigen;
  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. const auto dyn_n = n == Eigen::Dynamic ? H.rows() : n;
  571. const auto dyn_m = m == Eigen::Dynamic ? A.rows() : m;
  572. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  573. const auto dyn_nn = nn == Eigen::Dynamic ? dyn_n+dyn_m : nn;
  574. if(dyn_m == 0)
  575. {
  576. return igl::min_quad_with_fixed<Scalar,n,Hpd>(H,f,k,bc);
  577. }
  578. // min_x ½ xᵀ H x + xᵀ f subject to A x = b and x(k) = bc(k)
  579. // let zᵀ = [xᵀ λᵀ]
  580. // min_z ½ zᵀ [H Aᵀ;A 0] z + zᵀ [f;-b] z(k) = bc(k)
  581. const auto make_HH = [&]()
  582. {
  583. // Windows can't remember that nn is const.
  584. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  585. Eigen::Matrix<Scalar,nn,nn> HH =
  586. Eigen::Matrix<Scalar,nn,nn>::Zero(dyn_nn,dyn_nn);
  587. HH.topLeftCorner(dyn_n,dyn_n) = H;
  588. HH.bottomLeftCorner(dyn_m,dyn_n) = A;
  589. HH.topRightCorner(dyn_n,dyn_m) = A.transpose();
  590. return HH;
  591. };
  592. const Eigen::Matrix<Scalar,nn,nn> HH = make_HH();
  593. const auto make_ff = [&]()
  594. {
  595. // Windows can't remember that nn is const.
  596. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  597. Eigen::Matrix<Scalar,nn,1> ff(dyn_nn);
  598. ff.head(dyn_n) = f;
  599. ff.tail(dyn_m) = -b;
  600. return ff;
  601. };
  602. const Eigen::Matrix<Scalar,nn,1> ff = make_ff();
  603. const auto make_kk = [&]()
  604. {
  605. // Windows can't remember that nn is const.
  606. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  607. Eigen::Array<bool,nn,1> kk =
  608. Eigen::Array<bool,nn,1>::Constant(dyn_nn,1,false);
  609. kk.head(dyn_n) = k;
  610. return kk;
  611. };
  612. const Eigen::Array<bool,nn,1> kk = make_kk();
  613. const auto make_bcbc= [&]()
  614. {
  615. // Windows can't remember that nn is const.
  616. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  617. Eigen::Matrix<Scalar,nn,1> bcbc(dyn_nn);
  618. bcbc.head(dyn_n) = bc;
  619. return bcbc;
  620. };
  621. const Eigen::Matrix<Scalar,nn,1> bcbc = make_bcbc();
  622. const Eigen::Matrix<Scalar,nn,1> xx =
  623. min_quad_with_fixed<Scalar,nn,false>(HH,ff,kk,bcbc);
  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. assert(H.isApprox(H.transpose(),1e-7));
  634. assert(H.rows() == H.cols());
  635. assert(H.rows() == f.size());
  636. assert(H.rows() == k.size());
  637. assert(H.rows() == bc.size());
  638. const auto kcount = k.count();
  639. // Everything fixed
  640. if(kcount == (Eigen::Dynamic?H.rows():n))
  641. {
  642. return bc;
  643. }
  644. // Nothing fixed
  645. if(kcount == 0)
  646. {
  647. // avoid function call
  648. typedef Eigen::Matrix<Scalar,n,n> MatrixSn;
  649. typedef typename
  650. std::conditional<Hpd,Eigen::LLT<MatrixSn>,Eigen::CompleteOrthogonalDecomposition<MatrixSn>>::type
  651. Solver;
  652. return Solver(H).solve(-f);
  653. }
  654. // All-but-one fixed
  655. if( (Eigen::Dynamic?H.rows():n)-kcount == 1)
  656. {
  657. // which one is not fixed?
  658. int u = -1;
  659. for(int i=0;i<k.size();i++){ if(!k(i)){ u=i; break; } }
  660. assert(u>=0);
  661. // min ½ x(u) Huu x(u) + x(u)(fu + H(u,k)bc(k))
  662. // Huu x(u) = -(fu + H(u,k) bc(k))
  663. // x(u) = (-fu + ∑ -Huj bcj)/Huu
  664. Eigen::Matrix<Scalar,n,1> x = bc;
  665. x(u) = -f(u);
  666. for(int i=0;i<k.size();i++){ if(i!=u){ x(u)-=bc(i)*H(i,u); } }
  667. x(u) /= H(u,u);
  668. return x;
  669. }
  670. // Alec: Is there a smart template way to do this?
  671. // jdumas: I guess you could do a templated for-loop starting from 16, and
  672. // dispatching to the appropriate templated function when the argument matches
  673. // (with a fallback to the dynamic version). Cf this example:
  674. // https://gist.github.com/disconnect3d/13c2d035bb31b244df14
  675. switch(kcount)
  676. {
  677. case 0: assert(false && "Handled above."); return Eigen::Matrix<Scalar,n,1>();
  678. // % Matlibberish for generating these case statements:
  679. // 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
  680. case 1:
  681. {
  682. const bool D = (n-1<=0)||(1>=n)||(n>16);
  683. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:1,Hpd>(H,f,k,bc);
  684. }
  685. case 2:
  686. {
  687. const bool D = (n-2<=0)||(2>=n)||(n>16);
  688. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:2,Hpd>(H,f,k,bc);
  689. }
  690. case 3:
  691. {
  692. const bool D = (n-3<=0)||(3>=n)||(n>16);
  693. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:3,Hpd>(H,f,k,bc);
  694. }
  695. case 4:
  696. {
  697. const bool D = (n-4<=0)||(4>=n)||(n>16);
  698. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:4,Hpd>(H,f,k,bc);
  699. }
  700. case 5:
  701. {
  702. const bool D = (n-5<=0)||(5>=n)||(n>16);
  703. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:5,Hpd>(H,f,k,bc);
  704. }
  705. case 6:
  706. {
  707. const bool D = (n-6<=0)||(6>=n)||(n>16);
  708. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:6,Hpd>(H,f,k,bc);
  709. }
  710. case 7:
  711. {
  712. const bool D = (n-7<=0)||(7>=n)||(n>16);
  713. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:7,Hpd>(H,f,k,bc);
  714. }
  715. case 8:
  716. {
  717. const bool D = (n-8<=0)||(8>=n)||(n>16);
  718. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:8,Hpd>(H,f,k,bc);
  719. }
  720. case 9:
  721. {
  722. const bool D = (n-9<=0)||(9>=n)||(n>16);
  723. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:9,Hpd>(H,f,k,bc);
  724. }
  725. case 10:
  726. {
  727. const bool D = (n-10<=0)||(10>=n)||(n>16);
  728. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:10,Hpd>(H,f,k,bc);
  729. }
  730. case 11:
  731. {
  732. const bool D = (n-11<=0)||(11>=n)||(n>16);
  733. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:11,Hpd>(H,f,k,bc);
  734. }
  735. case 12:
  736. {
  737. const bool D = (n-12<=0)||(12>=n)||(n>16);
  738. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:12,Hpd>(H,f,k,bc);
  739. }
  740. case 13:
  741. {
  742. const bool D = (n-13<=0)||(13>=n)||(n>16);
  743. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:13,Hpd>(H,f,k,bc);
  744. }
  745. case 14:
  746. {
  747. const bool D = (n-14<=0)||(14>=n)||(n>16);
  748. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:14,Hpd>(H,f,k,bc);
  749. }
  750. case 15:
  751. {
  752. const bool D = (n-15<=0)||(15>=n)||(n>16);
  753. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:15,Hpd>(H,f,k,bc);
  754. }
  755. case 16:
  756. {
  757. const bool D = (n-16<=0)||(16>=n)||(n>16);
  758. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:16,Hpd>(H,f,k,bc);
  759. }
  760. default:
  761. return min_quad_with_fixed<Scalar,Eigen::Dynamic,Eigen::Dynamic,Hpd>(H,f,k,bc);
  762. }
  763. }
  764. template <typename Scalar, int n, int kcount, bool Hpd>
  765. IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
  766. const Eigen::Matrix<Scalar,n,n> & H,
  767. const Eigen::Matrix<Scalar,n,1> & f,
  768. const Eigen::Array<bool,n,1> & k,
  769. const Eigen::Matrix<Scalar,n,1> & bc)
  770. {
  771. // 0 and n should be handle outside this function
  772. static_assert(kcount==Eigen::Dynamic || kcount>0 ,"");
  773. static_assert(kcount==Eigen::Dynamic || kcount<n ,"");
  774. const int ucount = n==Eigen::Dynamic ? Eigen::Dynamic : n-kcount;
  775. static_assert(kcount==Eigen::Dynamic || ucount+kcount == n ,"");
  776. static_assert((n==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
  777. static_assert((kcount==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
  778. assert((n==Eigen::Dynamic) || n == H.rows());
  779. assert((kcount==Eigen::Dynamic) || kcount == k.count());
  780. typedef Eigen::Matrix<Scalar,ucount,ucount> MatrixSuu;
  781. typedef Eigen::Matrix<Scalar,ucount,kcount> MatrixSuk;
  782. typedef Eigen::Matrix<Scalar,n,1> VectorSn;
  783. typedef Eigen::Matrix<Scalar,ucount,1> VectorSu;
  784. typedef Eigen::Matrix<Scalar,kcount,1> VectorSk;
  785. const auto dyn_n = n==Eigen::Dynamic ? H.rows() : n;
  786. const auto dyn_kcount = kcount==Eigen::Dynamic ? k.count() : kcount;
  787. const auto dyn_ucount = ucount==Eigen::Dynamic ? dyn_n- dyn_kcount : ucount;
  788. // For ucount==2 or kcount==2 this calls the coefficient initiliazer rather
  789. // than the size initilizer, but I guess that's ok.
  790. MatrixSuu Huu(dyn_ucount,dyn_ucount);
  791. MatrixSuk Huk(dyn_ucount,dyn_kcount);
  792. VectorSu mrhs(dyn_ucount);
  793. VectorSk bck(dyn_kcount);
  794. {
  795. int ui = 0;
  796. int ki = 0;
  797. for(int i = 0;i<dyn_n;i++)
  798. {
  799. if(k(i))
  800. {
  801. bck(ki) = bc(i);
  802. ki++;
  803. }else
  804. {
  805. mrhs(ui) = f(i);
  806. int uj = 0;
  807. int kj = 0;
  808. for(int j = 0;j<dyn_n;j++)
  809. {
  810. if(k(j))
  811. {
  812. Huk(ui,kj) = H(i,j);
  813. kj++;
  814. }else
  815. {
  816. Huu(ui,uj) = H(i,j);
  817. uj++;
  818. }
  819. }
  820. ui++;
  821. }
  822. }
  823. }
  824. mrhs += Huk * bck;
  825. typedef typename
  826. std::conditional<Hpd,
  827. Eigen::LLT<MatrixSuu>,
  828. // LDLT should be faster for indefinite problems but already found some
  829. // cases where it was too inaccurate when called via quadprog_primal.
  830. // Ideally this function takes LLT,LDLT, or
  831. // CompleteOrthogonalDecomposition as a template parameter. "template
  832. // template" parameters did work because LLT,LDLT have different number of
  833. // template parameters from CompleteOrthogonalDecomposition. Perhaps
  834. // there's a way to take advantage of LLT and LDLT's default template
  835. // parameters (I couldn't figure out how).
  836. Eigen::CompleteOrthogonalDecomposition<MatrixSuu>>::type
  837. Solver;
  838. VectorSu xu = Solver(Huu).solve(-mrhs);
  839. VectorSn x(dyn_n);
  840. {
  841. int ui = 0;
  842. int ki = 0;
  843. for(int i = 0;i<dyn_n;i++)
  844. {
  845. if(k(i))
  846. {
  847. x(i) = bck(ki);
  848. ki++;
  849. }else
  850. {
  851. x(i) = xu(ui);
  852. ui++;
  853. }
  854. }
  855. }
  856. return x;
  857. }