eigs.cpp 4.7 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139
  1. #include "eigs.h"
  2. #include "../sort.h"
  3. #include "../IGL_ASSERT.h"
  4. #include "../placeholders.h"
  5. #include <Spectra/SymGEigsShiftSolver.h>
  6. template <
  7. typename EigsScalar,
  8. typename DerivedU,
  9. typename DerivedS,
  10. typename Solver>
  11. IGL_INLINE bool igl::spectra::eigs(
  12. const Eigen::SparseMatrix<EigsScalar> & A,
  13. const Eigen::SparseMatrix<EigsScalar> & B,
  14. const int k,
  15. const igl::EigsType type,
  16. Eigen::PlainObjectBase<DerivedU> & U,
  17. Eigen::PlainObjectBase<DerivedS> & S)
  18. {
  19. IGL_ASSERT(k > 0 && "k should be positive");
  20. IGL_ASSERT(k < A.rows() && "k should be less than size of A");
  21. IGL_ASSERT(type == igl::EIGS_TYPE_SM && "Only SM supported");
  22. // This seems like a hack. For the "eigs: grid" test this is necessary to get
  23. // at least 1e-4 error for the first 5 eigen values. It's annoying that this
  24. // means that the zero modes become O(sigma) and this is now rather large.
  25. //
  26. // I wonder if this is an issue with SparseLU and if UMFPACK would be better.
  27. //
  28. // Ideally this value would be 0.
  29. const EigsScalar sigma = 1e-8;
  30. return igl::spectra::eigs(A,B,k,sigma,U,S);
  31. }
  32. template <
  33. typename EigsScalar,
  34. typename DerivedU,
  35. typename DerivedS,
  36. typename Solver>
  37. IGL_INLINE bool igl::spectra::eigs(
  38. const Eigen::SparseMatrix<EigsScalar> & A,
  39. const Eigen::SparseMatrix<EigsScalar> & B,
  40. const int k,
  41. const EigsScalar sigma,
  42. Eigen::PlainObjectBase<DerivedU> & U,
  43. Eigen::PlainObjectBase<DerivedS> & S)
  44. {
  45. IGL_ASSERT(k > 0 && "k should be positive");
  46. IGL_ASSERT(k < A.rows() && "k should be less than size of A");
  47. class SparseMatProd
  48. {
  49. public:
  50. using Scalar = EigsScalar;
  51. const Eigen::SparseMatrix<Scalar> & m_B;
  52. SparseMatProd(const Eigen::SparseMatrix<Scalar> & B) : m_B(B) {}
  53. int rows() const { return m_B.rows(); }
  54. int cols() const { return m_B.cols(); }
  55. void perform_op(const Scalar *x_in, Scalar *y_out) const
  56. {
  57. typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1> VectorXS;
  58. Eigen::Map<const VectorXS> x(x_in, m_B.cols());
  59. Eigen::Map< VectorXS> y(y_out, m_B.rows());
  60. y = m_B * x;
  61. }
  62. };
  63. // Solver must expose .compute(A) and .solve(x)
  64. class ShiftInvert
  65. {
  66. public:
  67. using Scalar = EigsScalar;
  68. private:
  69. const Eigen::SparseMatrix<Scalar> & m_A;
  70. const Eigen::SparseMatrix<Scalar> & m_B;
  71. Scalar m_sigma;
  72. Solver m_solver;
  73. public:
  74. bool m_solver_is_successfully_factorized;
  75. ShiftInvert(
  76. const Eigen::SparseMatrix<Scalar>& A,
  77. const Eigen::SparseMatrix<Scalar>& B,
  78. const Scalar sigma):
  79. m_A(A), m_B(B)
  80. {
  81. IGL_ASSERT(m_A.rows() == m_A.cols() && "A must be square");
  82. IGL_ASSERT(m_B.rows() == m_B.cols() && "B must be square");
  83. IGL_ASSERT(m_A.rows() == m_B.cols() && "A and B must have the same size");
  84. set_shift(sigma, true);
  85. }
  86. void set_shift(const Scalar & sigma, const bool force = false)
  87. {
  88. if(sigma == m_sigma && !force)
  89. {
  90. return;
  91. }
  92. m_sigma = sigma;
  93. const Eigen::SparseMatrix<Scalar> C = m_A + m_sigma * m_B;
  94. m_solver.compute(C);
  95. m_solver_is_successfully_factorized = (m_solver.info() == Eigen::Success);
  96. }
  97. int rows() const { return m_A.rows(); }
  98. int cols() const { return m_A.cols(); }
  99. void perform_op(const Scalar* x_in,Scalar* y_out) const
  100. {
  101. typedef Eigen::Matrix<Scalar, Eigen::Dynamic, 1> VectorXS;
  102. Eigen::Map<const VectorXS>x(x_in, m_A.cols());
  103. Eigen::Map<VectorXS>y(y_out, m_A.rows());
  104. y = m_solver.solve(x);
  105. }
  106. };
  107. SparseMatProd Bop(B);
  108. ShiftInvert op(A, B, sigma);
  109. if(!op.m_solver_is_successfully_factorized)
  110. {
  111. return false;
  112. }
  113. Spectra::SymGEigsShiftSolver<ShiftInvert, SparseMatProd, Spectra::GEigsMode::ShiftInvert> geigs(op, Bop, k, 2*k, sigma);
  114. geigs.init();
  115. geigs.compute(Spectra::SortRule::LargestMagn);
  116. if (geigs.info() != Spectra::CompInfo::Successful)
  117. {
  118. return false;
  119. }
  120. U = geigs.eigenvectors().template cast<typename DerivedU::Scalar>();
  121. S = geigs.eigenvalues().template cast<typename DerivedS::Scalar>();
  122. Eigen::VectorXi I;
  123. igl::sort( Eigen::VectorXd(S), 1, false, S, I);
  124. U = U(igl::placeholders::all,I).eval();
  125. return true;
  126. }
  127. #ifdef IGL_STATIC_LIBRARY
  128. // Explicit template instantiation
  129. // generated by autoexplicit.sh
  130. template bool igl::spectra::eigs<double, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1>, Eigen::SparseLU<Eigen::SparseMatrix<double, 0, int>, Eigen::COLAMDOrdering<int> > >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::SparseMatrix<double, 0, int> const&, int , igl::EigsType, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  131. #endif