min_quad_with_fixed.h 9.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332
  1. #ifndef IGL_MIN_QUAD_WITH_FIXED_H
  2. #define IGL_MIN_QUAD_WITH_FIXED_H
  3. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  4. #include <Eigen/Core>
  5. #include <Eigen/Dense>
  6. #include <Eigen/Sparse>
  7. #include <Eigen/SparseExtra>
  8. namespace igl
  9. {
  10. template <typename T>
  11. struct min_quad_with_fixed_data;
  12. // MIN_QUAD_WITH_FIXED Minimize quadratic energy Z'*A*Z + Z'*B + C with
  13. // constraints that Z(known) = Y, optionally also subject to the constraints
  14. // Aeq*Z = Beq
  15. //
  16. // Templates:
  17. // T should be a eigen matrix primitive type like int or double
  18. // Inputs:
  19. // A n by n matrix of quadratic coefficients
  20. // B n by 1 column of linear coefficients
  21. // known list of indices to known rows in Z
  22. // Y list of fixed values corresponding to known rows in Z
  23. // Optional:
  24. // Aeq m by n list of linear equality constraint coefficients
  25. // Beq m by 1 list of linear equality constraint constant values
  26. // pd flag specifying whether A(unknown,unknown) is positive definite
  27. // Outputs:
  28. // data factorization struct with all necessary information to solve
  29. // using min_quad_with_fixed_solve
  30. // Returns true on success, false on error
  31. template <typename T>
  32. inline bool min_quad_with_fixed_precompute(
  33. const Eigen::SparseMatrix<T>& A,
  34. const Eigen::Matrix<int,Eigen::Dynamic,1> & known,
  35. const Eigen::SparseMatrix<T>& Aeq,
  36. const bool pd,
  37. min_quad_with_fixed_data<T> & data
  38. );
  39. // Solves a system previously factored using min_quad_with_fixed_precompute
  40. // Inputs:
  41. // data factorization struct with all necessary precomputation to solve
  42. // Outputs:
  43. // Z n by cols solution
  44. // Returns true on success, false on error
  45. template <typename T>
  46. inline bool min_quad_with_fixed_solve(
  47. const min_quad_with_fixed_data<T> & data,
  48. const Eigen::Matrix<T,Eigen::Dynamic,1> & B,
  49. const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> & Y,
  50. const Eigen::Matrix<T,Eigen::Dynamic,1> & Beq,
  51. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> & Z);
  52. }
  53. // Implementation
  54. #include <Eigen/SparseExtra>
  55. #include <cassert>
  56. #include <cstdio>
  57. #include <iostream>
  58. #include "slice.h"
  59. #include "is_symmetric.h"
  60. #include "find.h"
  61. #include "sparse.h"
  62. #include "repmat.h"
  63. #include "lu_lagrange.h"
  64. #include "full.h"
  65. template <typename T>
  66. struct igl::min_quad_with_fixed_data
  67. {
  68. // Size of original system: number of unknowns + number of knowns
  69. int n;
  70. // Whether A(unknown,unknown) is positive definite
  71. bool Auu_pd;
  72. // Whether A(unknown,unknown) is symmetric
  73. bool Auu_sym;
  74. // Indices of known variables
  75. Eigen::Matrix<int,Eigen::Dynamic,1> known;
  76. // Indices of unknown variables
  77. Eigen::Matrix<int,Eigen::Dynamic,1> unknown;
  78. // Indices of lagrange variables
  79. Eigen::Matrix<int,Eigen::Dynamic,1> lagrange;
  80. // Indices of unknown variable followed by Indices of lagrange variables
  81. Eigen::Matrix<int,Eigen::Dynamic,1> unknown_lagrange;
  82. // Matrix multiplied against Y when constructing right hand side
  83. Eigen::SparseMatrix<T> preY;
  84. // Tells whether system is sparse
  85. bool sparse;
  86. // Lower triangle of LU decomposition of final system matrix
  87. Eigen::SparseMatrix<T> L;
  88. // Upper triangle of LU decomposition of final system matrix
  89. Eigen::SparseMatrix<T> U;
  90. // Dense LU factorization
  91. Eigen::FullPivLU<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> > lu;
  92. };
  93. template <typename T>
  94. inline bool igl::min_quad_with_fixed_precompute(
  95. const Eigen::SparseMatrix<T>& A,
  96. const Eigen::Matrix<int,Eigen::Dynamic,1> & known,
  97. const Eigen::SparseMatrix<T>& Aeq,
  98. const bool pd,
  99. igl::min_quad_with_fixed_data<T> & data
  100. )
  101. {
  102. // number of rows
  103. int n = A.rows();
  104. // cache problem size
  105. data.n = n;
  106. int neq = Aeq.rows();
  107. // default is to have 0 linear equality constraints
  108. if(Aeq.size() != 0)
  109. {
  110. assert(n == Aeq.cols());
  111. }
  112. assert(A.rows() == n);
  113. assert(A.cols() == n);
  114. // number of known rows
  115. int kr = known.size();
  116. assert(kr == 0 || known.minCoeff() >= 0);
  117. assert(kr == 0 || known.maxCoeff() < n);
  118. assert(neq <= n);
  119. // cache known
  120. data.known = known;
  121. // get list of unknown indices
  122. data.unknown.resize(n-kr);
  123. std::vector<bool> unknown_mask;
  124. unknown_mask.resize(n,true);
  125. for(int i = 0;i<kr;i++)
  126. {
  127. unknown_mask[known(i)] = false;
  128. }
  129. int u = 0;
  130. for(int i = 0;i<n;i++)
  131. {
  132. if(unknown_mask[i])
  133. {
  134. data.unknown(u) = i;
  135. u++;
  136. }
  137. }
  138. // get list of lagrange multiplier indices
  139. data.lagrange.resize(neq);
  140. for(int i = 0;i<neq;i++)
  141. {
  142. data.lagrange(i) = n + i;
  143. }
  144. // cache unknown followed by lagrange indices
  145. data.unknown_lagrange.resize(data.unknown.size()+data.lagrange.size());
  146. data.unknown_lagrange << data.unknown, data.lagrange;
  147. Eigen::SparseMatrix<T> Auu;
  148. igl::slice(A,data.unknown,data.unknown,Auu);
  149. // determine if A(unknown,unknown) is symmetric and/or positive definite
  150. data.Auu_sym = igl::is_symmetric(Auu);
  151. // Positive definiteness is *not* determined, rather it is given as a
  152. // parameter
  153. data.Auu_pd = pd;
  154. // Append lagrange multiplier quadratic terms
  155. Eigen::SparseMatrix<T> new_A;
  156. Eigen::Matrix<int,Eigen::Dynamic,1> AI;
  157. Eigen::Matrix<int,Eigen::Dynamic,1> AJ;
  158. Eigen::Matrix<T,Eigen::Dynamic,1> AV;
  159. igl::find(A,AI,AJ,AV);
  160. Eigen::Matrix<int,Eigen::Dynamic,1> AeqI(0,1);
  161. Eigen::Matrix<int,Eigen::Dynamic,1> AeqJ(0,1);
  162. Eigen::Matrix<T,Eigen::Dynamic,1> AeqV(0,1);
  163. if(neq > 0)
  164. {
  165. igl::find(Aeq,AeqI,AeqJ,AeqV);
  166. }
  167. Eigen::Matrix<int,Eigen::Dynamic,1> new_AI(AV.size()+AeqV.size()*2);
  168. Eigen::Matrix<int,Eigen::Dynamic,1> new_AJ(AV.size()+AeqV.size()*2);
  169. Eigen::Matrix<T,Eigen::Dynamic,1> new_AV(AV.size()+AeqV.size()*2);
  170. new_AI << AI, (AeqI.array()+n).matrix(), AeqJ;
  171. new_AJ << AJ, AeqJ, (AeqI.array()+n).matrix();
  172. new_AV << AV, AeqV, AeqV;
  173. igl::sparse(new_AI,new_AJ,new_AV,n+neq,n+neq,new_A);
  174. // precompute RHS builders
  175. if(kr > 0)
  176. {
  177. Eigen::SparseMatrix<T> Aulk;
  178. igl::slice(new_A,data.unknown_lagrange,data.known,Aulk);
  179. Eigen::SparseMatrix<T> Akul;
  180. igl::slice(new_A,data.known,data.unknown_lagrange,Akul);
  181. //// This doesn't work!!!
  182. //data.preY = Aulk + Akul.transpose();
  183. Eigen::SparseMatrix<T> AkulT = Akul.transpose();
  184. data.preY = Aulk + AkulT;
  185. }else
  186. {
  187. data.preY.resize(data.unknown_lagrange.size(),0);
  188. }
  189. // Create factorization
  190. if(data.Auu_sym && data.Auu_pd)
  191. {
  192. data.sparse = true;
  193. Eigen::SparseMatrix<T> Aequ(0,0);
  194. if(neq>0)
  195. {
  196. Eigen::Matrix<int,Eigen::Dynamic,1> Aeq_rows;
  197. Aeq_rows.resize(neq);
  198. for(int i = 0;i<neq;i++)
  199. {
  200. Aeq_rows(i) = i;
  201. }
  202. igl::slice(Aeq,Aeq_rows,data.unknown,Aequ);
  203. }
  204. // Get transpose of Aequ
  205. Eigen::SparseMatrix<T> AequT = Aequ.transpose();
  206. // Compute LU decomposition
  207. bool lu_success = igl::lu_lagrange(Auu,AequT,data.L,data.U);
  208. if(!lu_success)
  209. {
  210. return false;
  211. }
  212. }else
  213. {
  214. Eigen::SparseMatrix<T> NA;
  215. igl::slice(new_A,data.unknown_lagrange,data.unknown_lagrange,NA);
  216. // Build LU decomposition of NA
  217. data.sparse = false;
  218. fprintf(stderr,
  219. "Warning: min_quad_with_fixed_precompute() resorting to dense LU\n");
  220. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NAfull;
  221. igl::full(NA,NAfull);
  222. data.lu =
  223. Eigen::FullPivLU<Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> >(NAfull);
  224. if(!data.lu.isInvertible())
  225. {
  226. fprintf(stderr,
  227. "Error: min_quad_with_fixed_precompute() LU not invertible\n");
  228. return false;
  229. }
  230. }
  231. return true;
  232. }
  233. template <typename T>
  234. inline bool igl::min_quad_with_fixed_solve(
  235. const igl::min_quad_with_fixed_data<T> & data,
  236. const Eigen::Matrix<T,Eigen::Dynamic,1> & B,
  237. const Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> & Y,
  238. const Eigen::Matrix<T,Eigen::Dynamic,1> & Beq,
  239. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> & Z)
  240. {
  241. // number of known rows
  242. int kr = data.known.size();
  243. if(kr!=0)
  244. {
  245. assert(kr == Y.rows());
  246. }
  247. // number of columns to solve
  248. int cols = Y.cols();
  249. // number of lagrange multipliers aka linear equality constraints
  250. int neq = data.lagrange.size();
  251. // append lagrange multiplier rhs's
  252. Eigen::Matrix<T,Eigen::Dynamic,1> BBeq(B.size() + Beq.size());
  253. BBeq << B, (Beq*-2.0);
  254. // Build right hand side
  255. Eigen::Matrix<T,Eigen::Dynamic,1> BBequl;
  256. igl::slice(BBeq,data.unknown_lagrange,BBequl);
  257. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> BBequlcols;
  258. igl::repmat(BBequl,1,cols,BBequlcols);
  259. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  260. if(kr == 0)
  261. {
  262. NB = BBequlcols;
  263. }else
  264. {
  265. NB = data.preY * Y + BBequlcols;
  266. }
  267. // resize output
  268. Z.resize(data.n,cols);
  269. // Set known values
  270. for(int i = 0;i < kr;i++)
  271. {
  272. for(int j = 0;j < cols;j++)
  273. {
  274. Z(data.known(i),j) = Y(i,j);
  275. }
  276. }
  277. //std::cout<<"NB=["<<std::endl<<NB<<std::endl<<"];"<<std::endl;
  278. if(data.sparse)
  279. {
  280. //std::cout<<"data.LIJV=["<<std::endl;print_ijv(data.L,1);std::cout<<std::endl<<"];"<<
  281. // std::endl<<"data.L=sparse(data.LIJV(:,1),data.LIJV(:,2),data.LIJV(:,3),"<<
  282. // data.L.rows()<<","<<data.L.cols()<<");"<<std::endl;
  283. //std::cout<<"data.UIJV=["<<std::endl;print_ijv(data.U,1);std::cout<<std::endl<<"];"<<
  284. // std::endl<<"data.U=sparse(data.UIJV(:,1),data.UIJV(:,2),data.UIJV(:,3),"<<
  285. // data.U.rows()<<","<<data.U.cols()<<");"<<std::endl;
  286. data.L.template triangularView<Eigen::Lower>().solveInPlace(NB);
  287. data.U.template triangularView<Eigen::Upper>().solveInPlace(NB);
  288. }else
  289. {
  290. NB = data.lu.solve(NB);
  291. }
  292. // Now NB contains sol/-0.5
  293. NB *= -0.5;
  294. // Now NB contains solution
  295. // Place solution in Z
  296. for(int i = 0;i<(NB.rows()-neq);i++)
  297. {
  298. for(int j = 0;j<NB.cols();j++)
  299. {
  300. Z(data.unknown_lagrange(i),j) = NB(i,j);
  301. }
  302. }
  303. return true;
  304. }
  305. #endif