min_quad_with_fixed.impl.h 28 KB

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