A conjugate gradient solver for sparse self-adjoint problems. More...
#include <ConjugateGradient.h>
Public Member Functions | |
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 |
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 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); SparseMatrix<double> A(n,n); // fill A and b ConjugateGradient<SparseMatrix<double> > 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. Here is a step by step execution example starting with a random guess and printing the evolution of the estimated error: *
x = VectorXd::Random(n); cg.setMaxIterations(1); int i = 0; do { x = cg.solveWithGuess(b,x); std::cout << i << " : " << cg.error() << std::endl; ++i; } while (cg.info()!=Success && i<100);
Note that such a step by step excution is slightly slower.
ConjugateGradient | ( | ) | [inline] |
Default constructor.
Referenced by ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::solveWithGuess().
ConjugateGradient | ( | const MatrixType & | A | ) | [inline] |
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().
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & analyzePattern | ( | const MatrixType & | A | ) | [inline, inherited] |
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 Eigen::Success.
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & compute | ( | const MatrixType & | A | ) | [inline, inherited] |
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 Eigen::Success.
RealScalar error | ( | ) | const [inline, inherited] |
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & factorize | ( | const MatrixType & | A | ) | [inline, inherited] |
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 Eigen::Success.
ComputationInfo info | ( | ) | const [inline, inherited] |
int iterations | ( | ) | const [inline, inherited] |
int maxIterations | ( | ) | const [inline, inherited] |
Preconditioner& preconditioner | ( | ) | [inline, inherited] |
const Preconditioner& preconditioner | ( | ) | const [inline, inherited] |
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setMaxIterations | ( | int | maxIters | ) | [inline, inherited] |
Sets the max number of iterations
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setTolerance | ( | RealScalar | tolerance | ) | [inline, inherited] |
Sets the tolerance threshold used by the stopping criteria
References IterativeSolverBase< Derived >::tolerance().
const internal::solve_retval<ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > , Rhs> solve | ( | const MatrixBase< Rhs > & | b | ) | const [inline, inherited] |
const internal::sparse_solve_retval<IterativeSolverBase, Rhs> solve | ( | const SparseMatrixBase< Rhs > & | b | ) | const [inline, inherited] |
References EigenBase< Derived >::derived(), and SparseMatrixBase< Derived >::rows().
const internal::solve_retval_with_guess<ConjugateGradient, Rhs, Guess> solveWithGuess | ( | const MatrixBase< Rhs > & | b, |
const Guess & | x0 | ||
) | const [inline] |
References ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::ConjugateGradient().
RealScalar tolerance | ( | ) | const [inline, inherited] |