min_quad_with_fixed.impl.h 28 KB

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