Public Types | Public Member Functions | Private Types | List of all members
Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > Class Template Reference

A conjugate gradient solver for sparse (or dense) least-square problems. More...

#include <LeastSquareConjugateGradient.h>

Inheritance diagram for Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >:
Inheritance graph
[legend]

Public Types

typedef _MatrixType MatrixType
 
typedef _Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
- Public Types inherited from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >
enum  
 
typedef internal::traits< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >::MatrixType MatrixType
 
typedef internal::traits< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >::Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, Dest &x) const
 
template<typename Rhs , typename Dest >
void _solve_with_guess_impl (const Rhs &b, Dest &x) const
 
 LeastSquaresConjugateGradient ()
 
template<typename MatrixDerived >
 LeastSquaresConjugateGradient (const EigenBase< MatrixDerived > &A)
 
 ~LeastSquaresConjugateGradient ()
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >
void _solve_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const
 
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & analyzePattern (const EigenBase< MatrixDerived > &A)
 
Index cols () const
 
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & compute (const EigenBase< MatrixDerived > &A)
 
RealScalar error () const
 
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & factorize (const EigenBase< MatrixDerived > &A)
 
ComputationInfo info () const
 
Index iterations () const
 
 IterativeSolverBase ()
 
 IterativeSolverBase (const EigenBase< MatrixDerived > &A)
 
Index maxIterations () const
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
Index rows () const
 
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & setMaxIterations (Index maxIters)
 
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & setTolerance (const RealScalar &tolerance)
 
const SolveWithGuess< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs, Guess > solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
RealScalar tolerance () const
 
 ~IterativeSolverBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & derived ()
 
const LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & derived () const
 
const Solve< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 ~SparseSolverBase ()
 

Private Types

typedef IterativeSolverBase< LeastSquaresConjugateGradientBase
 

Additional Inherited Members

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

Detailed Description

template<typename _MatrixType, typename _Preconditioner>
class Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >

A conjugate gradient solver for sparse (or dense) least-square problems.

This class allows to solve for A x = b linear problems using an iterative conjugate gradient algorithm. The matrix A can be non symmetric and rectangular, but the matrix A' A should be positive-definite to guaranty stability. Otherwise, the SparseLU or SparseQR classes might be preferable. 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.
_Preconditionerthe type of the preconditioner. Default is LeastSquareDiagonalPreconditioner

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 m=1000000, n = 10000;
VectorXd x(n), b(m);
SparseMatrix<double> A(m,n);
// fill A and b
LeastSquaresConjugateGradient<SparseMatrix<double> > lscg;
lscg.compute(A);
x = lscg.solve(b);
std::cout << "#iterations: " << lscg.iterations() << std::endl;
std::cout << "estimated error: " << lscg.error() << std::endl;
// update b, and solve again
x = lscg.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.

See also
class ConjugateGradient, SparseLU, SparseQR

Definition at line 98 of file LeastSquareConjugateGradient.h.

Member Typedef Documentation

template<typename _MatrixType, typename _Preconditioner>
typedef IterativeSolverBase<LeastSquaresConjugateGradient> Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Base
private

Definition at line 151 of file LeastSquareConjugateGradient.h.

template<typename _MatrixType, typename _Preconditioner>
typedef _MatrixType Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::MatrixType

Definition at line 158 of file LeastSquareConjugateGradient.h.

template<typename _MatrixType, typename _Preconditioner>
typedef _Preconditioner Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Preconditioner

Definition at line 161 of file LeastSquareConjugateGradient.h.

template<typename _MatrixType, typename _Preconditioner>
typedef MatrixType::RealScalar Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::RealScalar

Definition at line 160 of file LeastSquareConjugateGradient.h.

template<typename _MatrixType, typename _Preconditioner>
typedef MatrixType::Scalar Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Scalar

Definition at line 159 of file LeastSquareConjugateGradient.h.

Constructor & Destructor Documentation

template<typename _MatrixType, typename _Preconditioner>
Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::LeastSquaresConjugateGradient ( )
inline

Default constructor.

Definition at line 166 of file LeastSquareConjugateGradient.h.

template<typename _MatrixType, typename _Preconditioner>
template<typename MatrixDerived >
Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::LeastSquaresConjugateGradient ( 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.

Definition at line 179 of file LeastSquareConjugateGradient.h.

template<typename _MatrixType, typename _Preconditioner>
Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::~LeastSquaresConjugateGradient ( )
inline

Definition at line 181 of file LeastSquareConjugateGradient.h.

Member Function Documentation

template<typename _MatrixType, typename _Preconditioner>
template<typename Rhs , typename Dest >
void Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::_solve_impl ( const MatrixBase< Rhs > &  b,
Dest &  x 
) const
inline

Definition at line 206 of file LeastSquareConjugateGradient.h.

template<typename _MatrixType, typename _Preconditioner>
template<typename Rhs , typename Dest >
void Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::_solve_with_guess_impl ( const Rhs &  b,
Dest &  x 
) const
inline

Definition at line 185 of file LeastSquareConjugateGradient.h.


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


gtsam
Author(s):
autogenerated on Sat May 8 2021 02:52:52