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
-
This class follows the sparse solver concept .
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;
SparseMatrix<double> A(m,n);
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;
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
|
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | analyzePattern (const EigenBase< MatrixDerived > &A) |
|
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 |
|
| LeastSquaresConjugateGradient () |
|
template<typename MatrixDerived > |
| LeastSquaresConjugateGradient (const EigenBase< MatrixDerived > &A) |
|
Index | maxIterations () const |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () const |
|
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | setMaxIterations (Index maxIters) |
|
LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | setTolerance (const RealScalar &tolerance) |
|
const Solve< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
|
const Solve< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
|
const SolveWithGuess< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
|
RealScalar | tolerance () const |
|