min_quad_with_fixed.cpp 15 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515
  1. #include "min_quad_with_fixed.h"
  2. #include "slice.h"
  3. #include "is_symmetric.h"
  4. #include "find.h"
  5. #include "sparse.h"
  6. #include "repmat.h"
  7. #include "lu_lagrange.h"
  8. #include "full.h"
  9. #include "matlab_format.h"
  10. #include "EPS.h"
  11. #include "cat.h"
  12. //#include <Eigen/SparseExtra>
  13. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  14. #include <iostream>
  15. #include <unsupported/Eigen/SparseExtra>
  16. #include <cassert>
  17. #include <cstdio>
  18. #include <iostream>
  19. template <typename T, typename Derivedknown>
  20. IGL_INLINE bool igl::min_quad_with_fixed_precompute(
  21. const Eigen::SparseMatrix<T>& A,
  22. const Eigen::PlainObjectBase<Derivedknown> & known,
  23. const Eigen::SparseMatrix<T>& Aeq,
  24. const bool pd,
  25. min_quad_with_fixed_data<T> & data
  26. )
  27. {
  28. #define MIN_QUAD_WITH_FIXED_CPP_DEBUG
  29. using namespace Eigen;
  30. using namespace std;
  31. using namespace igl;
  32. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  33. cout<<" pre"<<endl;
  34. #endif
  35. // number of rows
  36. int n = A.rows();
  37. // cache problem size
  38. data.n = n;
  39. int neq = Aeq.rows();
  40. // default is to have 0 linear equality constraints
  41. if(Aeq.size() != 0)
  42. {
  43. assert(n == Aeq.cols());
  44. }
  45. assert(A.rows() == n);
  46. assert(A.cols() == n);
  47. // number of known rows
  48. int kr = known.size();
  49. assert(kr == 0 || known.minCoeff() >= 0);
  50. assert(kr == 0 || known.maxCoeff() < n);
  51. assert(neq <= n);
  52. // cache known
  53. data.known = known;
  54. // get list of unknown indices
  55. data.unknown.resize(n-kr);
  56. std::vector<bool> unknown_mask;
  57. unknown_mask.resize(n,true);
  58. for(int i = 0;i<kr;i++)
  59. {
  60. unknown_mask[known(i)] = false;
  61. }
  62. int u = 0;
  63. for(int i = 0;i<n;i++)
  64. {
  65. if(unknown_mask[i])
  66. {
  67. data.unknown(u) = i;
  68. u++;
  69. }
  70. }
  71. // get list of lagrange multiplier indices
  72. data.lagrange.resize(neq);
  73. for(int i = 0;i<neq;i++)
  74. {
  75. data.lagrange(i) = n + i;
  76. }
  77. // cache unknown followed by lagrange indices
  78. data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
  79. data.unknown_lagrange << data.unknown, data.lagrange;
  80. SparseMatrix<T> Auu;
  81. slice(A,data.unknown,data.unknown,Auu);
  82. // Positive definiteness is *not* determined, rather it is given as a
  83. // parameter
  84. data.Auu_pd = pd;
  85. if(data.Auu_pd)
  86. {
  87. // PD implies symmetric
  88. data.Auu_sym = true;
  89. assert(is_symmetric(Auu,EPS<double>()));
  90. }else
  91. {
  92. // determine if A(unknown,unknown) is symmetric and/or positive definite
  93. data.Auu_sym = is_symmetric(Auu,EPS<double>());
  94. }
  95. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  96. cout<<" qr"<<endl;
  97. #endif
  98. // Determine number of linearly independent constraints
  99. int nc = 0;
  100. if(neq>0)
  101. {
  102. // QR decomposition to determine row rank in Aequ
  103. slice(Aeq,data.unknown,2,data.Aequ);
  104. assert(data.Aequ.rows() == neq);
  105. assert(data.Aequ.cols() == data.unknown.size());
  106. data.AeqTQR.compute(data.Aequ.transpose().eval());
  107. switch(data.AeqTQR.info())
  108. {
  109. case Eigen::Success:
  110. break;
  111. case Eigen::NumericalIssue:
  112. cerr<<"Error: Numerical issue."<<endl;
  113. return false;
  114. case Eigen::InvalidInput:
  115. cerr<<"Error: Invalid input."<<endl;
  116. return false;
  117. default:
  118. cerr<<"Error: Other."<<endl;
  119. return false;
  120. }
  121. nc = data.AeqTQR.rank();
  122. assert(nc<=neq);
  123. data.Aeq_li = nc == neq;
  124. }else
  125. {
  126. data.Aeq_li = true;
  127. }
  128. if(data.Aeq_li)
  129. {
  130. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  131. cout<<" Aeq_li=true"<<endl;
  132. #endif
  133. // Append lagrange multiplier quadratic terms
  134. SparseMatrix<T> new_A;
  135. SparseMatrix<T> AeqT = Aeq.transpose();
  136. SparseMatrix<T> Z(neq,neq);
  137. // This is a bit slower. But why isn't cat fast?
  138. new_A = cat(1, cat(2, A, AeqT ),
  139. cat(2, Aeq, Z ));
  140. // precompute RHS builders
  141. if(kr > 0)
  142. {
  143. SparseMatrix<T> Aulk,Akul;
  144. // Slow
  145. slice(new_A,data.unknown_lagrange,data.known,Aulk);
  146. //// This doesn't work!!!
  147. //data.preY = Aulk + Akul.transpose();
  148. // Slow
  149. if(data.Auu_sym)
  150. {
  151. data.preY = Aulk*2;
  152. }else
  153. {
  154. slice(new_A,data.known,data.unknown_lagrange,Akul);
  155. SparseMatrix<T> AkulT = Akul.transpose();
  156. data.preY = Aulk + AkulT;
  157. }
  158. }else
  159. {
  160. data.preY.resize(data.unknown_lagrange.size(),0);
  161. }
  162. // Positive definite and no equality constraints (Postive definiteness
  163. // implies symmetric)
  164. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  165. cout<<" factorize"<<endl;
  166. #endif
  167. if(data.Auu_pd && neq == 0)
  168. {
  169. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  170. cout<<" llt"<<endl;
  171. #endif
  172. data.llt.compute(Auu);
  173. switch(data.llt.info())
  174. {
  175. case Eigen::Success:
  176. break;
  177. case Eigen::NumericalIssue:
  178. cerr<<"Error: Numerical issue."<<endl;
  179. return false;
  180. default:
  181. cerr<<"Error: Other."<<endl;
  182. return false;
  183. }
  184. data.solver_type = min_quad_with_fixed_data<T>::LLT;
  185. }else
  186. {
  187. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  188. cout<<" ldlt"<<endl;
  189. #endif
  190. // Either not PD or there are equality constraints
  191. SparseMatrix<T> NA;
  192. slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  193. data.NA = NA;
  194. // Ideally we'd use LDLT but Eigen doesn't support positive semi-definite
  195. // matrices:
  196. // http://forum.kde.org/viewtopic.php?f=74&t=106962&p=291990#p291990
  197. if(data.Auu_sym && false)
  198. {
  199. data.ldlt.compute(NA);
  200. switch(data.ldlt.info())
  201. {
  202. case Eigen::Success:
  203. break;
  204. case Eigen::NumericalIssue:
  205. cerr<<"Error: Numerical issue."<<endl;
  206. return false;
  207. default:
  208. cerr<<"Error: Other."<<endl;
  209. return false;
  210. }
  211. data.solver_type = min_quad_with_fixed_data<T>::LDLT;
  212. }else
  213. {
  214. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  215. cout<<" lu"<<endl;
  216. #endif
  217. // Resort to LU
  218. // Bottleneck >1/2
  219. data.lu.compute(NA);
  220. //std::cout<<"NA=["<<std::endl<<NA<<std::endl<<"];"<<std::endl;
  221. switch(data.lu.info())
  222. {
  223. case Eigen::Success:
  224. break;
  225. case Eigen::NumericalIssue:
  226. cerr<<"Error: Numerical issue."<<endl;
  227. return false;
  228. case Eigen::InvalidInput:
  229. cerr<<"Error: Invalid Input."<<endl;
  230. return false;
  231. default:
  232. cerr<<"Error: Other."<<endl;
  233. return false;
  234. }
  235. data.solver_type = min_quad_with_fixed_data<T>::LU;
  236. }
  237. }
  238. }else
  239. {
  240. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  241. cout<<" Aeq_li=false"<<endl;
  242. #endif
  243. data.neq = neq;
  244. const int nu = data.unknown.size();
  245. //cout<<"nu: "<<nu<<endl;
  246. //cout<<"neq: "<<neq<<endl;
  247. //cout<<"nc: "<<nc<<endl;
  248. //cout<<" matrixR"<<endl;
  249. SparseMatrix<T> AeqTR,AeqTQ;
  250. AeqTR = data.AeqTQR.matrixR();
  251. // This shouldn't be necessary
  252. AeqTR.prune(0.0);
  253. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  254. cout<<" matrixQ"<<endl;
  255. #endif
  256. // THIS IS ESSENTIALLY DENSE AND THIS IS BY FAR THE BOTTLENECK
  257. // http://forum.kde.org/viewtopic.php?f=74&t=117500
  258. AeqTQ = data.AeqTQR.matrixQ();
  259. // This shouldn't be necessary
  260. AeqTQ.prune(0.0);
  261. //cout<<"AeqTQ: "<<AeqTQ.rows()<<" "<<AeqTQ.cols()<<endl;
  262. //cout<<matlab_format(AeqTQ,"AeqTQ")<<endl;
  263. //cout<<" perms"<<endl;
  264. SparseMatrix<double> I(neq,neq);
  265. I.setIdentity();
  266. data.AeqTE = data.AeqTQR.colsPermutation() * I;
  267. data.AeqTET = data.AeqTQR.colsPermutation().transpose() * I;
  268. assert(AeqTR.rows() == neq);
  269. assert(AeqTQ.rows() == nu);
  270. assert(AeqTQ.cols() == nu);
  271. assert(AeqTR.cols() == neq);
  272. //cout<<" slice"<<endl;
  273. data.AeqTQ1 = AeqTQ.topLeftCorner(nu,nc);
  274. data.AeqTQ1T = data.AeqTQ1.transpose().eval();
  275. // ALREADY TRIM (Not 100% sure about this)
  276. data.AeqTR1 = AeqTR.topLeftCorner(nc,nc);
  277. data.AeqTR1T = data.AeqTR1.transpose().eval();
  278. //cout<<"AeqTR1T.size() "<<data.AeqTR1T.rows()<<" "<<data.AeqTR1T.cols()<<endl;
  279. // Null space
  280. data.AeqTQ2 = AeqTQ.bottomRightCorner(nu,nu-nc);
  281. data.AeqTQ2T = data.AeqTQ2.transpose().eval();
  282. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  283. cout<<" proj"<<endl;
  284. #endif
  285. // Projected hessian
  286. SparseMatrix<T> QRAuu = data.AeqTQ2T * Auu * data.AeqTQ2;
  287. {
  288. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  289. cout<<" factorize"<<endl;
  290. #endif
  291. // QRAuu should always be PD
  292. data.llt.compute(QRAuu);
  293. switch(data.llt.info())
  294. {
  295. case Eigen::Success:
  296. break;
  297. case Eigen::NumericalIssue:
  298. cerr<<"Error: Numerical issue."<<endl;
  299. return false;
  300. default:
  301. cerr<<"Error: Other."<<endl;
  302. return false;
  303. }
  304. data.solver_type = min_quad_with_fixed_data<T>::QR_LLT;
  305. }
  306. #ifdef MIN_QUAD_WITH_FIXED_CPP_DEBUG
  307. cout<<" smash"<<endl;
  308. #endif
  309. // Known value multiplier
  310. SparseMatrix<T> Auk;
  311. slice(A,data.unknown,data.known,Auk);
  312. SparseMatrix<T> Aku;
  313. slice(A,data.known,data.unknown,Aku);
  314. SparseMatrix<T> AkuT = Aku.transpose();
  315. data.preY = Auk + AkuT;
  316. // Needed during solve
  317. data.Auu = Auu;
  318. slice(Aeq,data.known,2,data.Aeqk);
  319. assert(data.Aeqk.rows() == neq);
  320. assert(data.Aeqk.cols() == data.known.size());
  321. }
  322. return true;
  323. }
  324. template <
  325. typename T,
  326. typename DerivedB,
  327. typename DerivedY,
  328. typename DerivedBeq,
  329. typename DerivedZ,
  330. typename Derivedsol>
  331. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  332. const min_quad_with_fixed_data<T> & data,
  333. const Eigen::PlainObjectBase<DerivedB> & B,
  334. const Eigen::PlainObjectBase<DerivedY> & Y,
  335. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  336. Eigen::PlainObjectBase<DerivedZ> & Z,
  337. Eigen::PlainObjectBase<Derivedsol> & sol)
  338. {
  339. using namespace std;
  340. using namespace Eigen;
  341. using namespace igl;
  342. // number of known rows
  343. int kr = data.known.size();
  344. if(kr!=0)
  345. {
  346. assert(kr == Y.rows());
  347. }
  348. // number of columns to solve
  349. int cols = Y.cols();
  350. assert(B.cols() == 1);
  351. assert(Beq.size() == 0 || Beq.cols() == 1);
  352. // resize output
  353. Z.resize(data.n,cols);
  354. // Set known values
  355. for(int i = 0;i < kr;i++)
  356. {
  357. for(int j = 0;j < cols;j++)
  358. {
  359. Z(data.known(i),j) = Y(i,j);
  360. }
  361. }
  362. if(data.Aeq_li)
  363. {
  364. // number of lagrange multipliers aka linear equality constraints
  365. int neq = data.lagrange.size();
  366. // append lagrange multiplier rhs's
  367. Eigen::Matrix<T,Eigen::Dynamic,1> BBeq(B.size() + Beq.size());
  368. BBeq << B, (Beq*-2.0);
  369. // Build right hand side
  370. Eigen::Matrix<T,Eigen::Dynamic,1> BBequl;
  371. igl::slice(BBeq,data.unknown_lagrange,BBequl);
  372. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> BBequlcols;
  373. igl::repmat(BBequl,1,cols,BBequlcols);
  374. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  375. if(kr == 0)
  376. {
  377. NB = BBequlcols;
  378. }else
  379. {
  380. NB = data.preY * Y + BBequlcols;
  381. }
  382. //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
  383. //cout<<matlab_format(NB,"NB")<<endl;
  384. switch(data.solver_type)
  385. {
  386. case igl::min_quad_with_fixed_data<T>::LLT:
  387. sol = data.llt.solve(NB);
  388. break;
  389. case igl::min_quad_with_fixed_data<T>::LDLT:
  390. sol = data.ldlt.solve(NB);
  391. break;
  392. case igl::min_quad_with_fixed_data<T>::LU:
  393. // Not a bottleneck
  394. sol = data.lu.solve(NB);
  395. break;
  396. default:
  397. cerr<<"Error: invalid solver type"<<endl;
  398. return false;
  399. }
  400. //std::cout<<"sol=["<<std::endl<<sol<<std::endl<<"];"<<std::endl;
  401. // Now sol contains sol/-0.5
  402. sol *= -0.5;
  403. // Now sol contains solution
  404. // Place solution in Z
  405. for(int i = 0;i<(sol.rows()-neq);i++)
  406. {
  407. for(int j = 0;j<sol.cols();j++)
  408. {
  409. Z(data.unknown_lagrange(i),j) = sol(i,j);
  410. }
  411. }
  412. }else
  413. {
  414. assert(data.solver_type == min_quad_with_fixed_data<T>::QR_LLT);
  415. PlainObjectBase<DerivedBeq> eff_Beq;
  416. // Adjust Aeq rhs to include known parts
  417. eff_Beq =
  418. //data.AeqTQR.colsPermutation().transpose() * (-data.Aeqk * Y + Beq);
  419. data.AeqTET * (-data.Aeqk * Y + Beq);
  420. // Where did this -0.5 come from? Probably the same place as above.
  421. PlainObjectBase<DerivedB> Bu;
  422. slice(B,data.unknown,Bu);
  423. PlainObjectBase<DerivedB> NB;
  424. NB = -0.5*(Bu + data.preY * Y);
  425. // Trim eff_Beq
  426. const int nc = data.AeqTQR.rank();
  427. const int neq = Beq.rows();
  428. eff_Beq = eff_Beq.topLeftCorner(nc,1).eval();
  429. data.AeqTR1T.template triangularView<Lower>().solveInPlace(eff_Beq);
  430. // Now eff_Beq = (data.AeqTR1T \ (data.AeqTET * (-data.Aeqk * Y + Beq)))
  431. PlainObjectBase<DerivedB> lambda_0;
  432. lambda_0 = data.AeqTQ1 * eff_Beq;
  433. //cout<<matlab_format(lambda_0,"lambda_0")<<endl;
  434. PlainObjectBase<DerivedB> QRB;
  435. QRB = -data.AeqTQ2T * (data.Auu * lambda_0) + data.AeqTQ2T * NB;
  436. PlainObjectBase<Derivedsol> lambda;
  437. lambda = data.llt.solve(QRB);
  438. // prepare output
  439. PlainObjectBase<Derivedsol> solu;
  440. solu = data.AeqTQ2 * lambda + lambda_0;
  441. // http://www.math.uh.edu/~rohop/fall_06/Chapter3.pdf
  442. PlainObjectBase<Derivedsol> solLambda;
  443. {
  444. PlainObjectBase<Derivedsol> temp1,temp2;
  445. temp1 = (data.AeqTQ1T * NB - data.AeqTQ1T * data.Auu * solu);
  446. data.AeqTR1.template triangularView<Upper>().solveInPlace(temp1);
  447. //cout<<matlab_format(temp1,"temp1")<<endl;
  448. temp2 = PlainObjectBase<Derivedsol>::Zero(neq,1);
  449. temp2.topLeftCorner(nc,1) = temp1;
  450. //solLambda = data.AeqTQR.colsPermutation() * temp2;
  451. solLambda = data.AeqTE * temp2;
  452. }
  453. // sol is [Z(unknown);Lambda]
  454. assert(data.unknown.size() == solu.rows());
  455. assert(cols == solu.cols());
  456. assert(data.neq == neq);
  457. assert(data.neq == solLambda.rows());
  458. assert(cols == solLambda.cols());
  459. sol.resize(data.unknown.size()+data.neq,cols);
  460. sol.block(0,0,solu.rows(),solu.cols()) = solu;
  461. sol.block(solu.rows(),0,solLambda.rows(),solLambda.cols()) = solLambda;
  462. for(int u = 0;u<data.unknown.size();u++)
  463. {
  464. for(int j = 0;j<Z.cols();j++)
  465. {
  466. Z(data.unknown(u),j) = solu(u,j);
  467. }
  468. }
  469. }
  470. return true;
  471. }
  472. template <
  473. typename T,
  474. typename DerivedB,
  475. typename DerivedY,
  476. typename DerivedBeq,
  477. typename DerivedZ>
  478. IGL_INLINE bool igl::min_quad_with_fixed_solve(
  479. const min_quad_with_fixed_data<T> & data,
  480. const Eigen::PlainObjectBase<DerivedB> & B,
  481. const Eigen::PlainObjectBase<DerivedY> & Y,
  482. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  483. Eigen::PlainObjectBase<DerivedZ> & Z)
  484. {
  485. Eigen::PlainObjectBase<DerivedZ> sol;
  486. return min_quad_with_fixed_solve(data,B,Y,Beq,Z,sol);
  487. }
  488. #ifndef IGL_HEADER_ONLY
  489. // Explicit template specialization
  490. template bool igl::min_quad_with_fixed_solve<double, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(igl::min_quad_with_fixed_data<double> const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> > const&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  491. template bool igl::min_quad_with_fixed_precompute<double, Eigen::Matrix<int, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::PlainObjectBase<Eigen::Matrix<int, -1, 1, 0, -1, 1> > const&, Eigen::SparseMatrix<double, 0, int> const&, bool, igl::min_quad_with_fixed_data<double>&);
  492. #endif