Eigen  3.3.90 (mercurial changeset 700bc7d4d76f)
 All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Eigen::BiCGSTAB< _MatrixType, _Preconditioner > Class Template Reference

Detailed Description

template<typename _MatrixType, typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
class Eigen::BiCGSTAB< _MatrixType, _Preconditioner >

A bi conjugate gradient stabilized solver for sparse square problems.

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

Template Parameters
_MatrixTypethe type of the sparse matrix A, can be a dense or a sparse matrix.
_Preconditionerthe type of the preconditioner. Default is DiagonalPreconditioner

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.

The tolerance corresponds to the relative residual error: |Ax-b|/|b|

Performance: when using sparse matrices, best performance is achied for a row-major sparse matrix format. Moreover, in this case multi-threading can be exploited if the user code is compiled with OpenMP enabled. See Eigen and multi-threading for details.

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 ... */
BiCGSTAB<SparseMatrix<double> > solver;
solver.compute(A);
x = solver.solve(b);
std::cout << "#iterations: " << solver.iterations() << std::endl;
std::cout << "estimated error: " << solver.error() << std::endl;
/* ... update b ... */
x = solver.solve(b); // solve again

By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method.

BiCGSTAB can also be used in a matrix-free context, see the following example .

See Also
class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
+ Inheritance diagram for Eigen::BiCGSTAB< _MatrixType, _Preconditioner >:

Public Member Functions

 BiCGSTAB ()
 
template<typename MatrixDerived >
 BiCGSTAB (const EigenBase< MatrixDerived > &A)
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >
BiCGSTAB< _MatrixType,
_Preconditioner > & 
analyzePattern (const EigenBase< MatrixDerived > &A)
 
BiCGSTAB< _MatrixType,
_Preconditioner > & 
compute (const EigenBase< MatrixDerived > &A)
 
RealScalar error () const
 
BiCGSTAB< _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
 
BiCGSTAB< _MatrixType,
_Preconditioner > & 
setMaxIterations (Index maxIters)
 
BiCGSTAB< _MatrixType,
_Preconditioner > & 
setTolerance (const RealScalar &tolerance)
 
const SolveWithGuess< BiCGSTAB
< _MatrixType, _Preconditioner >
, Rhs, Guess > 
solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
RealScalar tolerance () const
 
- Public Member Functions inherited from Eigen::SparseSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >
const Solve< BiCGSTAB
< _MatrixType, _Preconditioner >
, Rhs > 
solve (const MatrixBase< Rhs > &b) const
 
const Solve< BiCGSTAB
< _MatrixType, _Preconditioner >
, Rhs > 
solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 

Constructor & Destructor Documentation

template<typename _MatrixType, typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
Eigen::BiCGSTAB< _MatrixType, _Preconditioner >::BiCGSTAB ( )
inline

Default constructor.

template<typename _MatrixType, typename _Preconditioner = DiagonalPreconditioner<typename _MatrixType::Scalar>>
template<typename MatrixDerived >
Eigen::BiCGSTAB< _MatrixType, _Preconditioner >::BiCGSTAB ( const EigenBase< MatrixDerived > &  A)
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().

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.

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