Eigen  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::ConjugateGradient< MatrixType_, UpLo_, Preconditioner_ > Class Template Reference

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
MatrixType_the type of the matrix A, can be a dense or a sparse matrix.
UpLo_the 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.
Preconditioner_the type of the preconditioner. Default is DiagonalPreconditioner

This class follows the sparse solver concept .

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);
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
Matrix< double, Dynamic, 1 > VectorXd
DynamicĂ—1 vector of type double.
Definition: Matrix.h:501

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
+ Inheritance diagram for Eigen::ConjugateGradient< MatrixType_, UpLo_, Preconditioner_ >:

Public Member Functions

 ConjugateGradient ()
 
template<typename MatrixDerived >
 ConjugateGradient (const EigenBase< MatrixDerived > &A)
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< ConjugateGradient< MatrixType_, UpLo_, Preconditioner_ > >
ConjugateGradient< MatrixType_, UpLo_, Preconditioner_ > & analyzePattern (const EigenBase< MatrixDerived > &A)
 
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
 
Preconditioner & preconditioner ()
 
const Preconditioner & preconditioner () const
 
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
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
template<typename Rhs >
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
template<typename Rhs >
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 

Constructor & Destructor Documentation

◆ ConjugateGradient() [1/2]

template<typename MatrixType_ , int UpLo_, typename Preconditioner_ >
Eigen::ConjugateGradient< MatrixType_, UpLo_, Preconditioner_ >::ConjugateGradient ( )
inline

Default constructor.

◆ 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.

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