12 #ifndef EIGEN_MINRES_H_
13 #define EIGEN_MINRES_H_
29 template<
typename MatrixType,
typename Rhs,
typename Dest,
typename Preconditioner>
31 void minres(
const MatrixType& mat,
const Rhs& rhs, Dest& x,
32 const Preconditioner& precond,
int& iters,
33 typename Dest::RealScalar& tol_error)
36 typedef typename Dest::RealScalar RealScalar;
37 typedef typename Dest::Scalar Scalar;
38 typedef Matrix<Scalar,Dynamic,1> VectorType;
41 const RealScalar rhsNorm2(rhs.squaredNorm());
51 const int maxIters(iters);
52 const int N(mat.cols());
53 const RealScalar threshold2(tol_error*tol_error*rhsNorm2);
57 VectorType v( VectorType::Zero(N) );
58 VectorType v_new(rhs-mat*x);
59 RealScalar residualNorm2(v_new.squaredNorm());
61 VectorType w_new(precond.solve(v_new));
63 RealScalar beta_new2(v_new.dot(w_new));
64 eigen_assert(beta_new2 >= 0.0 &&
"PRECONDITIONER IS NOT POSITIVE DEFINITE");
65 RealScalar beta_new(sqrt(beta_new2));
66 const RealScalar beta_one(beta_new);
71 RealScalar c_old(1.0);
73 RealScalar s_old(0.0);
75 VectorType p_old(VectorType::Zero(N));
80 while ( iters < maxIters )
92 const RealScalar beta(beta_new);
98 v_new.noalias() = mat*w - beta*v_old;
99 const RealScalar alpha = v_new.dot(w);
101 w_new = precond.solve(v_new);
102 beta_new2 = v_new.dot(w_new);
103 eigen_assert(beta_new2 >= 0.0 &&
"PRECONDITIONER IS NOT POSITIVE DEFINITE");
104 beta_new = sqrt(beta_new2);
109 const RealScalar r2 =s*alpha+c*c_old*beta;
110 const RealScalar r3 =s_old*beta;
111 const RealScalar r1_hat=c*alpha-c_old*s*beta;
112 const RealScalar r1 =sqrt( std::pow(r1_hat,2) + std::pow(beta_new,2) );
122 p.noalias()=(w-r2*p_old-r3*p_oold) /r1;
123 x += beta_one*c*eta*p;
127 residualNorm2 *= s*s;
129 if ( residualNorm2 < threshold2)
140 tol_error = std::sqrt(residualNorm2 / rhsNorm2);
145 template<
typename _MatrixType,
int _UpLo=Lower,
146 typename _Preconditioner = IdentityPreconditioner>
152 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
153 struct traits<
MINRES<_MatrixType,_UpLo,_Preconditioner> >
155 typedef _MatrixType MatrixType;
156 typedef _Preconditioner Preconditioner;
197 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
198 class MINRES :
public IterativeSolverBase<MINRES<_MatrixType,_UpLo,_Preconditioner> >
201 typedef IterativeSolverBase<MINRES> Base;
202 using Base::mp_matrix;
204 using Base::m_iterations;
206 using Base::m_isInitialized;
208 typedef _MatrixType MatrixType;
209 typedef typename MatrixType::Scalar Scalar;
210 typedef typename MatrixType::Index Index;
211 typedef typename MatrixType::RealScalar RealScalar;
212 typedef _Preconditioner Preconditioner;
231 template<
typename MatrixDerived>
232 explicit MINRES(
const EigenBase<MatrixDerived>& A) : Base(A.derived()) {}
242 template<
typename Rhs,
typename Guess>
243 inline const internal::solve_retval_with_guess<MINRES, Rhs, Guess>
246 eigen_assert(m_isInitialized &&
"MINRES is not initialized.");
247 eigen_assert(Base::rows()==b.rows()
248 &&
"MINRES::solve(): invalid number of rows of the right hand side matrix b");
249 return internal::solve_retval_with_guess
250 <
MINRES, Rhs, Guess>(*
this, b.derived(), x0);
254 template<
typename Rhs,
typename Dest>
255 void _solveWithGuess(
const Rhs& b, Dest& x)
const
257 typedef typename internal::conditional<UpLo==(Lower|Upper),
259 SparseSelfAdjointView<const MatrixType, UpLo>
260 >::type MatrixWrapperType;
262 m_iterations = Base::maxIterations();
263 m_error = Base::m_tolerance;
265 for(
int j=0; j<b.cols(); ++j)
267 m_iterations = Base::maxIterations();
268 m_error = Base::m_tolerance;
270 typename Dest::ColXpr xj(x,j);
271 internal::minres(MatrixWrapperType(*mp_matrix), b.col(j), xj,
272 Base::m_preconditioner, m_iterations, m_error);
275 m_isInitialized =
true;
276 m_info = m_error <= Base::m_tolerance ? Success : NoConvergence;
280 template<
typename Rhs,
typename Dest>
281 void _solve(
const Rhs& b, Dest& x)
const
284 _solveWithGuess(b,x);
293 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner,
typename Rhs>
294 struct solve_retval<MINRES<_MatrixType,_UpLo,_Preconditioner>, Rhs>
295 : solve_retval_base<MINRES<_MatrixType,_UpLo,_Preconditioner>, Rhs>
297 typedef MINRES<_MatrixType,_UpLo,_Preconditioner> Dec;
298 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
300 template<typename Dest>
void evalTo(Dest& dst)
const
302 dec()._solve(rhs(),dst);
310 #endif // EIGEN_MINRES_H
const internal::solve_retval_with_guess< MINRES, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: MINRES.h:244
~MINRES()
Definition: MINRES.h:235
A minimal residual solver for sparse symmetric problems.
Definition: MINRES.h:148
MINRES(const EigenBase< MatrixDerived > &A)
Definition: MINRES.h:232
MINRES()
Definition: MINRES.h:219