min_quad_with_fixed.h 5.9 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179
  1. // This file is part of libigl, a simple c++ geometry processing library.
  2. //
  3. // Copyright (C) 2014 Alec Jacobson <alecjacobson@gmail.com>
  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. // Known Bugs: rows of Aeq **should probably** be linearly independent.
  23. // During precomputation, the rows of a Aeq are checked via QR. But in case
  24. // they're not then resulting probably will no longer be sparse: it will be
  25. // slow.
  26. //
  27. // MIN_QUAD_WITH_FIXED Minimize a quadratic energy of the form
  28. //
  29. // trace( 0.5*Z'*A*Z + Z'*B + constant )
  30. //
  31. // subject to
  32. //
  33. // Z(known,:) = Y, and
  34. // Aeq*Z = Beq
  35. //
  36. // Templates:
  37. // T should be a eigen matrix primitive type like int or double
  38. // Inputs:
  39. // A n by n matrix of quadratic coefficients
  40. // known list of indices to known rows in Z
  41. // Y list of fixed values corresponding to known rows in Z
  42. // Aeq m by n list of linear equality constraint coefficients
  43. // pd flag specifying whether A(unknown,unknown) is positive definite
  44. // Outputs:
  45. // data factorization struct with all necessary information to solve
  46. // using min_quad_with_fixed_solve
  47. // Returns true on success, false on error
  48. //
  49. // Benchmark: For a harmonic solve on a mesh with 325K facets, matlab 2.2
  50. // secs, igl/min_quad_with_fixed.h 7.1 secs
  51. //
  52. template <typename T, typename Derivedknown>
  53. IGL_INLINE bool min_quad_with_fixed_precompute(
  54. const Eigen::SparseMatrix<T>& A,
  55. const Eigen::MatrixBase<Derivedknown> & known,
  56. const Eigen::SparseMatrix<T>& Aeq,
  57. const bool pd,
  58. min_quad_with_fixed_data<T> & data
  59. );
  60. // Solves a system previously factored using min_quad_with_fixed_precompute
  61. //
  62. // Template:
  63. // T type of sparse matrix (e.g. double)
  64. // DerivedY type of Y (e.g. derived from VectorXd or MatrixXd)
  65. // DerivedZ type of Z (e.g. derived from VectorXd or MatrixXd)
  66. // Inputs:
  67. // data factorization struct with all necessary precomputation to solve
  68. // B n by k column of linear coefficients
  69. // Y b by k list of constant fixed values
  70. // Beq m by k list of linear equality constraint constant values
  71. // Outputs:
  72. // Z n by k solution
  73. // sol #unknowns+#lagrange by k solution to linear system
  74. // Returns true on success, false on error
  75. template <
  76. typename T,
  77. typename DerivedB,
  78. typename DerivedY,
  79. typename DerivedBeq,
  80. typename DerivedZ,
  81. typename Derivedsol>
  82. IGL_INLINE bool min_quad_with_fixed_solve(
  83. const min_quad_with_fixed_data<T> & data,
  84. const Eigen::MatrixBase<DerivedB> & B,
  85. const Eigen::MatrixBase<DerivedY> & Y,
  86. const Eigen::MatrixBase<DerivedBeq> & Beq,
  87. Eigen::PlainObjectBase<DerivedZ> & Z,
  88. Eigen::PlainObjectBase<Derivedsol> & sol);
  89. // Wrapper without sol
  90. template <
  91. typename T,
  92. typename DerivedB,
  93. typename DerivedY,
  94. typename DerivedBeq,
  95. typename DerivedZ>
  96. IGL_INLINE bool min_quad_with_fixed_solve(
  97. const min_quad_with_fixed_data<T> & data,
  98. const Eigen::MatrixBase<DerivedB> & B,
  99. const Eigen::MatrixBase<DerivedY> & Y,
  100. const Eigen::MatrixBase<DerivedBeq> & Beq,
  101. Eigen::PlainObjectBase<DerivedZ> & Z);
  102. template <
  103. typename T,
  104. typename Derivedknown,
  105. typename DerivedB,
  106. typename DerivedY,
  107. typename DerivedBeq,
  108. typename DerivedZ>
  109. IGL_INLINE bool min_quad_with_fixed(
  110. const Eigen::SparseMatrix<T>& A,
  111. const Eigen::MatrixBase<DerivedB> & B,
  112. const Eigen::MatrixBase<Derivedknown> & known,
  113. const Eigen::MatrixBase<DerivedY> & Y,
  114. const Eigen::SparseMatrix<T>& Aeq,
  115. const Eigen::MatrixBase<DerivedBeq> & Beq,
  116. const bool pd,
  117. Eigen::PlainObjectBase<DerivedZ> & Z);
  118. }
  119. template <typename T>
  120. struct igl::min_quad_with_fixed_data
  121. {
  122. // Size of original system: number of unknowns + number of knowns
  123. int n;
  124. // Whether A(unknown,unknown) is positive definite
  125. bool Auu_pd;
  126. // Whether A(unknown,unknown) is symmetric
  127. bool Auu_sym;
  128. // Indices of known variables
  129. Eigen::VectorXi known;
  130. // Indices of unknown variables
  131. Eigen::VectorXi unknown;
  132. // Indices of lagrange variables
  133. Eigen::VectorXi lagrange;
  134. // Indices of unknown variable followed by Indices of lagrange variables
  135. Eigen::VectorXi unknown_lagrange;
  136. // Matrix multiplied against Y when constructing right hand side
  137. Eigen::SparseMatrix<T> preY;
  138. enum SolverType
  139. {
  140. LLT = 0,
  141. LDLT = 1,
  142. LU = 2,
  143. QR_LLT = 3,
  144. NUM_SOLVER_TYPES = 4
  145. } solver_type;
  146. // Solvers
  147. Eigen::SimplicialLLT <Eigen::SparseMatrix<T > > llt;
  148. Eigen::SimplicialLDLT<Eigen::SparseMatrix<T > > ldlt;
  149. Eigen::SparseLU<Eigen::SparseMatrix<T, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > lu;
  150. // QR factorization
  151. // Are rows of Aeq linearly independent?
  152. bool Aeq_li;
  153. // Columns of Aeq corresponding to unknowns
  154. int neq;
  155. Eigen::SparseQR<Eigen::SparseMatrix<T>, Eigen::COLAMDOrdering<int> > AeqTQR;
  156. Eigen::SparseMatrix<T> Aeqk;
  157. Eigen::SparseMatrix<T> Aequ;
  158. Eigen::SparseMatrix<T> Auu;
  159. Eigen::SparseMatrix<T> AeqTQ1;
  160. Eigen::SparseMatrix<T> AeqTQ1T;
  161. Eigen::SparseMatrix<T> AeqTQ2;
  162. Eigen::SparseMatrix<T> AeqTQ2T;
  163. Eigen::SparseMatrix<T> AeqTR1;
  164. Eigen::SparseMatrix<T> AeqTR1T;
  165. Eigen::SparseMatrix<T> AeqTE;
  166. Eigen::SparseMatrix<T> AeqTET;
  167. // Debug
  168. Eigen::SparseMatrix<T> NA;
  169. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  170. };
  171. #ifndef IGL_STATIC_LIBRARY
  172. # include "min_quad_with_fixed.cpp"
  173. #endif
  174. #endif