WPILibC++ 2023.4.3
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > Class Template Reference

A conjugate gradient solver for sparse (or dense) self-adjoint problems. More...

#include </home/runner/work/allwpilib/allwpilib/wpimath/src/main/native/thirdparty/eigen/include/Eigen/src/IterativeLinearSolvers/ConjugateGradient.h>

Inheritance diagram for Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >:
Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > > Eigen::SparseSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > > Eigen::internal::noncopyable

Public Types

enum  { UpLo = _UpLo }
 
typedef _MatrixType MatrixType
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::RealScalar RealScalar
 
typedef _Preconditioner Preconditioner
 
- Public Types inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
enum  
 
typedef internal::traits< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::MatrixType MatrixType
 
typedef internal::traits< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >::Preconditioner Preconditioner
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 
typedef MatrixType::RealScalar RealScalar
 

Public Member Functions

 ConjugateGradient ()
 Default constructor. More...
 
template<typename MatrixDerived >
 ConjugateGradient (const EigenBase< MatrixDerived > &A)
 Initialize the solver with matrix A for further Ax=b solving. More...
 
 ~ConjugateGradient ()
 
template<typename Rhs , typename Dest >
void _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
 IterativeSolverBase ()
 Default constructor. More...
 
 IterativeSolverBase (const EigenBase< MatrixDerived > &A)
 Initialize the solver with matrix A for further Ax=b solving. More...
 
 ~IterativeSolverBase ()
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & analyzePattern (const EigenBase< MatrixDerived > &A)
 Initializes the iterative solver for the sparsity pattern of the matrix A for further solving Ax=b problems. More...
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & factorize (const EigenBase< MatrixDerived > &A)
 Initializes the iterative solver with the numerical values of the matrix A for further solving Ax=b problems. More...
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & compute (const EigenBase< MatrixDerived > &A)
 Initializes the iterative solver with the matrix A for further solving Ax=b problems. More...
 
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)
 Sets the tolerance threshold used by the stopping criteria. More...
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
Index maxIterations () const
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setMaxIterations (Index maxIters)
 Sets the max number of iterations. More...
 
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
 
- Public Member Functions inherited from Eigen::SparseSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
 SparseSolverBase ()
 Default constructor. More...
 
 ~SparseSolverBase ()
 
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & derived ()
 
const ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & derived () const
 
const Solve< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 

Additional Inherited Members

- Protected Types inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
typedef SparseSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > > Base
 
typedef internal::generic_matrix_wrapper< MatrixTypeMatrixWrapper
 
typedef MatrixWrapper::ActualMatrixType ActualMatrixType
 
- Protected Member Functions inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
void init ()
 
const ActualMatrixTypematrix () const
 
void grab (const InputType &A)
 
- Protected Attributes inherited from Eigen::IterativeSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
MatrixWrapper m_matrixWrapper
 
Preconditioner m_preconditioner
 
Index m_maxIterations
 
RealScalar m_tolerance
 
RealScalar m_error
 
Index m_iterations
 
ComputationInfo m_info
 
bool m_analysisIsOk
 
bool m_factorizationIsOk
 
bool m_isInitialized
 
- Protected Attributes inherited from Eigen::SparseSolverBase< ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > >
bool m_isInitialized
 

Detailed Description

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
_MatrixTypethe type of the matrix A, can be a dense or a sparse matrix.
_UpLothe 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.
_Preconditionerthe 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);
// fill A and b
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);
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
View matrix as a lower triangular matrix.
Definition: Constants.h:209
@ Upper
View matrix as an upper triangular matrix.
Definition: Constants.h:211
b
Definition: data.h:44

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

Member Typedef Documentation

◆ MatrixType

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
typedef _MatrixType Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::MatrixType

◆ Preconditioner

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
typedef _Preconditioner Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Preconditioner

◆ RealScalar

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
typedef MatrixType::RealScalar Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::RealScalar

◆ Scalar

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
typedef MatrixType::Scalar Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::Scalar

Member Enumeration Documentation

◆ anonymous enum

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
anonymous enum
Enumerator
UpLo 

Constructor & Destructor Documentation

◆ ConjugateGradient() [1/2]

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::ConjugateGradient ( )
inline

Default constructor.

◆ ConjugateGradient() [2/2]

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
template<typename MatrixDerived >
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::ConjugateGradient ( const EigenBase< MatrixDerived > &  A)
inlineexplicit

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().

Warning
this class stores a reference to the matrix A as well as some precomputed values that depend on it. Therefore, if A is changed this class becomes invalid. Call compute() to update it with the new matrix A, or modify a copy of A.

◆ ~ConjugateGradient()

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::~ConjugateGradient ( )
inline

Member Function Documentation

◆ _solve_vector_with_guess_impl()

template<typename _MatrixType , int _UpLo, typename _Preconditioner >
template<typename Rhs , typename Dest >
void Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >::_solve_vector_with_guess_impl ( const Rhs &  b,
Dest &  x 
) const
inline

The documentation for this class was generated from the following file: