min_quad_with_fixed.h 8.0 KB

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