min_quad_with_fixed.h 5.9 KB

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