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;
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>
155 typedef _MatrixType MatrixType;
156 typedef _Preconditioner Preconditioner;
197 template<
typename _MatrixType,
int _UpLo,
typename _Preconditioner>
202 using Base::mp_matrix;
204 using Base::m_iterations;
206 using Base::m_isInitialized;
208 typedef _MatrixType MatrixType;
210 typedef typename MatrixType::Index Index;
211 typedef typename MatrixType::RealScalar RealScalar;
212 typedef _Preconditioner Preconditioner;
231 template<
typename MatrixDerived>
242 template<
typename Rhs,
typename 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");
250 <
MINRES, Rhs, Guess>(*
this, b.derived(), x0);
254 template<
typename Rhs,
typename Dest>
255 void _solveWithGuess(
const Rhs& b, Dest& x)
const 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;
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>
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 Definition: gtest_unittest.cc:5031
Definition: ForwardDeclarations.h:124
iterative scaling algorithm to equilibrate rows and column norms in matrices
Definition: TestIMU_Common.h:87
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
Definition: SparseSelfAdjointView.h:49
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:26
MINRES(const EigenBase< MatrixDerived > &A)
Initialize the solver with matrix A for further Ax=b solving.
Definition: MINRES.h:232
View matrix as an upper triangular matrix.
Definition: Constants.h:169
MINRES()
Default constructor.
Definition: MINRES.h:219
A minimal residual solver for sparse symmetric problems.
Definition: MINRES.h:148
Definition: SparseSolve.h:86
Computation was successful.
Definition: Constants.h:376
Definition: BandTriangularSolver.h:13
~MINRES()
Destructor.
Definition: MINRES.h:235
Definition: ForwardDeclarations.h:125
View matrix as a lower triangular matrix.
Definition: Constants.h:167
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:21
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
Iterative procedure did not converge.
Definition: Constants.h:380
Definition: ForwardDeclarations.h:17
const internal::solve_retval_with_guess< MINRES, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition: MINRES.h:244
double Scalar
Common scalar type.
Definition: FlexibleKalmanBase.h:48