min_quad_with_fixed.h 5.4 KB

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