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

A conjugate gradient solver for sparse (or dense) self-adjoint problems. More...

#include <ConjugateGradient.h>

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

Public Types

enum  { UpLo = _UpLo }
 
typedef _MatrixType MatrixType
 
typedef _Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
- Public Types inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
enum  
 
typedef internal::traits< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::MatrixType MatrixType
 
typedef internal::traits< ConjugateGradient< _MatrixType, _UpLo, _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
 
 ConjugateGradient ()
 
template<typename MatrixDerived >
 ConjugateGradient (const EigenBase< MatrixDerived > &A)
 
 ~ConjugateGradient ()
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _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
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & analyzePattern (const EigenBase< MatrixDerived > &A)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & compute (const EigenBase< MatrixDerived > &A)
 
RealScalar error () const
 
ConjugateGradient< _MatrixType, _UpLo, _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
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setMaxIterations (Index maxIters)
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setTolerance (const RealScalar &tolerance)
 
const SolveWithGuess< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >, Rhs, Guess > solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
RealScalar tolerance () const
 
 ~IterativeSolverBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & derived ()
 
const ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & derived () const
 
const Solve< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 ~SparseSolverBase ()
 

Private Types

typedef IterativeSolverBase< ConjugateGradientBase
 

Additional Inherited Members

- Protected Types inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
typedef MatrixWrapper::ActualMatrixType ActualMatrixType
 
typedef SparseSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > > Base
 
typedef internal::generic_matrix_wrapper< MatrixTypeMatrixWrapper
 
- Protected Member Functions inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
void grab (const InputType &A)
 
void init ()
 
const ActualMatrixTypematrix () const
 
- Protected Attributes inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _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< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
bool m_isInitialized
 

Detailed Description

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
class Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >

A conjugate gradient solver for sparse (or dense) self-adjoint problems.

This class allows to solve for A.x = b linear problems using an iterative conjugate gradient algorithm. The matrix A must be selfadjoint. The matrix A and the vectors x and b can be either dense or sparse.

Template Parameters
_MatrixTypethe type of the matrix A, can be a dense or a sparse matrix.
_UpLothe triangular part that will be used for the computations. It can be Lower, Upper, or Lower|Upper in which the full matrix entries will be considered. Default is Lower, best performance is Lower|Upper.
_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.

The tolerance corresponds to the relative residual error: |Ax-b|/|b|

Performance: Even though the default value of _UpLo is Lower, significantly higher performance is achieved when using a complete matrix and Lower|Upper as the _UpLo template parameter. Moreover, in this case multi-threading can be exploited if the user code is compiled with OpenMP enabled. See Eigen and multi-threading for details.

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
ConjugateGradient<SparseMatrix<double>, Lower|Upper> cg;
cg.compute(A);
x = cg.solve(b);
std::cout << "#iterations: " << cg.iterations() << std::endl;
std::cout << "estimated error: " << cg.error() << std::endl;
// update b, and solve again
x = cg.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.

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

See also
class LeastSquaresConjugateGradient, class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner

Definition at line 97 of file ConjugateGradient.h.

Member Typedef Documentation

◆ Base

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
typedef IterativeSolverBase<ConjugateGradient> Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Base
private

Definition at line 160 of file ConjugateGradient.h.

◆ MatrixType

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
typedef _MatrixType Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::MatrixType

Definition at line 167 of file ConjugateGradient.h.

◆ Preconditioner

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
typedef _Preconditioner Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Preconditioner

Definition at line 170 of file ConjugateGradient.h.

◆ RealScalar

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
typedef MatrixType::RealScalar Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::RealScalar

Definition at line 169 of file ConjugateGradient.h.

◆ Scalar

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
typedef MatrixType::Scalar Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Scalar

Definition at line 168 of file ConjugateGradient.h.

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
anonymous enum
Enumerator
UpLo 

Definition at line 172 of file ConjugateGradient.h.

Constructor & Destructor Documentation

◆ ConjugateGradient() [1/2]

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::ConjugateGradient ( )
inline

Default constructor.

Definition at line 179 of file ConjugateGradient.h.

◆ ConjugateGradient() [2/2]

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
template<typename MatrixDerived >
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::ConjugateGradient ( 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 192 of file ConjugateGradient.h.

◆ ~ConjugateGradient()

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::~ConjugateGradient ( )
inline

Definition at line 194 of file ConjugateGradient.h.

Member Function Documentation

◆ _solve_vector_with_guess_impl()

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

Definition at line 198 of file ConjugateGradient.h.


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


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