min_quad_with_fixed.cpp 14 KB

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