min_quad_with_fixed.h 4.4 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134
  1. #ifndef IGL_MIN_QUAD_WITH_FIXED_H
  2. #define IGL_MIN_QUAD_WITH_FIXED_H
  3. #include "igl_inline.h"
  4. #define EIGEN_YES_I_KNOW_SPARSE_MODULE_IS_NOT_STABLE_YET
  5. #include <Eigen/Core>
  6. #include <Eigen/Dense>
  7. #include <Eigen/Sparse>
  8. //#include <Eigen/SparseExtra>
  9. // Bug in unsupported/Eigen/SparseExtra needs iostream first
  10. #include <iostream>
  11. #include <unsupported/Eigen/SparseExtra>
  12. namespace igl
  13. {
  14. template <typename T>
  15. struct min_quad_with_fixed_data;
  16. // Known Bugs: rows of Aeq **must** be linearly independent. Should be using
  17. // QR decomposition otherwise:
  18. // http://www.okstate.edu/sas/v8/sashtml/ormp/chap5/sect32.htm
  19. //
  20. // MIN_QUAD_WITH_FIXED Minimize quadratic energy Z'*A*Z + Z'*B + C with
  21. // constraints that Z(known) = Y, optionally also subject to the constraints
  22. // Aeq*Z = Beq
  23. //
  24. // Templates:
  25. // T should be a eigen matrix primitive type like int or double
  26. // Inputs:
  27. // A n by n matrix of quadratic coefficients
  28. // known list of indices to known rows in Z
  29. // Y list of fixed values corresponding to known rows in Z
  30. // Aeq m by n list of linear equality constraint coefficients
  31. // pd flag specifying whether A(unknown,unknown) is positive definite
  32. // Outputs:
  33. // data factorization struct with all necessary information to solve
  34. // using min_quad_with_fixed_solve
  35. // Returns true on success, false on error
  36. //
  37. // Benchmark: For a harmonic solve on a mesh with 325K facets, matlab 2.2
  38. // secs, igl/min_quad_with_fixed.h 7.1 secs
  39. //
  40. template <typename T, typename Derivedknown>
  41. IGL_INLINE bool min_quad_with_fixed_precompute(
  42. const Eigen::SparseMatrix<T>& A,
  43. const Eigen::PlainObjectBase<Derivedknown> & known,
  44. const Eigen::SparseMatrix<T>& Aeq,
  45. const bool pd,
  46. min_quad_with_fixed_data<T> & data
  47. );
  48. // Solves a system previously factored using min_quad_with_fixed_precompute
  49. //
  50. // Template:
  51. // T type of sparse matrix (e.g. double)
  52. // DerivedY type of Y (e.g. derived from VectorXd or MatrixXd)
  53. // DerivedZ type of Z (e.g. derived from VectorXd or MatrixXd)
  54. // Inputs:
  55. // data factorization struct with all necessary precomputation to solve
  56. // B n by 1 column of linear coefficients
  57. // Beq m by 1 list of linear equality constraint constant values
  58. // Outputs:
  59. // Z n by cols solution
  60. // sol #unknowns+#lagrange by cols solution to linear system
  61. // Returns true on success, false on error
  62. template <
  63. typename T,
  64. typename DerivedB,
  65. typename DerivedY,
  66. typename DerivedBeq,
  67. typename DerivedZ,
  68. typename Derivedsol>
  69. IGL_INLINE bool min_quad_with_fixed_solve(
  70. const min_quad_with_fixed_data<T> & data,
  71. const Eigen::PlainObjectBase<DerivedB> & B,
  72. const Eigen::PlainObjectBase<DerivedY> & Y,
  73. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  74. Eigen::PlainObjectBase<DerivedZ> & Z,
  75. Eigen::PlainObjectBase<Derivedsol> & sol);
  76. // Wrapper without sol
  77. template <
  78. typename T,
  79. typename DerivedB,
  80. typename DerivedY,
  81. typename DerivedBeq,
  82. typename DerivedZ>
  83. IGL_INLINE bool min_quad_with_fixed_solve(
  84. const min_quad_with_fixed_data<T> & data,
  85. const Eigen::PlainObjectBase<DerivedB> & B,
  86. const Eigen::PlainObjectBase<DerivedY> & Y,
  87. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  88. Eigen::PlainObjectBase<DerivedZ> & Z);
  89. }
  90. template <typename T>
  91. struct igl::min_quad_with_fixed_data
  92. {
  93. // Size of original system: number of unknowns + number of knowns
  94. int n;
  95. // Whether A(unknown,unknown) is positive definite
  96. bool Auu_pd;
  97. // Whether A(unknown,unknown) is symmetric
  98. bool Auu_sym;
  99. // Indices of known variables
  100. Eigen::VectorXi known;
  101. // Indices of unknown variables
  102. Eigen::VectorXi unknown;
  103. // Indices of lagrange variables
  104. Eigen::VectorXi lagrange;
  105. // Indices of unknown variable followed by Indices of lagrange variables
  106. Eigen::VectorXi unknown_lagrange;
  107. // Matrix multiplied against Y when constructing right hand side
  108. Eigen::SparseMatrix<T> preY;
  109. enum SolverType
  110. {
  111. LLT = 0,
  112. LDLT = 1,
  113. LU = 2,
  114. NUM_SOLVER_TYPES = 3
  115. } solver_type;
  116. // Solvers
  117. Eigen::SimplicialLLT <Eigen::SparseMatrix<T > > llt;
  118. Eigen::SimplicialLDLT<Eigen::SparseMatrix<T > > ldlt;
  119. Eigen::SparseLU<Eigen::SparseMatrix<T, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > lu;
  120. // Debug
  121. Eigen::SparseMatrix<T> NA;
  122. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  123. };
  124. #ifdef IGL_HEADER_ONLY
  125. # include "min_quad_with_fixed.cpp"
  126. #endif
  127. #endif