00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012 #ifndef EIGEN_MINRES_H_
00013 #define EIGEN_MINRES_H_
00014
00015
00016 namespace Eigen {
00017
00018 namespace internal {
00019
00029 template<typename MatrixType, typename Rhs, typename Dest, typename Preconditioner>
00030 EIGEN_DONT_INLINE
00031 void minres(const MatrixType& mat, const Rhs& rhs, Dest& x,
00032 const Preconditioner& precond, int& iters,
00033 typename Dest::RealScalar& tol_error)
00034 {
00035 using std::sqrt;
00036 typedef typename Dest::RealScalar RealScalar;
00037 typedef typename Dest::Scalar Scalar;
00038 typedef Matrix<Scalar,Dynamic,1> VectorType;
00039
00040
00041 const int maxIters(iters);
00042 const int N(mat.cols());
00043 const RealScalar rhsNorm2(rhs.squaredNorm());
00044 const RealScalar threshold2(tol_error*tol_error*rhsNorm2);
00045
00046
00047
00048 VectorType v( VectorType::Zero(N) );
00049 VectorType v_new(rhs-mat*x);
00050 RealScalar residualNorm2(v_new.squaredNorm());
00051
00052 VectorType w_new(precond.solve(v_new));
00053
00054 RealScalar beta_new2(v_new.dot(w_new));
00055 eigen_assert(beta_new2 >= 0 && "PRECONDITIONER IS NOT POSITIVE DEFINITE");
00056 RealScalar beta_new(sqrt(beta_new2));
00057 const RealScalar beta_one(beta_new);
00058 v_new /= beta_new;
00059 w_new /= beta_new;
00060
00061 RealScalar c(1.0);
00062 RealScalar c_old(1.0);
00063 RealScalar s(0.0);
00064 RealScalar s_old(0.0);
00065
00066 VectorType p_old(VectorType::Zero(N));
00067 VectorType p(p_old);
00068 RealScalar eta(1.0);
00069
00070 iters = 0;
00071 while ( iters < maxIters ){
00072
00073
00074
00075
00076
00077
00078
00079
00080
00081
00082
00083 const RealScalar beta(beta_new);
00084
00085 const VectorType v_old(v);
00086 v = v_new;
00087
00088 const VectorType w(w_new);
00089 v_new.noalias() = mat*w - beta*v_old;
00090 const RealScalar alpha = v_new.dot(w);
00091 v_new -= alpha*v;
00092 w_new = precond.solve(v_new);
00093 beta_new2 = v_new.dot(w_new);
00094 eigen_assert(beta_new2 >= 0 && "PRECONDITIONER IS NOT POSITIVE DEFINITE");
00095 beta_new = sqrt(beta_new2);
00096 v_new /= beta_new;
00097 w_new /= beta_new;
00098
00099
00100 const RealScalar r2 =s*alpha+c*c_old*beta;
00101 const RealScalar r3 =s_old*beta;
00102 const RealScalar r1_hat=c*alpha-c_old*s*beta;
00103 const RealScalar r1 =sqrt( std::pow(r1_hat,2) + std::pow(beta_new,2) );
00104 c_old = c;
00105 s_old = s;
00106 c=r1_hat/r1;
00107 s=beta_new/r1;
00108
00109
00110
00111 const VectorType p_oold(p_old);
00112 p_old = p;
00113 p.noalias()=(w-r2*p_old-r3*p_oold) /r1;
00114 x += beta_one*c*eta*p;
00115 residualNorm2 *= s*s;
00116
00117 if ( residualNorm2 < threshold2){
00118 break;
00119 }
00120
00121 eta=-s*eta;
00122 iters++;
00123 }
00124 tol_error = std::sqrt(residualNorm2 / rhsNorm2);
00125 }
00126
00127 }
00128
00129 template< typename _MatrixType, int _UpLo=Lower,
00130 typename _Preconditioner = IdentityPreconditioner>
00131
00132 class MINRES;
00133
00134 namespace internal {
00135
00136 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
00137 struct traits<MINRES<_MatrixType,_UpLo,_Preconditioner> >
00138 {
00139 typedef _MatrixType MatrixType;
00140 typedef _Preconditioner Preconditioner;
00141 };
00142
00143 }
00144
00194 template< typename _MatrixType, int _UpLo, typename _Preconditioner>
00195 class MINRES : public IterativeSolverBase<MINRES<_MatrixType,_UpLo,_Preconditioner> >
00196 {
00197
00198 typedef IterativeSolverBase<MINRES> Base;
00199 using Base::mp_matrix;
00200 using Base::m_error;
00201 using Base::m_iterations;
00202 using Base::m_info;
00203 using Base::m_isInitialized;
00204 public:
00205 typedef _MatrixType MatrixType;
00206 typedef typename MatrixType::Scalar Scalar;
00207 typedef typename MatrixType::Index Index;
00208 typedef typename MatrixType::RealScalar RealScalar;
00209 typedef _Preconditioner Preconditioner;
00210
00211 enum {UpLo = _UpLo};
00212
00213 public:
00214
00216 MINRES() : Base() {}
00217
00228 MINRES(const MatrixType& A) : Base(A) {}
00229
00231 ~MINRES(){}
00232
00238 template<typename Rhs,typename Guess>
00239 inline const internal::solve_retval_with_guess<MINRES, Rhs, Guess>
00240 solveWithGuess(const MatrixBase<Rhs>& b, const Guess& x0) const
00241 {
00242 eigen_assert(m_isInitialized && "MINRES is not initialized.");
00243 eigen_assert(Base::rows()==b.rows()
00244 && "MINRES::solve(): invalid number of rows of the right hand side matrix b");
00245 return internal::solve_retval_with_guess
00246 <MINRES, Rhs, Guess>(*this, b.derived(), x0);
00247 }
00248
00250 template<typename Rhs,typename Dest>
00251 void _solveWithGuess(const Rhs& b, Dest& x) const
00252 {
00253 m_iterations = Base::maxIterations();
00254 m_error = Base::m_tolerance;
00255
00256 for(int j=0; j<b.cols(); ++j)
00257 {
00258 m_iterations = Base::maxIterations();
00259 m_error = Base::m_tolerance;
00260
00261 typename Dest::ColXpr xj(x,j);
00262 internal::minres(mp_matrix->template selfadjointView<UpLo>(), b.col(j), xj,
00263 Base::m_preconditioner, m_iterations, m_error);
00264 }
00265
00266 m_isInitialized = true;
00267 m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
00268 }
00269
00271 template<typename Rhs,typename Dest>
00272 void _solve(const Rhs& b, Dest& x) const
00273 {
00274 x.setZero();
00275 _solveWithGuess(b,x);
00276 }
00277
00278 protected:
00279
00280 };
00281
00282 namespace internal {
00283
00284 template<typename _MatrixType, int _UpLo, typename _Preconditioner, typename Rhs>
00285 struct solve_retval<MINRES<_MatrixType,_UpLo,_Preconditioner>, Rhs>
00286 : solve_retval_base<MINRES<_MatrixType,_UpLo,_Preconditioner>, Rhs>
00287 {
00288 typedef MINRES<_MatrixType,_UpLo,_Preconditioner> Dec;
00289 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
00290
00291 template<typename Dest> void evalTo(Dest& dst) const
00292 {
00293 dec()._solve(rhs(),dst);
00294 }
00295 };
00296
00297 }
00298
00299 }
00300
00301 #endif // EIGEN_MINRES_H
00302