min_quad_with_fixed.h 8.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Alec Jacobson <[email protected]>
  4. //
  5. // This Source Code Form is subject to the terms of the Mozilla Public License
  6. // v. 2.0. If a copy of the MPL was not distributed with this file, You can
  7. // obtain one at http://mozilla.org/MPL/2.0/.
  8. #ifndef IGL_MIN_QUAD_WITH_FIXED_H
  9. #define IGL_MIN_QUAD_WITH_FIXED_H
  10. #include "igl_inline.h"
  11. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  12. #include <Eigen/Core>
  13. #include <Eigen/Dense>
  14. #include <Eigen/Sparse>
  15. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  16. #include <iostream>
  17. #include <unsupported/Eigen/SparseExtra>
  18. namespace igl
  19. {
  20. template <typename T>
  21. struct min_quad_with_fixed_data;
  22. /// Minimize a convex quadratic energy subject to fixed value and linear
  23. /// equality constraints. Problems of the form
  24. ///
  25. /// trace( 0.5*Z'*A*Z + Z'*B + constant )
  26. ///
  27. /// subject to
  28. ///
  29. /// Z(known,:) = Y, and
  30. /// Aeq*Z = Beq
  31. ///
  32. /// @tparam T should be a eigen matrix primitive type like int or double
  33. /// @param[in] A n by n matrix of quadratic coefficients
  34. /// @param[in] known list of indices to known rows in Z
  35. /// @param[in] Y list of fixed values corresponding to known rows in Z
  36. /// @param[in] Aeq m by n list of linear equality constraint coefficients
  37. /// @param[in] pd flag specifying whether A(unknown,unknown) is positive definite
  38. /// @param[in,out] data factorization struct with all necessary information to solve
  39. /// using min_quad_with_fixed_solve
  40. /// @return true on success, false on error
  41. ///
  42. /// \pre rows of Aeq **should probably** be linearly independent.
  43. /// During precomputation, the rows of a Aeq are checked via QR. But in case
  44. /// they're not then resulting probably will no longer be sparse: it will be
  45. /// slow.
  46. ///
  47. template <typename T, typename Derivedknown>
  48. IGL_INLINE bool min_quad_with_fixed_precompute(
  49. const Eigen::SparseMatrix<T>& A,
  50. const Eigen::MatrixBase<Derivedknown> & known,
  51. const Eigen::SparseMatrix<T>& Aeq,
  52. const bool pd,
  53. min_quad_with_fixed_data<T> & data
  54. );
  55. /// Solves a system previously factored using min_quad_with_fixed_precompute
  56. ///
  57. /// @tparam T type of sparse matrix (e.g. double)
  58. /// @tparam DerivedY type of Y (e.g. derived from VectorXd or MatrixXd)
  59. /// @tparam DerivedZ type of Z (e.g. derived from VectorXd or MatrixXd)
  60. /// @param[in] data factorization struct with all necessary precomputation to solve
  61. /// @param[in] B n by k column of linear coefficients
  62. /// @param[in] Y b by k list of constant fixed values
  63. /// @param[in] Beq m by k list of linear equality constraint constant values
  64. /// @param[out] Z n by k solution
  65. /// @param[out] sol #unknowns+#lagrange by k solution to linear system
  66. /// Returns true on success, false on error
  67. template <
  68. typename T,
  69. typename DerivedB,
  70. typename DerivedY,
  71. typename DerivedBeq,
  72. typename DerivedZ,
  73. typename Derivedsol>
  74. IGL_INLINE bool min_quad_with_fixed_solve(
  75. const min_quad_with_fixed_data<T> & data,
  76. const Eigen::MatrixBase<DerivedB> & B,
  77. const Eigen::MatrixBase<DerivedY> & Y,
  78. const Eigen::MatrixBase<DerivedBeq> & Beq,
  79. Eigen::PlainObjectBase<DerivedZ> & Z,
  80. Eigen::PlainObjectBase<Derivedsol> & sol);
  81. /// \overload
  82. template <
  83. typename T,
  84. typename DerivedB,
  85. typename DerivedY,
  86. typename DerivedBeq,
  87. typename DerivedZ>
  88. IGL_INLINE bool min_quad_with_fixed_solve(
  89. const min_quad_with_fixed_data<T> & data,
  90. const Eigen::MatrixBase<DerivedB> & B,
  91. const Eigen::MatrixBase<DerivedY> & Y,
  92. const Eigen::MatrixBase<DerivedBeq> & Beq,
  93. Eigen::PlainObjectBase<DerivedZ> & Z);
  94. /// \overload
  95. /// \brief Minimize convex quadratic energy subject to fixed value and linear
  96. /// equality constraints. Without prefactorization.
  97. template <
  98. typename T,
  99. typename Derivedknown,
  100. typename DerivedB,
  101. typename DerivedY,
  102. typename DerivedBeq,
  103. typename DerivedZ>
  104. IGL_INLINE bool min_quad_with_fixed(
  105. const Eigen::SparseMatrix<T>& A,
  106. const Eigen::MatrixBase<DerivedB> & B,
  107. const Eigen::MatrixBase<Derivedknown> & known,
  108. const Eigen::MatrixBase<DerivedY> & Y,
  109. const Eigen::SparseMatrix<T>& Aeq,
  110. const Eigen::MatrixBase<DerivedBeq> & Beq,
  111. const bool pd,
  112. Eigen::PlainObjectBase<DerivedZ> & Z);
  113. /// Dense version optimized for very small, known at compile time sizes. Still
  114. /// works for Eigen::Dynamic (and then everything needs to be Dynamic).
  115. ///
  116. /// min_x ½ xᵀ H x + xᵀ f
  117. /// subject to
  118. /// A x = b
  119. /// x(i) = bc(i) iff k(i)==true
  120. ///
  121. /// @tparam Scalar (e.g., double)
  122. /// @tparam n #H or Eigen::Dynamic if not known at compile time
  123. /// @tparam m #A or Eigen::Dynamic if not known at compile time
  124. /// @tparam Hpd whether H is positive definite (LLT used) or not (QR used)
  125. /// @param[in] H #H by #H quadratic coefficients (only lower triangle used)
  126. /// @param[in] f #H linear coefficients
  127. /// @param[in] k #H list of flags whether to fix value
  128. /// @param[in] bc #H value to fix to (if !k(i) then bc(i) is ignored)
  129. /// @param[in] A #A by #H list of linear equality constraint coefficients, must be
  130. /// linearly independent (with self and fixed value constraints)
  131. /// @param[in] b #A list of linear equality right-hand sides
  132. /// @return #H-long solution x
  133. template <typename Scalar, int n, int m, bool Hpd=true>
  134. IGL_INLINE Eigen::Matrix<Scalar,n,1> min_quad_with_fixed(
  135. const Eigen::Matrix<Scalar,n,n> & H,
  136. const Eigen::Matrix<Scalar,n,1> & f,
  137. const Eigen::Array<bool,n,1> & k,
  138. const Eigen::Matrix<Scalar,n,1> & bc,
  139. const Eigen::Matrix<Scalar,m,n> & A,
  140. const Eigen::Matrix<Scalar,m,1> & b);
  141. /// \overload
  142. template <typename Scalar, int n, bool Hpd=true>
  143. IGL_INLINE Eigen::Matrix<Scalar,n,1> min_quad_with_fixed(
  144. const Eigen::Matrix<Scalar,n,n> & H,
  145. const Eigen::Matrix<Scalar,n,1> & f,
  146. const Eigen::Array<bool,n,1> & k,
  147. const Eigen::Matrix<Scalar,n,1> & bc);
  148. /// \overload
  149. ///
  150. /// \brief Special wrapper where the number of constrained values (i.e., true values
  151. /// in k) is exposed as a template parameter. Not intended to be called
  152. /// directly. The overhead of calling the overloads above is already minimal.
  153. template <typename Scalar, int n, int kcount, bool Hpd/*no default*/>
  154. IGL_INLINE Eigen::Matrix<Scalar,n,1> min_quad_with_fixed(
  155. const Eigen::Matrix<Scalar,n,n> & H,
  156. const Eigen::Matrix<Scalar,n,1> & f,
  157. const Eigen::Array<bool,n,1> & k,
  158. const Eigen::Matrix<Scalar,n,1> & bc);
  159. }
  160. /// Parameters and precomputed values for min_quad_with_fixed
  161. template <typename T>
  162. struct igl::min_quad_with_fixed_data
  163. {
  164. /// Size of original system: number of unknowns + number of knowns
  165. int n;
  166. /// Whether A(unknown,unknown) is positive definite
  167. bool Auu_pd;
  168. /// Whether A(unknown,unknown) is symmetric
  169. bool Auu_sym;
  170. /// Indices of known variables
  171. Eigen::VectorXi known;
  172. /// Indices of unknown variables
  173. Eigen::VectorXi unknown;
  174. /// Indices of lagrange variables
  175. Eigen::VectorXi lagrange;
  176. /// Indices of unknown variable followed by Indices of lagrange variables
  177. Eigen::VectorXi unknown_lagrange;
  178. /// Matrix multiplied against Y when constructing right hand side
  179. Eigen::SparseMatrix<T> preY;
  180. /// Type of solver used
  181. enum SolverType
  182. {
  183. LLT = 0,
  184. LDLT = 1,
  185. LU = 2,
  186. QR_LLT = 3,
  187. NUM_SOLVER_TYPES = 4
  188. } solver_type;
  189. /// Solver data (factorization)
  190. Eigen::SimplicialLLT <Eigen::SparseMatrix<T > > llt;
  191. Eigen::SimplicialLDLT<Eigen::SparseMatrix<T > > ldlt;
  192. Eigen::SparseLU<Eigen::SparseMatrix<T, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > lu;
  193. /// QR factorization
  194. /// Are rows of Aeq linearly independent?
  195. bool Aeq_li;
  196. /// Columns of Aeq corresponding to unknowns
  197. int neq;
  198. Eigen::SparseQR<Eigen::SparseMatrix<T>, Eigen::COLAMDOrdering<int> > AeqTQR;
  199. Eigen::SparseMatrix<T> Aeqk;
  200. Eigen::SparseMatrix<T> Aequ;
  201. Eigen::SparseMatrix<T> Auu;
  202. Eigen::SparseMatrix<T> AeqTQ1;
  203. Eigen::SparseMatrix<T> AeqTQ1T;
  204. Eigen::SparseMatrix<T> AeqTQ2;
  205. Eigen::SparseMatrix<T> AeqTQ2T;
  206. Eigen::SparseMatrix<T> AeqTR1;
  207. Eigen::SparseMatrix<T> AeqTR1T;
  208. Eigen::SparseMatrix<T> AeqTE;
  209. Eigen::SparseMatrix<T> AeqTET;
  210. /// @private Debug
  211. Eigen::SparseMatrix<T> NA;
  212. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  213. };
  214. #ifndef IGL_STATIC_LIBRARY
  215. # include "min_quad_with_fixed.impl.h"
  216. #endif
  217. #endif