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.
_MatrixType | the type of the matrix A, can be a dense or a sparse matrix. |
_Preconditioner | the type of the preconditioner. Default is LeastSquareDiagonalPreconditioner |
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; 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.
LeastSquaresConjugateGradient () | |
template<typename MatrixDerived > | |
LeastSquaresConjugateGradient (const EigenBase< MatrixDerived > &A) | |
Public Member Functions inherited from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > > | |
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 |
IterativeSolverBase () | |
IterativeSolverBase (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 SolveWithGuess< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
RealScalar | tolerance () const |
Public Member Functions inherited from Eigen::SparseSolverBase< Derived > | |
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 |
SparseSolverBase () | |
| inline |
Default constructor.
| 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().
© Eigen.
Licensed under the MPL2 License.
https://eigen.tuxfamily.org/dox/classEigen_1_1LeastSquaresConjugateGradient.html