$treeview $search $mathjax
Eigen  3.2.5
$projectbrief
$projectbrief
$searchbox

ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > Class Template Reference
[IterativeLinearSolvers module]

A conjugate gradient solver for sparse self-adjoint problems. More...

Inheritance diagram for ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >:

List of all members.

Public Member Functions

ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
analyzePattern (const MatrixType &A)
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
compute (const MatrixType &A)
 ConjugateGradient (const MatrixType &A)
 ConjugateGradient ()
RealScalar error () const
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
factorize (const MatrixType &A)
ComputationInfo info () const
int iterations () const
int maxIterations () const
const Preconditioner & preconditioner () const
Preconditioner & preconditioner ()
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
setMaxIterations (int maxIters)
ConjugateGradient< _MatrixType,
_UpLo, _Preconditioner > & 
setTolerance (const RealScalar &tolerance)
const
internal::sparse_solve_retval
< IterativeSolverBase, Rhs > 
solve (const SparseMatrixBase< Rhs > &b) const
const internal::solve_retval
< ConjugateGradient
< _MatrixType, _UpLo,
_Preconditioner >, Rhs > 
solve (const MatrixBase< Rhs > &b) const
template<typename Rhs , typename Guess >
const
internal::solve_retval_with_guess
< ConjugateGradient, Rhs,
Guess > 
solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
RealScalar tolerance () const

Detailed Description

template<typename _MatrixType, int _UpLo, typename _Preconditioner>
class Eigen::ConjugateGradient< _MatrixType, _UpLo, _Preconditioner >

A conjugate gradient solver for sparse self-adjoint problems.

This class allows to solve for A.x = b sparse linear problems using a conjugate gradient algorithm. The sparse matrix A must be selfadjoint. The vectors x and b can be either dense or sparse.

Template Parameters:
_MatrixType the type of the matrix A, can be a dense or a sparse matrix.
_UpLo the triangular part that will be used for the computations. It can be Lower, Upper, or Lower|Upper in which the full matrix entries will be considered. Default is Lower.
_Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner

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

Constructor & Destructor Documentation

ConjugateGradient (  )  [inline]

Default constructor.

ConjugateGradient ( const MatrixType &  A  )  [inline]

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.

Member Function Documentation

ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & analyzePattern ( const MatrixType &  A  )  [inline, inherited]

Initializes the iterative solver for the sparcity pattern of the matrix A for further solving Ax=b problems.

Currently, this function mostly call analyzePattern on the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products.

References Eigen::Success.

ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & compute ( const MatrixType &  A  )  [inline, inherited]

Initializes the iterative solver with the matrix A for further solving Ax=b problems.

Currently, this function mostly initialized/compute the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products.

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.

References Eigen::Success.

RealScalar error (  )  const [inline, inherited]
Returns:
the tolerance error reached during the last solve
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & factorize ( const MatrixType &  A  )  [inline, inherited]

Initializes the iterative solver with the numerical values of the matrix A for further solving Ax=b problems.

Currently, this function mostly call factorize on the preconditioner.

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.

References Eigen::Success.

ComputationInfo info (  )  const [inline, inherited]
Returns:
Success if the iterations converged, and NoConvergence otherwise.
int iterations (  )  const [inline, inherited]
Returns:
the number of iterations performed during the last solve
int maxIterations (  )  const [inline, inherited]
Returns:
the max number of iterations
const Preconditioner& preconditioner (  )  const [inline, inherited]
Returns:
a read-only reference to the preconditioner.
Preconditioner& preconditioner (  )  [inline, inherited]
Returns:
a read-write reference to the preconditioner for custom configuration.
ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setMaxIterations ( int  maxIters  )  [inline, inherited]

Sets the max number of iterations

ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > & setTolerance ( const RealScalar &  tolerance  )  [inline, inherited]

Sets the tolerance threshold used by the stopping criteria

const internal::sparse_solve_retval<IterativeSolverBase, Rhs> solve ( const SparseMatrixBase< Rhs > &  b  )  const [inline, inherited]
Returns:
the solution x of $ A x = b $ using the current decomposition of A.
See also:
compute()

References EigenBase< Derived >::derived(), and SparseMatrixBase< Derived >::rows().

const internal::solve_retval<ConjugateGradient< _MatrixType, _UpLo, _Preconditioner > , Rhs> solve ( const MatrixBase< Rhs > &  b  )  const [inline, inherited]
Returns:
the solution x of $ A x = b $ using the current decomposition of A.
See also:
compute()
const internal::solve_retval_with_guess<ConjugateGradient, Rhs, Guess> solveWithGuess ( const MatrixBase< Rhs > &  b,
const Guess &  x0 
) const [inline]
Returns:
the solution x of $ A x = b $ using the current decomposition of A x0 as an initial solution.
See also:
compute()
RealScalar tolerance (  )  const [inline, inherited]
Returns:
the tolerance threshold used by the stopping criteria

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