template<typename _MatrixType, int _UpLo, typename _Preconditioner>
ConjugateGradient class
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
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 |
Contents
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.
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);
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.
Base classes
-
template<typename Derived>class IterativeSolverBase
- Base class for linear iterative solvers.
Constructors, destructors, conversion operators
- ConjugateGradient()
-
template<typename MatrixDerived>ConjugateGradient(const EigenBase<MatrixDerived>& A) explicit
Function documentation
template<typename _MatrixType, int _UpLo, typename _Preconditioner>
Eigen:: ConjugateGradient<_MatrixType, _UpLo, _Preconditioner>:: ConjugateGradient()
Default constructor.
template<typename _MatrixType, int _UpLo, typename _Preconditioner>
template<typename MatrixDerived>
Eigen:: ConjugateGradient<_MatrixType, _UpLo, _Preconditioner>:: ConjugateGradient(const EigenBase<MatrixDerived>& A) explicit
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().