A conjugate gradient solver for sparse self-adjoint problems. More...
#include <ConjugateGradient.h>
Public Types | |
enum | { UpLo = _UpLo } |
typedef MatrixType::Index | Index |
typedef _MatrixType | MatrixType |
typedef _Preconditioner | Preconditioner |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType::Scalar | Scalar |
Public Member Functions | |
template<typename Rhs , typename Dest > | |
void | _solve (const Rhs &b, Dest &x) const |
template<typename Rhs , typename Dest > | |
void | _solveWithGuess (const Rhs &b, Dest &x) const |
ConjugateGradient () | |
ConjugateGradient (const MatrixType &A) | |
template<typename Rhs , typename Guess > | |
const internal::solve_retval_with_guess < ConjugateGradient, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
~ConjugateGradient () | |
Private Types | |
typedef IterativeSolverBase < ConjugateGradient > | Base |
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:
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.
Definition at line 145 of file ConjugateGradient.h.
typedef IterativeSolverBase<ConjugateGradient> Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Base [private] |
Definition at line 147 of file ConjugateGradient.h.
typedef MatrixType::Index Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Index |
Reimplemented from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >.
Definition at line 156 of file ConjugateGradient.h.
typedef _MatrixType Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::MatrixType |
Reimplemented from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >.
Definition at line 154 of file ConjugateGradient.h.
typedef _Preconditioner Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Preconditioner |
Reimplemented from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >.
Definition at line 158 of file ConjugateGradient.h.
typedef MatrixType::RealScalar Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::RealScalar |
Reimplemented from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >.
Definition at line 157 of file ConjugateGradient.h.
typedef MatrixType::Scalar Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Scalar |
Reimplemented from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >.
Definition at line 155 of file ConjugateGradient.h.
anonymous enum |
Definition at line 160 of file ConjugateGradient.h.
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::ConjugateGradient | ( | ) | [inline] |
Default constructor.
Definition at line 167 of file ConjugateGradient.h.
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::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().
Definition at line 179 of file ConjugateGradient.h.
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::~ConjugateGradient | ( | ) | [inline] |
Definition at line 181 of file ConjugateGradient.h.
void Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::_solve | ( | const Rhs & | b, |
Dest & | x | ||
) | const [inline] |
Definition at line 225 of file ConjugateGradient.h.
void Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::_solveWithGuess | ( | const Rhs & | b, |
Dest & | x | ||
) | const [inline] |
Definition at line 201 of file ConjugateGradient.h.
const internal::solve_retval_with_guess<ConjugateGradient, Rhs, Guess> Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::solveWithGuess | ( | const MatrixBase< Rhs > & | b, |
const Guess & | x0 | ||
) | const [inline] |
Definition at line 190 of file ConjugateGradient.h.