MOAB
4.9.3pre
|
A conjugate gradient solver for sparse (or dense) least-square problems. More...
#include <LeastSquareConjugateGradient.h>
Public Types | |
typedef _MatrixType | MatrixType |
typedef MatrixType::Scalar | Scalar |
typedef MatrixType::RealScalar | RealScalar |
typedef _Preconditioner | Preconditioner |
Public Member Functions | |
LeastSquaresConjugateGradient () | |
template<typename MatrixDerived > | |
LeastSquaresConjugateGradient (const EigenBase< MatrixDerived > &A) | |
~LeastSquaresConjugateGradient () | |
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 |
Private Types | |
typedef IterativeSolverBase < LeastSquaresConjugateGradient > | Base |
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 |
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.
Definition at line 149 of file LeastSquareConjugateGradient.h.
typedef IterativeSolverBase<LeastSquaresConjugateGradient> Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Base [private] |
Reimplemented from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >.
Definition at line 151 of file LeastSquareConjugateGradient.h.
typedef _MatrixType Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::MatrixType |
Reimplemented from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >.
Definition at line 158 of file LeastSquareConjugateGradient.h.
typedef _Preconditioner Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Preconditioner |
Reimplemented from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >.
Definition at line 161 of file LeastSquareConjugateGradient.h.
typedef MatrixType::RealScalar Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::RealScalar |
Reimplemented from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >.
Definition at line 160 of file LeastSquareConjugateGradient.h.
typedef MatrixType::Scalar Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Scalar |
Reimplemented from Eigen::IterativeSolverBase< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > >.
Definition at line 159 of file LeastSquareConjugateGradient.h.
Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::LeastSquaresConjugateGradient | ( | ) | [inline] |
Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::LeastSquaresConjugateGradient | ( | const EigenBase< MatrixDerived > & | A | ) | [inline, explicit] |
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().
Definition at line 179 of file LeastSquareConjugateGradient.h.
: Base(A.derived()) {}
Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::~LeastSquaresConjugateGradient | ( | ) | [inline] |
Definition at line 181 of file LeastSquareConjugateGradient.h.
{}
void Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::_solve_impl | ( | const MatrixBase< Rhs > & | b, |
Dest & | x | ||
) | const [inline] |
Definition at line 206 of file LeastSquareConjugateGradient.h.
{ x.setZero(); _solve_with_guess_impl(b.derived(),x); }
void Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::_solve_with_guess_impl | ( | const Rhs & | b, |
Dest & | x | ||
) | const [inline] |
Definition at line 185 of file LeastSquareConjugateGradient.h.
{ m_iterations = Base::maxIterations(); m_error = Base::m_tolerance; for(Index j=0; j<b.cols(); ++j) { m_iterations = Base::maxIterations(); m_error = Base::m_tolerance; typename Dest::ColXpr xj(x,j); internal::least_square_conjugate_gradient(matrix(), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error); } m_isInitialized = true; m_info = m_error <= Base::m_tolerance ? Success : NoConvergence; }