Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > Class Template Reference

A conjugate gradient solver for sparse self-adjoint problems. More...

#include <ConjugateGradient.h>

+ Inheritance diagram for ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >:

Public Types

enum  { 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
void _solve_sparse (const Rhs &b, SparseMatrix< DestScalar, DestOptions, DestIndex > &dest) const
template<typename Rhs , typename Dest >
void _solveWithGuess (const Rhs &b, Dest &x) const
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
analyzePattern (const MatrixType &A)
Index cols () const
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
compute (const MatrixType &A)
 ConjugateGradient ()
 ConjugateGradient (const MatrixType &A)
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
derived ()
const ConjugateGradient
< _MatrixType, _UpLo,
_Preconditioner > & 
derived () const
RealScalar error () const
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
factorize (const MatrixType &A)
ComputationInfo info () const
int iterations () const
int maxIterations () const
Preconditionerpreconditioner ()
const Preconditionerpreconditioner () const
Index rows () const
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
setMaxIterations (int maxIters)
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
setTolerance (RealScalar tolerance)
const internal::solve_retval
< ConjugateGradient
< _MatrixType, _UpLo,
_Preconditioner >, Rhs > 
solve (const MatrixBase< Rhs > &b) const
const
internal::sparse_solve_retval
< IterativeSolverBase, Rhs > 
solve (const SparseMatrixBase< Rhs > &b) const
template<typename Rhs , typename Guess >
const
internal::solve_retval_with_guess
< ConjugateGradient, Rhs,
Guess > 
solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
RealScalar tolerance () const
 ~ConjugateGradient ()

Protected Member Functions

void init ()

Protected Attributes

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, int _UpLo, typename _Preconditioner>
class Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >

A conjugate gradient solver for sparse self-adjoint problems.

This class allows to solve for A.x = b sparse linear problems using a conjugate gradient algorithm. The sparse matrix A must be selfadjoint. 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.
_UpLothe triangular part that will be used for the computations. It can be Lower or Upper. Default is Lower.
_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
ConjugateGradient<SparseMatrix<double> > 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. 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

Member Typedef Documentation

typedef MatrixType::Index Index
typedef _MatrixType MatrixType
typedef _Preconditioner Preconditioner
typedef MatrixType::RealScalar RealScalar
typedef MatrixType::Scalar Scalar

Member Enumeration Documentation

anonymous enum
Enumerator:
UpLo 

Constructor & Destructor Documentation

ConjugateGradient ( )
inline
ConjugateGradient ( 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.
~ConjugateGradient ( )
inline

Member Function Documentation

void _solve ( const Rhs &  b,
Dest &  x 
) const
inline
void _solve_sparse ( const Rhs &  b,
SparseMatrix< DestScalar, DestOptions, DestIndex > &  dest 
) const
inlineinherited
void _solveWithGuess ( const Rhs &  b,
Dest &  x 
) const
inline
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & analyzePattern ( const MatrixType A)
inlineinherited

Initializes the iterative solver for the sparcity pattern of the matrix A for further solving Ax=b problems.

Currently, this function mostly call analyzePattern on the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products.

References IterativeSolverBase< Derived >::derived(), IterativeSolverBase< Derived >::m_analysisIsOk, IterativeSolverBase< Derived >::m_info, IterativeSolverBase< Derived >::m_isInitialized, IterativeSolverBase< Derived >::m_preconditioner, and Eigen::Success.

Index cols ( void  ) const
inlineinherited
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & compute ( const MatrixType A)
inlineinherited

Initializes the iterative solver with the matrix A for further solving Ax=b problems.

Currently, this function mostly initialized/compute the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products.

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.

References IterativeSolverBase< Derived >::derived(), IterativeSolverBase< Derived >::m_analysisIsOk, IterativeSolverBase< Derived >::m_factorizationIsOk, IterativeSolverBase< Derived >::m_info, IterativeSolverBase< Derived >::m_isInitialized, IterativeSolverBase< Derived >::m_preconditioner, IterativeSolverBase< Derived >::mp_matrix, and Eigen::Success.

ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & derived ( )
inlineinherited
const ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & derived ( ) const
inlineinherited
RealScalar error ( ) const
inlineinherited
Returns
the tolerance error reached during the last solve

References eigen_assert, IterativeSolverBase< Derived >::m_error, and IterativeSolverBase< Derived >::m_isInitialized.

ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & factorize ( const MatrixType A)
inlineinherited

Initializes the iterative solver with the numerical values of the matrix A for further solving Ax=b problems.

Currently, this function mostly call factorize on the preconditioner.

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.

References IterativeSolverBase< Derived >::derived(), eigen_assert, IterativeSolverBase< Derived >::m_analysisIsOk, IterativeSolverBase< Derived >::m_factorizationIsOk, IterativeSolverBase< Derived >::m_info, IterativeSolverBase< Derived >::m_preconditioner, IterativeSolverBase< Derived >::mp_matrix, and Eigen::Success.

ComputationInfo info ( ) const
inlineinherited
Returns
Success if the iterations converged, and NoConvergence otherwise.

References eigen_assert, IterativeSolverBase< Derived >::m_info, and IterativeSolverBase< Derived >::m_isInitialized.

void init ( )
inlineprotectedinherited
int iterations ( ) const
inlineinherited
Returns
the number of iterations performed during the last solve

References eigen_assert, IterativeSolverBase< Derived >::m_isInitialized, and IterativeSolverBase< Derived >::m_iterations.

int maxIterations ( ) const
inlineinherited
Preconditioner& preconditioner ( )
inlineinherited
Returns
a read-write reference to the preconditioner for custom configuration.

References IterativeSolverBase< Derived >::m_preconditioner.

const Preconditioner& preconditioner ( ) const
inlineinherited
Returns
a read-only reference to the preconditioner.

References IterativeSolverBase< Derived >::m_preconditioner.

Index rows ( void  ) const
inlineinherited
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setMaxIterations ( int  maxIters)
inlineinherited
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setTolerance ( RealScalar  tolerance)
inlineinherited
const internal::solve_retval<ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > , Rhs> solve ( const MatrixBase< Rhs > &  b) const
inlineinherited
Returns
the solution x of $ A x = b $ using the current decomposition of A.
See Also
compute()

References IterativeSolverBase< Derived >::derived(), eigen_assert, IterativeSolverBase< Derived >::m_isInitialized, and IterativeSolverBase< Derived >::rows().

const internal::sparse_solve_retval<IterativeSolverBase, Rhs> solve ( const SparseMatrixBase< Rhs > &  b) const
inlineinherited
const internal::solve_retval_with_guess<ConjugateGradient, Rhs, Guess> solveWithGuess ( const MatrixBase< Rhs > &  b,
const Guess &  x0 
) const
inline
RealScalar tolerance ( ) const
inlineinherited
Returns
the tolerance threshold used by the stopping criteria

References IterativeSolverBase< Derived >::m_tolerance.

Member Data Documentation

bool m_analysisIsOk
mutableprotectedinherited
RealScalar m_error
mutableprotectedinherited
bool m_factorizationIsOk
mutableprotectedinherited
ComputationInfo m_info
mutableprotectedinherited
bool m_isInitialized
mutableprotectedinherited
int m_iterations
mutableprotectedinherited
int m_maxIterations
protectedinherited
Preconditioner m_preconditioner
protectedinherited
RealScalar m_tolerance
protectedinherited
const MatrixType* mp_matrix
protectedinherited

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