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 MatrixType
 
typedef _Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
- Public Types inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
enum  
 
typedef internal::traits< GMRES< _MatrixType, _Preconditioner > >::MatrixType MatrixType
 
typedef internal::traits< GMRES< _MatrixType, _Preconditioner > >::Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

template<typename Rhs , typename Dest >
void _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const
 
Index get_restart ()
 
 GMRES ()
 
template<typename MatrixDerived >
 GMRES (const EigenBase< MatrixDerived > &A)
 
void set_restart (const Index restart)
 
 ~GMRES ()
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
void _solve_impl (const Rhs &b, Dest &x) const
 
void _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const
 
internal::enable_if< Rhs::ColsAtCompileTime!=1 &&DestDerived::ColsAtCompileTime!=1 >::type _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const
 
internal::enable_if< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 >::type _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const
 
GMRES< _MatrixType, _Preconditioner > & analyzePattern (const EigenBase< MatrixDerived > &A)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
GMRES< _MatrixType, _Preconditioner > & compute (const EigenBase< MatrixDerived > &A)
 
RealScalar error () const
 
GMRES< _MatrixType, _Preconditioner > & factorize (const EigenBase< MatrixDerived > &A)
 
ComputationInfo info () const
 
Index iterations () const
 
 IterativeSolverBase ()
 
 IterativeSolverBase (const EigenBase< MatrixDerived > &A)
 
Index maxIterations () const
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
GMRES< _MatrixType, _Preconditioner > & setMaxIterations (Index maxIters)
 
GMRES< _MatrixType, _Preconditioner > & setTolerance (const RealScalar &tolerance)
 
const SolveWithGuess< GMRES< _MatrixType, _Preconditioner >, Rhs, Guess > solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
RealScalar tolerance () const
 
 ~IterativeSolverBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< GMRES< _MatrixType, _Preconditioner > >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
GMRES< _MatrixType, _Preconditioner > & derived ()
 
const GMRES< _MatrixType, _Preconditioner > & derived () const
 
const Solve< GMRES< _MatrixType, _Preconditioner >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< GMRES< _MatrixType, _Preconditioner >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 ~SparseSolverBase ()
 

Private Types

typedef IterativeSolverBase< GMRESBase
 

Private Attributes

Index m_restart
 

Additional Inherited Members

- Protected Types inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
typedef MatrixWrapper::ActualMatrixType ActualMatrixType
 
typedef SparseSolverBase< GMRES< _MatrixType, _Preconditioner > > Base
 
typedef internal::generic_matrix_wrapper< MatrixTypeMatrixWrapper
 
- Protected Member Functions inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
void grab (const InputType &A)
 
void init ()
 
const ActualMatrixTypematrix () const
 
- Protected Attributes inherited from Eigen::IterativeSolverBase< GMRES< _MatrixType, _Preconditioner > >
bool m_analysisIsOk
 
RealScalar m_error
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
Index m_iterations
 
MatrixWrapper m_matrixWrapper
 
Index m_maxIterations
 
Preconditioner m_preconditioner
 
RealScalar m_tolerance
 
- Protected Attributes inherited from Eigen::SparseSolverBase< GMRES< _MatrixType, _Preconditioner > >
bool m_isInitialized
 

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);
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.

GMRES can also be used in a matrix-free context, see the following example .

See also
class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner

Definition at line 221 of file GMRES.h.

Member Typedef Documentation

◆ Base

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

Definition at line 271 of file GMRES.h.

◆ MatrixType

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

Definition at line 283 of file GMRES.h.

◆ Preconditioner

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

Definition at line 286 of file GMRES.h.

◆ RealScalar

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

Definition at line 285 of file GMRES.h.

◆ Scalar

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

Definition at line 284 of file GMRES.h.

Constructor & Destructor Documentation

◆ GMRES() [1/2]

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

Default constructor.

Definition at line 291 of file GMRES.h.

◆ GMRES() [2/2]

template<typename _MatrixType, typename _Preconditioner>
template<typename MatrixDerived >
Eigen::GMRES< _MatrixType, _Preconditioner >::GMRES ( const EigenBase< MatrixDerived > &  A)
inlineexplicit

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 304 of file GMRES.h.

◆ ~GMRES()

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

Definition at line 306 of file GMRES.h.

Member Function Documentation

◆ _solve_vector_with_guess_impl()

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

Definition at line 319 of file GMRES.h.

◆ get_restart()

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

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

Definition at line 310 of file GMRES.h.

◆ set_restart()

template<typename _MatrixType, typename _Preconditioner>
void Eigen::GMRES< _MatrixType, _Preconditioner >::set_restart ( const Index  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 315 of file GMRES.h.

Member Data Documentation

◆ m_restart

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

Definition at line 279 of file GMRES.h.


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


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:41:41