min_quad_with_fixed.cpp 16 KB

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