min_quad_with_fixed.impl.h 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887
  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,1> VectorXT;
  386. typedef Matrix<T,Dynamic,Dynamic> MatrixXT;
  387. // number of known rows
  388. int kr = data.known.size();
  389. if(kr!=0)
  390. {
  391. assert(kr == Y.rows());
  392. }
  393. // number of columns to solve
  394. int cols = Y.cols();
  395. assert(B.cols() == 1 || B.cols() == cols);
  396. assert(Beq.size() == 0 || Beq.cols() == 1 || Beq.cols() == cols);
  397. // resize output
  398. Z.resize(data.n,cols);
  399. // Set known values
  400. for(int i = 0;i < kr;i++)
  401. {
  402. for(int j = 0;j < cols;j++)
  403. {
  404. Z(data.known(i),j) = Y(i,j);
  405. }
  406. }
  407. if(data.Aeq_li)
  408. {
  409. // number of lagrange multipliers aka linear equality constraints
  410. int neq = data.lagrange.size();
  411. // append lagrange multiplier rhs's
  412. MatrixXT BBeq(B.rows() + Beq.rows(),cols);
  413. if(B.size() > 0)
  414. {
  415. BBeq.topLeftCorner(B.rows(),cols) = B.replicate(1,B.cols()==cols?1:cols);
  416. }
  417. if(Beq.size() > 0)
  418. {
  419. BBeq.bottomLeftCorner(Beq.rows(),cols) = -2.0*Beq.replicate(1,Beq.cols()==cols?1:cols);
  420. }
  421. // Build right hand side
  422. MatrixXT BBequlcols;
  423. igl::slice(BBeq,data.unknown_lagrange,1,BBequlcols);
  424. MatrixXT NB;
  425. if(kr == 0)
  426. {
  427. NB = BBequlcols;
  428. }else
  429. {
  430. NB = data.preY * Y + BBequlcols;
  431. }
  432. //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
  433. //cout<<matlab_format(NB,"NB")<<endl;
  434. switch(data.solver_type)
  435. {
  436. case igl::min_quad_with_fixed_data<T>::LLT:
  437. sol = data.llt.solve(NB);
  438. break;
  439. case igl::min_quad_with_fixed_data<T>::LDLT:
  440. sol = data.ldlt.solve(NB);
  441. break;
  442. case igl::min_quad_with_fixed_data<T>::LU:
  443. // Not a bottleneck
  444. sol = data.lu.solve(NB);
  445. break;
  446. default:
  447. cerr<<"Error: invalid solver type"<<endl;
  448. return false;
  449. }
  450. //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
  451. // Now sol contains sol/-0.5
  452. sol *= -0.5;
  453. // Now sol contains solution
  454. // Place solution in Z
  455. for(int i = 0;i<(sol.rows()-neq);i++)
  456. {
  457. for(int j = 0;j<sol.cols();j++)
  458. {
  459. Z(data.unknown_lagrange(i),j) = sol(i,j);
  460. }
  461. }
  462. }else
  463. {
  464. assert(data.solver_type == min_quad_with_fixed_data<T>::QR_LLT);
  465. MatrixXT eff_Beq;
  466. // Adjust Aeq rhs to include known parts
  467. eff_Beq =
  468. //data.AeqTQR.colsPermutation().transpose() * (-data.Aeqk * Y + Beq);
  469. data.AeqTET * (-data.Aeqk * Y + Beq.replicate(1,Beq.cols()==cols?1:cols));
  470. // Where did this -0.5 come from? Probably the same place as above.
  471. MatrixXT Bu;
  472. slice(B,data.unknown,1,Bu);
  473. MatrixXT NB;
  474. NB = -0.5*(Bu.replicate(1,B.cols()==cols?1:cols) + data.preY * Y);
  475. // Trim eff_Beq
  476. const int nc = data.AeqTQR.rank();
  477. const int neq = Beq.rows();
  478. eff_Beq = eff_Beq.topLeftCorner(nc,cols).eval();
  479. data.AeqTR1T.template triangularView<Lower>().solveInPlace(eff_Beq);
  480. // Now eff_Beq = (data.AeqTR1T \ (data.AeqTET * (-data.Aeqk * Y + Beq)))
  481. MatrixXT lambda_0;
  482. lambda_0 = data.AeqTQ1 * eff_Beq;
  483. //cout<<matlab_format(lambda_0,"lambda_0")<<endl;
  484. MatrixXT QRB;
  485. QRB = -data.AeqTQ2T * (data.Auu * lambda_0) + data.AeqTQ2T * NB;
  486. Derivedsol lambda;
  487. lambda = data.llt.solve(QRB);
  488. // prepare output
  489. Derivedsol solu;
  490. solu = data.AeqTQ2 * lambda + lambda_0;
  491. // http://www.math.uh.edu/~rohop/fall_06/Chapter3.pdf
  492. Derivedsol solLambda;
  493. {
  494. Derivedsol temp1,temp2;
  495. temp1 = (data.AeqTQ1T * NB - data.AeqTQ1T * data.Auu * solu);
  496. data.AeqTR1.template triangularView<Upper>().solveInPlace(temp1);
  497. //cout<<matlab_format(temp1,"temp1")<<endl;
  498. temp2 = Derivedsol::Zero(neq,cols);
  499. temp2.topLeftCorner(nc,cols) = temp1;
  500. //solLambda = data.AeqTQR.colsPermutation() * temp2;
  501. solLambda = data.AeqTE * temp2;
  502. }
  503. // sol is [Z(unknown);Lambda]
  504. assert(data.unknown.size() == solu.rows());
  505. assert(cols == solu.cols());
  506. assert(data.neq == neq);
  507. assert(data.neq == solLambda.rows());
  508. assert(cols == solLambda.cols());
  509. sol.resize(data.unknown.size()+data.neq,cols);
  510. sol.block(0,0,solu.rows(),solu.cols()) = solu;
  511. sol.block(solu.rows(),0,solLambda.rows(),solLambda.cols()) = solLambda;
  512. for(int u = 0;u<data.unknown.size();u++)
  513. {
  514. for(int j = 0;j<Z.cols();j++)
  515. {
  516. Z(data.unknown(u),j) = solu(u,j);
  517. }
  518. }
  519. }
  520. return true;
  521. }
  522. template <
  523. typename T,
  524. typename DerivedB,
  525. typename DerivedY,
  526. typename DerivedBeq,
  527. typename DerivedZ>
  528. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  529. const min_quad_with_fixed_data<T> & data,
  530. const Eigen::MatrixBase<DerivedB> & B,
  531. const Eigen::MatrixBase<DerivedY> & Y,
  532. const Eigen::MatrixBase<DerivedBeq> & Beq,
  533. Eigen::PlainObjectBase<DerivedZ> & Z)
  534. {
  535. Eigen::Matrix<typename DerivedZ::Scalar, Eigen::Dynamic, Eigen::Dynamic> sol;
  536. return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
  537. }
  538. template <
  539. typename T,
  540. typename Derivedknown,
  541. typename DerivedB,
  542. typename DerivedY,
  543. typename DerivedBeq,
  544. typename DerivedZ>
  545. IGL_INLINE bool igl::min_quad_with_fixed(
  546. const Eigen::SparseMatrix<T>& A,
  547. const Eigen::MatrixBase<DerivedB> & B,
  548. const Eigen::MatrixBase<Derivedknown> & known,
  549. const Eigen::MatrixBase<DerivedY> & Y,
  550. const Eigen::SparseMatrix<T>& Aeq,
  551. const Eigen::MatrixBase<DerivedBeq> & Beq,
  552. const bool pd,
  553. Eigen::PlainObjectBase<DerivedZ> & Z)
  554. {
  555. min_quad_with_fixed_data<T> data;
  556. if(!min_quad_with_fixed_precompute(A,known,Aeq,pd,data))
  557. {
  558. return false;
  559. }
  560. return min_quad_with_fixed_solve(data,B,Y,Beq,Z);
  561. }
  562. template <typename Scalar, int n, int m, bool Hpd>
  563. IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
  564. const Eigen::Matrix<Scalar,n,n> & H,
  565. const Eigen::Matrix<Scalar,n,1> & f,
  566. const Eigen::Array<bool,n,1> & k,
  567. const Eigen::Matrix<Scalar,n,1> & bc,
  568. const Eigen::Matrix<Scalar,m,n> & A,
  569. const Eigen::Matrix<Scalar,m,1> & b)
  570. {
  571. const auto dyn_n = n == Eigen::Dynamic ? H.rows() : n;
  572. const auto dyn_m = m == Eigen::Dynamic ? A.rows() : m;
  573. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  574. const auto dyn_nn = nn == Eigen::Dynamic ? dyn_n+dyn_m : nn;
  575. if(dyn_m == 0)
  576. {
  577. return igl::min_quad_with_fixed<Scalar,n,Hpd>(H,f,k,bc);
  578. }
  579. // min_x ½ xᵀ H x + xᵀ f subject to A x = b and x(k) = bc(k)
  580. // let zᵀ = [xᵀ λᵀ]
  581. // min_z ½ zᵀ [H Aᵀ;A 0] z + zᵀ [f;-b] z(k) = bc(k)
  582. const auto make_HH = [&]()
  583. {
  584. // Windows can't remember that nn is const.
  585. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  586. Eigen::Matrix<Scalar,nn,nn> HH =
  587. Eigen::Matrix<Scalar,nn,nn>::Zero(dyn_nn,dyn_nn);
  588. HH.topLeftCorner(dyn_n,dyn_n) = H;
  589. HH.bottomLeftCorner(dyn_m,dyn_n) = A;
  590. HH.topRightCorner(dyn_n,dyn_m) = A.transpose();
  591. return HH;
  592. };
  593. const Eigen::Matrix<Scalar,nn,nn> HH = make_HH();
  594. const auto make_ff = [&]()
  595. {
  596. // Windows can't remember that nn is const.
  597. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  598. Eigen::Matrix<Scalar,nn,1> ff(dyn_nn);
  599. ff.head(dyn_n) = f;
  600. ff.tail(dyn_m) = -b;
  601. return ff;
  602. };
  603. const Eigen::Matrix<Scalar,nn,1> ff = make_ff();
  604. const auto make_kk = [&]()
  605. {
  606. // Windows can't remember that nn is const.
  607. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  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. // Windows can't remember that nn is const.
  617. constexpr const int nn = n == Eigen::Dynamic ? Eigen::Dynamic : n+m;
  618. Eigen::Matrix<Scalar,nn,1> bcbc(dyn_nn);
  619. bcbc.head(dyn_n) = bc;
  620. return bcbc;
  621. };
  622. const Eigen::Matrix<Scalar,nn,1> bcbc = make_bcbc();
  623. const Eigen::Matrix<Scalar,nn,1> xx =
  624. min_quad_with_fixed<Scalar,nn,false>(HH,ff,kk,bcbc);
  625. return xx.head(dyn_n);
  626. }
  627. template <typename Scalar, int n, bool Hpd>
  628. IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
  629. const Eigen::Matrix<Scalar,n,n> & H,
  630. const Eigen::Matrix<Scalar,n,1> & f,
  631. const Eigen::Array<bool,n,1> & k,
  632. const Eigen::Matrix<Scalar,n,1> & bc)
  633. {
  634. assert(H.isApprox(H.transpose(),1e-7));
  635. assert(H.rows() == H.cols());
  636. assert(H.rows() == f.size());
  637. assert(H.rows() == k.size());
  638. assert(H.rows() == bc.size());
  639. const auto kcount = k.count();
  640. // Everything fixed
  641. if(kcount == (Eigen::Dynamic?H.rows():n))
  642. {
  643. return bc;
  644. }
  645. // Nothing fixed
  646. if(kcount == 0)
  647. {
  648. // avoid function call
  649. typedef Eigen::Matrix<Scalar,n,n> MatrixSn;
  650. typedef typename
  651. std::conditional<Hpd,Eigen::LLT<MatrixSn>,Eigen::CompleteOrthogonalDecomposition<MatrixSn>>::type
  652. Solver;
  653. return Solver(H).solve(-f);
  654. }
  655. // All-but-one fixed
  656. if( (Eigen::Dynamic?H.rows():n)-kcount == 1)
  657. {
  658. // which one is not fixed?
  659. int u = -1;
  660. for(int i=0;i<k.size();i++){ if(!k(i)){ u=i; break; } }
  661. assert(u>=0);
  662. // min ½ x(u) Huu x(u) + x(u)(fu + H(u,k)bc(k))
  663. // Huu x(u) = -(fu + H(u,k) bc(k))
  664. // x(u) = (-fu + ∑ -Huj bcj)/Huu
  665. Eigen::Matrix<Scalar,n,1> x = bc;
  666. x(u) = -f(u);
  667. for(int i=0;i<k.size();i++){ if(i!=u){ x(u)-=bc(i)*H(i,u); } }
  668. x(u) /= H(u,u);
  669. return x;
  670. }
  671. // Alec: Is there a smart template way to do this?
  672. // jdumas: I guess you could do a templated for-loop starting from 16, and
  673. // dispatching to the appropriate templated function when the argument matches
  674. // (with a fallback to the dynamic version). Cf this example:
  675. // https://gist.github.com/disconnect3d/13c2d035bb31b244df14
  676. switch(kcount)
  677. {
  678. case 0: assert(false && "Handled above."); return Eigen::Matrix<Scalar,n,1>();
  679. // % Matlibberish for generating these case statements:
  680. // 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
  681. case 1:
  682. {
  683. const bool D = (n-1<=0)||(1>=n)||(n>16);
  684. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:1,Hpd>(H,f,k,bc);
  685. }
  686. case 2:
  687. {
  688. const bool D = (n-2<=0)||(2>=n)||(n>16);
  689. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:2,Hpd>(H,f,k,bc);
  690. }
  691. case 3:
  692. {
  693. const bool D = (n-3<=0)||(3>=n)||(n>16);
  694. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:3,Hpd>(H,f,k,bc);
  695. }
  696. case 4:
  697. {
  698. const bool D = (n-4<=0)||(4>=n)||(n>16);
  699. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:4,Hpd>(H,f,k,bc);
  700. }
  701. case 5:
  702. {
  703. const bool D = (n-5<=0)||(5>=n)||(n>16);
  704. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:5,Hpd>(H,f,k,bc);
  705. }
  706. case 6:
  707. {
  708. const bool D = (n-6<=0)||(6>=n)||(n>16);
  709. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:6,Hpd>(H,f,k,bc);
  710. }
  711. case 7:
  712. {
  713. const bool D = (n-7<=0)||(7>=n)||(n>16);
  714. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:7,Hpd>(H,f,k,bc);
  715. }
  716. case 8:
  717. {
  718. const bool D = (n-8<=0)||(8>=n)||(n>16);
  719. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:8,Hpd>(H,f,k,bc);
  720. }
  721. case 9:
  722. {
  723. const bool D = (n-9<=0)||(9>=n)||(n>16);
  724. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:9,Hpd>(H,f,k,bc);
  725. }
  726. case 10:
  727. {
  728. const bool D = (n-10<=0)||(10>=n)||(n>16);
  729. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:10,Hpd>(H,f,k,bc);
  730. }
  731. case 11:
  732. {
  733. const bool D = (n-11<=0)||(11>=n)||(n>16);
  734. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:11,Hpd>(H,f,k,bc);
  735. }
  736. case 12:
  737. {
  738. const bool D = (n-12<=0)||(12>=n)||(n>16);
  739. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:12,Hpd>(H,f,k,bc);
  740. }
  741. case 13:
  742. {
  743. const bool D = (n-13<=0)||(13>=n)||(n>16);
  744. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:13,Hpd>(H,f,k,bc);
  745. }
  746. case 14:
  747. {
  748. const bool D = (n-14<=0)||(14>=n)||(n>16);
  749. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:14,Hpd>(H,f,k,bc);
  750. }
  751. case 15:
  752. {
  753. const bool D = (n-15<=0)||(15>=n)||(n>16);
  754. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:15,Hpd>(H,f,k,bc);
  755. }
  756. case 16:
  757. {
  758. const bool D = (n-16<=0)||(16>=n)||(n>16);
  759. return min_quad_with_fixed<Scalar,D?Eigen::Dynamic:n,D?Eigen::Dynamic:16,Hpd>(H,f,k,bc);
  760. }
  761. default:
  762. return min_quad_with_fixed<Scalar,Eigen::Dynamic,Eigen::Dynamic,Hpd>(H,f,k,bc);
  763. }
  764. }
  765. template <typename Scalar, int n, int kcount, bool Hpd>
  766. IGL_INLINE Eigen::Matrix<Scalar,n,1> igl::min_quad_with_fixed(
  767. const Eigen::Matrix<Scalar,n,n> & H,
  768. const Eigen::Matrix<Scalar,n,1> & f,
  769. const Eigen::Array<bool,n,1> & k,
  770. const Eigen::Matrix<Scalar,n,1> & bc)
  771. {
  772. // 0 and n should be handle outside this function
  773. static_assert(kcount==Eigen::Dynamic || kcount>0 ,"");
  774. static_assert(kcount==Eigen::Dynamic || kcount<n ,"");
  775. const int ucount = n==Eigen::Dynamic ? Eigen::Dynamic : n-kcount;
  776. static_assert(kcount==Eigen::Dynamic || ucount+kcount == n ,"");
  777. static_assert((n==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
  778. static_assert((kcount==Eigen::Dynamic) == (ucount==Eigen::Dynamic),"");
  779. assert((n==Eigen::Dynamic) || n == H.rows());
  780. assert((kcount==Eigen::Dynamic) || kcount == k.count());
  781. typedef Eigen::Matrix<Scalar,ucount,ucount> MatrixSuu;
  782. typedef Eigen::Matrix<Scalar,ucount,kcount> MatrixSuk;
  783. typedef Eigen::Matrix<Scalar,n,1> VectorSn;
  784. typedef Eigen::Matrix<Scalar,ucount,1> VectorSu;
  785. typedef Eigen::Matrix<Scalar,kcount,1> VectorSk;
  786. const auto dyn_n = n==Eigen::Dynamic ? H.rows() : n;
  787. const auto dyn_kcount = kcount==Eigen::Dynamic ? k.count() : kcount;
  788. const auto dyn_ucount = ucount==Eigen::Dynamic ? dyn_n- dyn_kcount : ucount;
  789. // For ucount==2 or kcount==2 this calls the coefficient initiliazer rather
  790. // than the size initilizer, but I guess that's ok.
  791. MatrixSuu Huu(dyn_ucount,dyn_ucount);
  792. MatrixSuk Huk(dyn_ucount,dyn_kcount);
  793. VectorSu mrhs(dyn_ucount);
  794. VectorSk bck(dyn_kcount);
  795. {
  796. int ui = 0;
  797. int ki = 0;
  798. for(int i = 0;i<dyn_n;i++)
  799. {
  800. if(k(i))
  801. {
  802. bck(ki) = bc(i);
  803. ki++;
  804. }else
  805. {
  806. mrhs(ui) = f(i);
  807. int uj = 0;
  808. int kj = 0;
  809. for(int j = 0;j<dyn_n;j++)
  810. {
  811. if(k(j))
  812. {
  813. Huk(ui,kj) = H(i,j);
  814. kj++;
  815. }else
  816. {
  817. Huu(ui,uj) = H(i,j);
  818. uj++;
  819. }
  820. }
  821. ui++;
  822. }
  823. }
  824. }
  825. mrhs += Huk * bck;
  826. typedef typename
  827. std::conditional<Hpd,
  828. Eigen::LLT<MatrixSuu>,
  829. // LDLT should be faster for indefinite problems but already found some
  830. // cases where it was too inaccurate when called via quadprog_primal.
  831. // Ideally this function takes LLT,LDLT, or
  832. // CompleteOrthogonalDecomposition as a template parameter. "template
  833. // template" parameters did work because LLT,LDLT have different number of
  834. // template parameters from CompleteOrthogonalDecomposition. Perhaps
  835. // there's a way to take advantage of LLT and LDLT's default template
  836. // parameters (I couldn't figure out how).
  837. Eigen::CompleteOrthogonalDecomposition<MatrixSuu>>::type
  838. Solver;
  839. VectorSu xu = Solver(Huu).solve(-mrhs);
  840. VectorSn x(dyn_n);
  841. {
  842. int ui = 0;
  843. int ki = 0;
  844. for(int i = 0;i<dyn_n;i++)
  845. {
  846. if(k(i))
  847. {
  848. x(i) = bck(ki);
  849. ki++;
  850. }else
  851. {
  852. x(i) = xu(ui);
  853. ui++;
  854. }
  855. }
  856. }
  857. return x;
  858. }