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
-
_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. |
_Preconditioner | the 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;
SparseMatrix<double> A(n,n);
ConjugateGradient<SparseMatrix<double> > cg;
x = cg.solve(b);
std::cout << "#iterations: " << cg.iterations() << std::endl;
std::cout << "estimated error: " << cg.error() << std::endl;
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 SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
|
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | analyzePattern (const EigenBase< InputDerived > &A) |
|
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | compute (const EigenBase< InputDerived > &A) |
|
| ConjugateGradient () |
|
template<typename MatrixDerived > |
| ConjugateGradient (const EigenBase< MatrixDerived > &A) |
|
RealScalar | error () const |
|
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | factorize (const EigenBase< InputDerived > &A) |
|
ComputationInfo | info () const |
|
int | iterations () const |
|
int | maxIterations () const |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () const |
|
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | setMaxIterations (int maxIters) |
|
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | setTolerance (const 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 |
|
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 EigenBase< Derived >::derived(), and Eigen::Success.