template<typename _MatrixType, typename _Preconditioner>
BiCGSTAB class
A bi conjugate gradient stabilized solver for sparse square problems.
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 |
Contents
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.
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.
Base classes
-
template<typename Derived>class IterativeSolverBase
- Base class for linear iterative solvers.
Constructors, destructors, conversion operators
Function documentation
template<typename _MatrixType, typename _Preconditioner>
Eigen:: BiCGSTAB<_MatrixType, _Preconditioner>:: BiCGSTAB()
Default constructor.
template<typename _MatrixType, typename _Preconditioner>
template<typename MatrixDerived>
Eigen:: BiCGSTAB<_MatrixType, _Preconditioner>:: BiCGSTAB(const EigenBase<MatrixDerived>& A) 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().