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

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

#include <GMRES.h>

Inheritance diagram for Eigen::GMRES< _MatrixType, _Preconditioner >:
Inheritance graph
[legend]

Public Types

typedef MatrixType::Index Index
 
typedef _MatrixType MatrixType
 
typedef _Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
- Public Types inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
typedef MatrixType::Index Index
 
typedef internal::traits< GMRES< _MatrixType, _Preconditioner > >::MatrixType MatrixType
 
typedef internal::traits< GMRES< _MatrixType, _Preconditioner > >::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
 
 ~GMRES ()
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
void _solve_sparse (const Rhs &b, SparseMatrix< DestScalar, DestOptions, DestIndex > &dest) const
 
GMRES< _MatrixType, _Preconditioner > & analyzePattern (const MatrixType &A)
 
Index cols () const
 
GMRES< _MatrixType, _Preconditioner > & compute (const MatrixType &A)
 
GMRES< _MatrixType, _Preconditioner > & derived ()
 
const GMRES< _MatrixType, _Preconditioner > & derived () const
 
RealScalar error () const
 
GMRES< _MatrixType, _Preconditioner > & factorize (const MatrixType &A)
 
ComputationInfo info () const
 
int iterations () const
 
 IterativeSolverBase ()
 
 IterativeSolverBase (const MatrixType &A)
 
int maxIterations () const
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
Index rows () const
 
GMRES< _MatrixType, _Preconditioner > & setMaxIterations (int maxIters)
 
GMRES< _MatrixType, _Preconditioner > & setTolerance (const RealScalar &tolerance)
 
const internal::solve_retval< GMRES< _MatrixType, _Preconditioner >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const internal::sparse_solve_retval< IterativeSolverBase, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
RealScalar tolerance () const
 
 ~IterativeSolverBase ()
 

Private Types

typedef IterativeSolverBase< GMRESBase
 

Private Attributes

int m_restart
 

Additional Inherited Members

- Protected Member Functions inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
void init ()
 
- Protected Attributes inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
bool m_analysisIsOk
 
RealScalar m_error
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
bool m_isInitialized
 
int m_iterations
 
int m_maxIterations
 
Preconditioner m_preconditioner
 
RealScalar m_tolerance
 
const MatrixTypemp_matrix
 

Detailed Description

template<typename _MatrixType, typename _Preconditioner>
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);
// 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

Definition at line 204 of file GMRES.h.

Member Typedef Documentation

template<typename _MatrixType , typename _Preconditioner >
typedef IterativeSolverBase<GMRES> Eigen::GMRES< _MatrixType, _Preconditioner >::Base
private

Definition at line 265 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::Index Eigen::GMRES< _MatrixType, _Preconditioner >::Index

Definition at line 278 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef _MatrixType Eigen::GMRES< _MatrixType, _Preconditioner >::MatrixType

Definition at line 276 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef _Preconditioner Eigen::GMRES< _MatrixType, _Preconditioner >::Preconditioner

Definition at line 280 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::RealScalar Eigen::GMRES< _MatrixType, _Preconditioner >::RealScalar

Definition at line 279 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::Scalar Eigen::GMRES< _MatrixType, _Preconditioner >::Scalar

Definition at line 277 of file GMRES.h.

Constructor & Destructor Documentation

template<typename _MatrixType , typename _Preconditioner >
Eigen::GMRES< _MatrixType, _Preconditioner >::GMRES ( )
inline

Default constructor.

Definition at line 285 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
Eigen::GMRES< _MatrixType, _Preconditioner >::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.

Definition at line 297 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
Eigen::GMRES< _MatrixType, _Preconditioner >::~GMRES ( )
inline

Definition at line 299 of file GMRES.h.

Member Function Documentation

template<typename _MatrixType , typename _Preconditioner >
template<typename Rhs , typename Dest >
void Eigen::GMRES< _MatrixType, _Preconditioner >::_solve ( const Rhs &  b,
Dest &  x 
) const
inline

Definition at line 348 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
template<typename Rhs , typename Dest >
void Eigen::GMRES< _MatrixType, _Preconditioner >::_solveWithGuess ( const Rhs &  b,
Dest &  x 
) const
inline

Definition at line 328 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
int Eigen::GMRES< _MatrixType, _Preconditioner >::get_restart ( )
inline

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

Definition at line 303 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
void Eigen::GMRES< _MatrixType, _Preconditioner >::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.

Definition at line 308 of file GMRES.h.

template<typename _MatrixType , typename _Preconditioner >
template<typename Rhs , typename Guess >
const internal::solve_retval_with_guess<GMRES, Rhs, Guess> Eigen::GMRES< _MatrixType, _Preconditioner >::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()

Definition at line 317 of file GMRES.h.

Member Data Documentation

template<typename _MatrixType , typename _Preconditioner >
int Eigen::GMRES< _MatrixType, _Preconditioner >::m_restart
private

Definition at line 273 of file GMRES.h.


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


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Mon Jun 10 2019 12:35:35