A minimal residual solver for sparse symmetric problems. More...
#include <MINRES.h>
Public Types | |
enum | { UpLo = _UpLo } |
typedef MatrixType::Index | Index |
typedef _MatrixType | MatrixType |
typedef _Preconditioner | Preconditioner |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType::Scalar | Scalar |
Public Member Functions | |
template<typename Rhs , typename Dest > | |
void | _solve (const Rhs &b, Dest &x) const |
template<typename Rhs , typename Dest > | |
void | _solveWithGuess (const Rhs &b, Dest &x) const |
MINRES () | |
MINRES (const MatrixType &A) | |
template<typename Rhs , typename Guess > | |
const internal::solve_retval_with_guess < MINRES, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
~MINRES () | |
Private Types | |
typedef IterativeSolverBase < MINRES > | Base |
A minimal residual solver for sparse symmetric problems.
This class allows to solve for A.x = b sparse linear problems using the MINRES algorithm of Paige and Saunders (1975). The sparse matrix A must be symmetric (possibly indefinite). The vectors x and b can be either dense or sparse.
_MatrixType | the type of the sparse matrix A, can be a dense or a sparse matrix. |
_UpLo | the triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower. |
_Preconditioner | the type of the preconditioner. Default is DiagonalPreconditioner |
The maximal number of iterations and tolerance value can be controlled via the setMaxIterations() and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations and NumTraits<Scalar>::epsilon() for the tolerance.
This class can be used as the direct solver classes. Here is a typical usage example:
int n = 10000; VectorXd x(n), b(n); SparseMatrix<double> A(n,n); // fill A and b MINRES<SparseMatrix<double> > mr; mr.compute(A); x = mr.solve(b); std::cout << "#iterations: " << mr.iterations() << std::endl; std::cout << "estimated error: " << mr.error() << std::endl; // update b, and solve again x = mr.solve(b);
By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method. Here is a step by step execution example starting with a random guess and printing the evolution of the estimated error: *
x = VectorXd::Random(n); mr.setMaxIterations(1); int i = 0; do { x = mr.solveWithGuess(b,x); std::cout << i << " : " << mr.error() << std::endl; ++i; } while (mr.info()!=Success && i<100);
Note that such a step by step excution is slightly slower.
typedef IterativeSolverBase<MINRES> Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::Base [private] |
typedef MatrixType::Index Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::Index |
Reimplemented from Eigen::IterativeSolverBase< MINRES< _MatrixType, _UpLo, _Preconditioner > >.
typedef _MatrixType Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::MatrixType |
Reimplemented from Eigen::IterativeSolverBase< MINRES< _MatrixType, _UpLo, _Preconditioner > >.
typedef _Preconditioner Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::Preconditioner |
Reimplemented from Eigen::IterativeSolverBase< MINRES< _MatrixType, _UpLo, _Preconditioner > >.
typedef MatrixType::RealScalar Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::RealScalar |
Reimplemented from Eigen::IterativeSolverBase< MINRES< _MatrixType, _UpLo, _Preconditioner > >.
typedef MatrixType::Scalar Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::Scalar |
Reimplemented from Eigen::IterativeSolverBase< MINRES< _MatrixType, _UpLo, _Preconditioner > >.
anonymous enum |
Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::MINRES | ( | ) | [inline] |
Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::MINRES | ( | const MatrixType & | A | ) | [inline] |
Initialize the solver with matrix A for further Ax=b
solving.
This constructor is a shortcut for the default constructor followed by a call to compute().
Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::~MINRES | ( | ) | [inline] |
void Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::_solve | ( | const Rhs & | b, |
Dest & | x | ||
) | const [inline] |
void Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::_solveWithGuess | ( | const Rhs & | b, |
Dest & | x | ||
) | const [inline] |
const internal::solve_retval_with_guess<MINRES, Rhs, Guess> Eigen::MINRES< _MatrixType, _UpLo, _Preconditioner >::solveWithGuess | ( | const MatrixBase< Rhs > & | b, |
const Guess & | x0 | ||
) | const [inline] |