eigs.cpp 5.3 KB

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