|
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | analyzePattern (const MatrixType &A) |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | compute (const MatrixType &A) |
| |
| | ConjugateGradient () |
| |
| | ConjugateGradient (const MatrixType &A) |
| |
| RealScalar | error () const |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | factorize (const MatrixType &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 (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 |
| |
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 sparse 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 or Upper. 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;
VectorXd x(n), b(n);
std::cout <<
"#iterations: " << cg.
iterations() << std::endl;
std::cout <<
"estimated error: " << cg.
error() << std::endl;
ConjugateGradient()
Definition ConjugateGradient.h:166
Derived & compute(const MatrixType &A)
Definition IterativeSolverBase.h:103
int iterations() const
Definition IterativeSolverBase.h:149
RealScalar error() const
Definition IterativeSolverBase.h:156
const internal::solve_retval< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition IterativeSolverBase.h:167
A versatible sparse matrix representation.
Definition SparseMatrix.h:87
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:
int i = 0;
do {
std::cout << i <<
" : " << cg.
error() << std::endl;
++i;
const internal::solve_retval_with_guess< ConjugateGradient, Rhs, Guess > solveWithGuess(const MatrixBase< Rhs > &b, const Guess &x0) const
Definition ConjugateGradient.h:189
static const CwiseNullaryOp< internal::scalar_random_op< Scalar >, Matrix< double, _Rows, _Cols, _Options, _MaxRows, _MaxCols > > Random()
Definition Random.h:97
Derived & setMaxIterations(int maxIters)
Definition IterativeSolverBase.h:142
ComputationInfo info() const
Definition IterativeSolverBase.h:190
@ Success
Definition Constants.h:369
Note that such a step by step excution is slightly slower.
- See also
- class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner