![]() |
Eigen
3.2.7
|
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.
_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:
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 .
Public Member Functions | |
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 |
|
inline |
Default constructor.
Referenced by ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::solveWithGuess().
|
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().
|
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 EigenBase< Derived >::derived(), and Eigen::Success.
|
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.
References EigenBase< Derived >::derived(), and Eigen::Success.
|
inlineinherited |
|
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.
References EigenBase< Derived >::derived(), and Eigen::Success.
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
|
inlineinherited |
Sets the max number of iterations
|
inlineinherited |
Sets the tolerance threshold used by the stopping criteria
References IterativeSolverBase< Derived >::tolerance().
|
inlineinherited |
|
inlineinherited |
References EigenBase< Derived >::derived(), and SparseMatrixBase< Derived >::rows().
|
inline |
References ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::ConjugateGradient().
|
inlineinherited |