eigs.cpp 4.5 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158
  1. #include "eigs.h"
  2. #include "cotmatrix.h"
  3. #include "sort.h"
  4. #include "slice.h"
  5. #include "massmatrix.h"
  6. #include <iostream>
  7. template <
  8. typename Atype,
  9. typename Btype,
  10. typename DerivedU,
  11. typename DerivedS>
  12. IGL_INLINE bool igl::eigs(
  13. const Eigen::SparseMatrix<Atype> & A,
  14. const Eigen::SparseMatrix<Btype> & iB,
  15. const size_t k,
  16. const EigsType type,
  17. Eigen::PlainObjectBase<DerivedU> & sU,
  18. Eigen::PlainObjectBase<DerivedS> & sS)
  19. {
  20. using namespace Eigen;
  21. using namespace std;
  22. const size_t n = A.rows();
  23. assert(A.cols() == n && "A should be square.");
  24. assert(iB.rows() == n && "B should be match A's dims.");
  25. assert(iB.cols() == n && "B should be square.");
  26. assert(type == EIGS_TYPE_SM && "Only low frequencies are supported");
  27. DerivedU U(n,k);
  28. DerivedS S(k,1);
  29. typedef Atype Scalar;
  30. typedef Eigen::Matrix<typename DerivedU::Scalar,DerivedU::RowsAtCompileTime,1> VectorXS;
  31. // Rescale B for better numerics
  32. const Scalar rescale = std::abs(iB.diagonal().maxCoeff());
  33. const Eigen::SparseMatrix<Btype> B = iB/rescale;
  34. Scalar tol = 1e-4;
  35. Scalar conv = 1e-14;
  36. int max_iter = 100;
  37. int i = 0;
  38. while(true)
  39. {
  40. // Random initial guess
  41. VectorXS y = VectorXS::Random(n,1);
  42. Scalar eff_sigma = 0;
  43. if(i>0)
  44. {
  45. eff_sigma = 1e-8+std::abs(S(i-1));
  46. }
  47. // whether to use rayleigh quotient method
  48. bool ray = false;
  49. Scalar err = std::numeric_limits<Scalar>::infinity();
  50. int iter;
  51. Scalar sigma = std::numeric_limits<Scalar>::infinity();
  52. VectorXS x;
  53. for(iter = 0;iter<max_iter;iter++)
  54. {
  55. if(i>0 && !ray)
  56. {
  57. // project-out existing modes
  58. for(int j = 0;j<i;j++)
  59. {
  60. const VectorXS u = U.col(j);
  61. y = (y - u*u.dot(B*y)/u.dot(B * u)).eval();
  62. }
  63. }
  64. // normalize
  65. x = y/sqrt(y.dot(B*y));
  66. // current guess at eigen value
  67. sigma = x.dot(A*x)/x.dot(B*x);
  68. //x *= sigma>0?1.:-1.;
  69. Scalar err_prev = err;
  70. err = (A*x-sigma*B*x).array().abs().maxCoeff();
  71. if(err<conv)
  72. {
  73. break;
  74. }
  75. if(ray || err<tol)
  76. {
  77. eff_sigma = sigma;
  78. ray = true;
  79. }
  80. Scalar tikhonov = std::abs(eff_sigma)<1e-12?1e-10:0;
  81. switch(type)
  82. {
  83. default:
  84. assert(false && "Not supported");
  85. break;
  86. case EIGS_TYPE_SM:
  87. {
  88. SimplicialLDLT<SparseMatrix<Scalar> > solver;
  89. const SparseMatrix<Scalar> C = A-eff_sigma*B+tikhonov*B;
  90. //mw.save(C,"C");
  91. //mw.save(eff_sigma,"eff_sigma");
  92. //mw.save(tikhonov,"tikhonov");
  93. solver.compute(C);
  94. switch(solver.info())
  95. {
  96. case Eigen::Success:
  97. break;
  98. case Eigen::NumericalIssue:
  99. cerr<<"Error: Numerical issue."<<endl;
  100. return false;
  101. default:
  102. cerr<<"Error: Other."<<endl;
  103. return false;
  104. }
  105. const VectorXS rhs = B*x;
  106. y = solver.solve(rhs);
  107. //mw.save(rhs,"rhs");
  108. //mw.save(y,"y");
  109. //mw.save(x,"x");
  110. //mw.write("eigs.mat");
  111. //if(i == 1)
  112. //return false;
  113. break;
  114. }
  115. }
  116. }
  117. if(iter == max_iter)
  118. {
  119. cerr<<"Failed to converge."<<endl;
  120. return false;
  121. }
  122. if(i==0 || (S.head(i).array()-sigma).abs().maxCoeff()>1e-14)
  123. {
  124. U.col(i) = x;
  125. S(i) = sigma;
  126. i++;
  127. if(i == k)
  128. {
  129. break;
  130. }
  131. }else
  132. {
  133. // restart with new random guess.
  134. cout<<"RESTART!"<<endl;
  135. }
  136. }
  137. // finally sort
  138. VectorXi I;
  139. igl::sort(S,1,false,sS,I);
  140. sU = igl::slice(U,I,2);
  141. sS /= rescale;
  142. sU /= sqrt(rescale);
  143. return true;
  144. }
  145. #ifdef IGL_STATIC_LIBRARY
  146. // Explicit template specialization
  147. template bool igl::eigs<double, double, Eigen::Matrix<double, -1, -1, 0, -1, -1>, Eigen::Matrix<double, -1, 1, 0, -1, 1> >(Eigen::SparseMatrix<double, 0, int> const&, Eigen::SparseMatrix<double, 0, int> const&, unsigned long, igl::EigsType, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, -1, 0, -1, -1> >&, Eigen::PlainObjectBase<Eigen::Matrix<double, -1, 1, 0, -1, 1> >&);
  148. #ifdef WIN32
  149. template bool igl::eigs<double, double, Eigen::Matrix<double,-1,-1,0,-1,-1>, Eigen::Matrix<double,-1,1,0,-1,1> >(Eigen::SparseMatrix<double,0,int> const &,Eigen::SparseMatrix<double,0,int> const &,unsigned long long, igl::EigsType, Eigen::PlainObjectBase< Eigen::Matrix<double,-1,-1,0,-1,-1> > &, Eigen::PlainObjectBase<Eigen::Matrix<double,-1,1,0,-1,1> > &);
  150. #endif
  151. #endif