min_quad_with_fixed.h 3.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113
  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. // MIN_QUAD_WITH_FIXED Minimize quadratic energy Z'*A*Z + Z'*B + C with
  17. // constraints that Z(known) = Y, optionally also subject to the constraints
  18. // Aeq*Z = Beq
  19. //
  20. // Templates:
  21. // T should be a eigen matrix primitive type like int or double
  22. // Inputs:
  23. // A n by n matrix of quadratic coefficients
  24. // known list of indices to known rows in Z
  25. // Y list of fixed values corresponding to known rows in Z
  26. // Aeq m by n list of linear equality constraint coefficients
  27. // pd flag specifying whether A(unknown,unknown) is positive definite
  28. // Outputs:
  29. // data factorization struct with all necessary information to solve
  30. // using min_quad_with_fixed_solve
  31. // Returns true on success, false on error
  32. //
  33. // Benchmark: For a harmonic solve on a mesh with 325K facets, matlab 2.2
  34. // secs, igl/min_quad_with_fixed.h 7.1 secs
  35. //
  36. template <typename T, typename Derivedknown>
  37. IGL_INLINE bool min_quad_with_fixed_precompute(
  38. const Eigen::SparseMatrix<T>& A,
  39. const Eigen::PlainObjectBase<Derivedknown> & known,
  40. const Eigen::SparseMatrix<T>& Aeq,
  41. const bool pd,
  42. min_quad_with_fixed_data<T> & data
  43. );
  44. // Solves a system previously factored using min_quad_with_fixed_precompute
  45. //
  46. // Template:
  47. // T type of sparse matrix (e.g. double)
  48. // DerivedY type of Y (e.g. derived from VectorXd or MatrixXd)
  49. // DerivedZ type of Z (e.g. derived from VectorXd or MatrixXd)
  50. // Inputs:
  51. // data factorization struct with all necessary precomputation to solve
  52. // B n by 1 column of linear coefficients
  53. // Beq m by 1 list of linear equality constraint constant values
  54. // Outputs:
  55. // Z n by cols solution
  56. // Returns true on success, false on error
  57. template <
  58. typename T,
  59. typename DerivedB,
  60. typename DerivedY,
  61. typename DerivedBeq,
  62. typename DerivedZ>
  63. IGL_INLINE bool min_quad_with_fixed_solve(
  64. const min_quad_with_fixed_data<T> & data,
  65. const Eigen::PlainObjectBase<DerivedB> & B,
  66. const Eigen::PlainObjectBase<DerivedY> & Y,
  67. const Eigen::PlainObjectBase<DerivedBeq> & Beq,
  68. Eigen::PlainObjectBase<DerivedZ> & Z);
  69. }
  70. template <typename T>
  71. struct igl::min_quad_with_fixed_data
  72. {
  73. // Size of original system: number of unknowns + number of knowns
  74. int n;
  75. // Whether A(unknown,unknown) is positive definite
  76. bool Auu_pd;
  77. // Whether A(unknown,unknown) is symmetric
  78. bool Auu_sym;
  79. // Indices of known variables
  80. Eigen::VectorXi known;
  81. // Indices of unknown variables
  82. Eigen::VectorXi unknown;
  83. // Indices of lagrange variables
  84. Eigen::VectorXi lagrange;
  85. // Indices of unknown variable followed by Indices of lagrange variables
  86. Eigen::VectorXi unknown_lagrange;
  87. // Matrix multiplied against Y when constructing right hand side
  88. Eigen::SparseMatrix<T> preY;
  89. enum SolverType
  90. {
  91. LLT = 0,
  92. LDLT = 1,
  93. LU = 2,
  94. NUM_SOLVER_TYPES = 3
  95. } solver_type;
  96. // Solvers
  97. Eigen::SimplicialLLT <Eigen::SparseMatrix<T > > llt;
  98. Eigen::SimplicialLDLT<Eigen::SparseMatrix<T > > ldlt;
  99. Eigen::SparseLU<Eigen::SparseMatrix<T, Eigen::ColMajor>, Eigen::COLAMDOrdering<int> > lu;
  100. // Debug
  101. Eigen::SparseMatrix<T> NA;
  102. Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic> NB;
  103. };
  104. #ifdef IGL_HEADER_ONLY
  105. # include "min_quad_with_fixed.cpp"
  106. #endif
  107. #endif