MOAB  4.9.3pre
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 >:
Collaboration diagram for Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >:

List of all members.

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

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 149 of file LeastSquareConjugateGradient.h.


Member Typedef Documentation

template<typename _MatrixType , typename _Preconditioner >
typedef IterativeSolverBase<LeastSquaresConjugateGradient> Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Base [private]
template<typename _MatrixType , typename _Preconditioner >
typedef _MatrixType Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::MatrixType
template<typename _MatrixType , typename _Preconditioner >
typedef _Preconditioner Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Preconditioner
template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::RealScalar Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::RealScalar
template<typename _MatrixType , typename _Preconditioner >
typedef MatrixType::Scalar Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >::Scalar

Constructor & Destructor Documentation

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

Default constructor.

Definition at line 166 of file LeastSquareConjugateGradient.h.

: Base() {}
template<typename _MatrixType , typename _Preconditioner >
template<typename MatrixDerived >
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().

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.

: Base(A.derived()) {}
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.

  {
    x.setZero();
    _solve_with_guess_impl(b.derived(),x);
  }
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]

The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines