eigs.cpp 5.0 KB

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