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 |
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);
std::cout <<
"#iterations: " << solver.
iterations() << std::endl;
std::cout <<
"estimated error: " << solver.
error() << std::endl;
BiCGSTAB()
Definition BiCGSTAB.h:177
int iterations() const
Definition IterativeSolverBase.h:154
RealScalar error() const
Definition IterativeSolverBase.h:161
const internal::solve_retval< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition IterativeSolverBase.h:172
Derived & compute(const EigenBase< InputDerived > &A)
Definition IterativeSolverBase.h:108
A versatible sparse matrix representation.
Definition SparseMatrix.h:87
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
|
| BiCGSTAB< _MatrixType, _Preconditioner > & | analyzePattern (const EigenBase< InputDerived > &A) |
| |
| | BiCGSTAB () |
| |
| template<typename MatrixDerived> |
| | BiCGSTAB (const EigenBase< MatrixDerived > &A) |
| |
| BiCGSTAB< _MatrixType, _Preconditioner > & | compute (const EigenBase< InputDerived > &A) |
| |
| RealScalar | error () const |
| |
| BiCGSTAB< _MatrixType, _Preconditioner > & | factorize (const EigenBase< InputDerived > &A) |
| |
| ComputationInfo | info () const |
| |
| int | iterations () const |
| |
| int | maxIterations () const |
| |
| Preconditioner & | preconditioner () |
| |
| const Preconditioner & | preconditioner () const |
| |
| BiCGSTAB< _MatrixType, _Preconditioner > & | setMaxIterations (int maxIters) |
| |
| BiCGSTAB< _MatrixType, _Preconditioner > & | setTolerance (const RealScalar &tolerance) |
| |
| const internal::solve_retval< BiCGSTAB< _MatrixType, _Preconditioner >, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| |
| const internal::sparse_solve_retval< IterativeSolverBase, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
| |
| template<typename Rhs, typename Guess> |
| const internal::solve_retval_with_guess< BiCGSTAB, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
| |
| RealScalar | tolerance () const |
| |
template<typename _MatrixType, typename _Preconditioner>
template<typename MatrixDerived>
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.