Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
Eigen::BiCGSTAB< MatrixType_, Preconditioner_ > Class Template Reference

Detailed Description

template<typename MatrixType_, typename Preconditioner_>
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
MatrixType_the type of the sparse matrix A, can be a dense or a sparse matrix.
Preconditioner_the 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);
/* ... fill A and b ... */
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
A bi conjugate gradient stabilized solver for sparse square problems.
Definition BiCGSTAB.h:153
RealScalar error() const
Definition IterativeSolverBase.h:270
Derived & compute(const EigenBase< MatrixDerived > &A)
Definition IterativeSolverBase.h:210
Index iterations() const
Definition IterativeSolverBase.h:262
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:186
A versatible sparse matrix representation.
Definition SparseUtil.h:47
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SparseSolverBase.h:84

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

◆ BiCGSTAB() [1/2]

template<typename MatrixType_ , typename Preconditioner_ >
Eigen::BiCGSTAB< MatrixType_, Preconditioner_ >::BiCGSTAB ( )
inline

Default constructor.

◆ BiCGSTAB() [2/2]

template<typename MatrixType_ , typename Preconditioner_ >
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: