|
| | ConjugateGradient () |
| |
| template<typename MatrixDerived > |
| | ConjugateGradient (const EigenBase< MatrixDerived > &A) |
| |
| | ~ConjugateGradient () |
| |
| template<typename Rhs , typename Dest > |
| void | _solve_with_guess_impl (const Rhs &b, Dest &x) const |
| |
| template<typename Rhs , typename Dest > |
| void | _solve_impl (const MatrixBase< Rhs > &b, Dest &x) const |
| |
| template<typename Rhs , typename DestDerived > |
| void | _solve_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const |
| |
| 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) |
| |
| Index | rows () const |
| |
| Index | cols () const |
| |
| 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 |
| |
| template<typename Rhs , typename Dest > |
| void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
| |
| ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | derived () |
| |
| const ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & | 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 _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 TopicMultiThreading for details.
This class can be used as the direct solver classes. Here is a typical usage example:
int n = 10000;
VectorXd x(n), b(n);
cg.compute(A);
x = cg.solve(b);
std::cout << "#iterations: " << cg.iterations() << std::endl;
std::cout << "estimated error: " << cg.error() << std::endl;
x = cg.solve(b);
A conjugate gradient solver for sparse (or dense) self-adjoint problems.
Definition ConjugateGradient.h:159
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
@ Lower
Definition Constants.h:204
@ Upper
Definition Constants.h:206
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
template<typename _MatrixType , int _UpLo, typename _Preconditioner >
template<typename Rhs , typename Dest >
199 {
202 enum {
203 TransposeInput = (!MatrixWrapper::MatrixFree)
205 && (!MatrixType::IsRowMajor)
206 && (!NumTraits<Scalar>::IsComplex)
207 };
208 typedef typename internal::conditional<TransposeInput,Transpose<const ActualMatrixType>,
ActualMatrixType const&>::type RowMajorWrapper;
211 RowMajorWrapper,
212 typename MatrixWrapper::template ConstSelfAdjointViewReturnType<UpLo>::Type
213 >::type SelfAdjointWrapper;
216
217 for(
Index j=0; j<
b.cols(); ++j)
218 {
221
222 typename Dest::ColXpr xj(x,j);
223 RowMajorWrapper row_mat(
matrix());
225 }
226
229 }
#define EIGEN_IMPLIES(a, b)
Definition Macros.h:902
#define EIGEN_STATIC_ASSERT(CONDITION, MSG)
Definition StaticAssert.h:124
ComputationInfo m_info
Definition IterativeSolverBase.h:388
RealScalar m_error
Definition IterativeSolverBase.h:386
Index m_iterations
Definition IterativeSolverBase.h:387
const ActualMatrixType & matrix() const
Definition IterativeSolverBase.h:369
internal::generic_matrix_wrapper< MatrixType > MatrixWrapper
Definition IterativeSolverBase.h:366
Index maxIterations() const
Definition IterativeSolverBase.h:281
MatrixWrapper::ActualMatrixType ActualMatrixType
Definition IterativeSolverBase.h:367
Preconditioner m_preconditioner
Definition IterativeSolverBase.h:381
bool m_isInitialized
Definition SparseSolverBase.h:119
RealScalar m_tolerance
Definition IterativeSolverBase.h:384
@ Success
Definition Constants.h:432
@ NoConvergence
Definition Constants.h:436
EIGEN_DONT_INLINE void conjugate_gradient(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error)
Definition ConjugateGradient.h:28
References Eigen::internal::conjugate_gradient(), EIGEN_IMPLIES, EIGEN_STATIC_ASSERT, Eigen::Lower, Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::m_error, Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::m_info, Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::m_isInitialized, Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::m_iterations, Eigen::IterativeSolverBase< Derived >::m_preconditioner, Eigen::IterativeSolverBase< Derived >::m_tolerance, Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::matrix(), Eigen::IterativeSolverBase< Derived >::maxIterations(), Eigen::NoConvergence, Eigen::Success, Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::UpLo, and Eigen::Upper.
Referenced by Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::_solve_impl().