Public Types | Public Member Functions | List of all members
GMRES< _MatrixType, _Preconditioner > Class Template Reference

A GMRES solver for sparse square problems. More...

#include <GMRES.h>

Public Types

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
 
int get_restart ()
 
 GMRES ()
 
 GMRES (const MatrixType &A)
 
void set_restart (const int restart)
 
template<typename Rhs , typename Guess >
const
internal::solve_retval_with_guess
< GMRES, Rhs, Guess > 
solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
- Public Member Functions inherited from IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
GMRES< _MatrixType,
_Preconditioner > & 
analyzePattern (const MatrixType &A)
 
GMRES< _MatrixType,
_Preconditioner > & 
compute (const MatrixType &A)
 
RealScalar error () const
 
GMRES< _MatrixType,
_Preconditioner > & 
factorize (const MatrixType &A)
 
ComputationInfo info () const
 
int iterations () const
 
 IterativeSolverBase (const MatrixType &A)
 
int maxIterations () const
 
const Preconditioner & preconditioner () const
 
Preconditioner & preconditioner ()
 
GMRES< _MatrixType,
_Preconditioner > & 
setMaxIterations (int maxIters)
 
GMRES< _MatrixType,
_Preconditioner > & 
setTolerance (RealScalar tolerance)
 
const
internal::sparse_solve_retval
< IterativeSolverBase, Rhs > 
solve (const SparseMatrixBase< Rhs > &b) const
 
const internal::solve_retval
< GMRES< _MatrixType,
_Preconditioner >, Rhs > 
solve (const MatrixBase< Rhs > &b) const
 
RealScalar tolerance () const
 

Detailed Description

template<typename _MatrixType, typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
class Eigen::GMRES< _MatrixType, _Preconditioner >

A GMRES solver for sparse square problems.

This class allows to solve for A.x = b sparse linear problems using a generalized minimal residual method. The vectors x and b can be either dense or sparse.

Template Parameters
_MatrixTypethe type of the sparse matrix A, can be a dense or a sparse matrix.
_Preconditionerthe 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
GMRES<SparseMatrix<double> > solver(A);
x = solver.solve(b);
std::cout << "#iterations: " << solver.iterations() << std::endl;
std::cout << "estimated error: " << solver.error() << std::endl;
// update b, and solve again
x = solver.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:

See Also
class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner

Constructor & Destructor Documentation

GMRES ( )
inline

Default constructor.

Referenced by GMRES< _MatrixType, _Preconditioner >::solveWithGuess().

GMRES ( 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().

Warning
this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if A is changed this class becomes invalid. Call compute() to update it with the new matrix A, or modify a copy of A.

Member Function Documentation

int get_restart ( )
inline

Get the number of iterations after that a restart is performed.

void set_restart ( const int  restart)
inline

Set the number of iterations after that a restart is performed.

Parameters
restartnumber of iterations for a restarti, default is 30.
const internal::solve_retval_with_guess<GMRES, Rhs, Guess> solveWithGuess ( const MatrixBase< Rhs > &  b,
const Guess &  x0 
) const
inline
Returns
the solution x of $ A x = b $ using the current decomposition of A x0 as an initial solution.
See Also
compute()

References GMRES< _MatrixType, _Preconditioner >::GMRES().


The documentation for this class was generated from the following file: