|
| | ConjugateGradient () |
| |
| template<typename MatrixDerived > |
| | ConjugateGradient (const EigenBase< MatrixDerived > &A) |
| |
| | ~ConjugateGradient () |
| |
| template<typename Rhs , typename Dest > |
| void | _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const |
| |
| | IterativeSolverBase () |
| |
| | IterativeSolverBase (const EigenBase< MatrixDerived > &A) |
| |
| | ~IterativeSolverBase () |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | analyzePattern (const EigenBase< MatrixDerived > &A) |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | factorize (const EigenBase< MatrixDerived > &A) |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | compute (const EigenBase< MatrixDerived > &A) |
| |
| EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
| |
| EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
| |
| RealScalar | tolerance () const |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | setTolerance (const RealScalar &tolerance) |
| |
| Preconditioner & | preconditioner () |
| |
| const Preconditioner & | preconditioner () const |
| |
| Index | maxIterations () const |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | setMaxIterations (Index maxIters) |
| |
| Index | iterations () const |
| |
| RealScalar | error () const |
| |
| const SolveWithGuess< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
| |
| ComputationInfo | info () const |
| |
| void | _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const |
| |
| internal::enable_if< Rhs::ColsAtCompileTime!=1 &&DestDerived::ColsAtCompileTime!=1 >::type | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const |
| |
| internal::enable_if< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 >::type | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const |
| |
| void | _solve_impl (const Rhs &b, Dest &x) const |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | derived () |
| |
| const ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | derived () const |
| |
| | SparseSolverBase () |
| |
| | ~SparseSolverBase () |
| |
| Derived & | derived () |
| |
| const Derived & | derived () const |
| |
| template<typename Rhs > |
| const Solve< Derived, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| |
| template<typename Rhs > |
| const Solve< Derived, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
| |
| template<typename Rhs , typename Dest > |
| void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
| |
template<typename _MatrixType,
int _UpLo, typename _Preconditioner>
class Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
This class allows to solve for A.x = b linear problems using an iterative conjugate gradient algorithm. The matrix A must be selfadjoint. The matrix A and 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, best performance is Lower|Upper. |
| _Preconditioner | the type of the preconditioner. Default is DiagonalPreconditioner |
\implsparsesolverconcept
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.
The tolerance corresponds to the relative residual error: |Ax-b|/|b|
Performance: Even though the default value of _UpLo is Lower, significantly higher performance is achieved when using a complete matrix and Lower|Upper as the _UpLo template parameter. Moreover, in this case multi-threading can be exploited if the user code is compiled with OpenMP enabled. See Eigen and multi-threading for details.
This class can be used as the direct solver classes. Here is a typical usage example:
std::cout << "#iterations: " << cg.iterations() << std::endl;
std::cout << "estimated error: " << cg.error() << std::endl;
int n
Definition BiCGSTAB_simple.cpp:1
Scalar * b
Definition benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition bench_gemm.cpp:48
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition ConjugateGradient.h:159
Derived & compute(const EigenBase< MatrixDerived > &A)
Definition IterativeSolverBase.h:238
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition gnuplot_common_settings.hh:12
@ Lower
Definition Constants.h:209
@ Upper
Definition Constants.h:211
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 LeastSquaresConjugateGradient, class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner