|
template<typename Rhs , typename Dest > |
void | _solve (const Rhs &b, Dest &x) const |
|
template<typename Rhs , typename Dest > |
void | _solveWithGuess (const Rhs &b, Dest &x) const |
|
| BiCGSTAB () |
|
| BiCGSTAB (const MatrixType &A) |
|
template<typename Rhs , typename Guess > |
const internal::solve_retval_with_guess< BiCGSTAB, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
|
| ~BiCGSTAB () |
|
void | _solve_sparse (const Rhs &b, SparseMatrix< DestScalar, DestOptions, DestIndex > &dest) const |
|
BiCGSTAB< _MatrixType, _Preconditioner > & | analyzePattern (const MatrixType &A) |
|
Index | cols () const |
|
BiCGSTAB< _MatrixType, _Preconditioner > & | compute (const MatrixType &A) |
|
BiCGSTAB< _MatrixType, _Preconditioner > & | derived () |
|
const BiCGSTAB< _MatrixType, _Preconditioner > & | derived () const |
|
RealScalar | error () const |
|
BiCGSTAB< _MatrixType, _Preconditioner > & | factorize (const MatrixType &A) |
|
ComputationInfo | info () const |
|
int | iterations () const |
|
| IterativeSolverBase () |
|
| IterativeSolverBase (const MatrixType &A) |
|
int | maxIterations () const |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () const |
|
Index | rows () 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 |
|
RealScalar | tolerance () const |
|
| ~IterativeSolverBase () |
|
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;
SparseMatrix<double> A(n,n);
BiCGSTAB<SparseMatrix<double> > solver;
solver(A);
std::cout << "#iterations: " << solver.iterations() << std::endl;
std::cout << "estimated error: " << solver.error() << std::endl;
By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method. Here is a step by step execution example starting with a random guess and printing the evolution of the estimated error:
solver.setMaxIterations(1);
int i = 0;
do {
x = solver.solveWithGuess(b,x);
std::cout << i << " : " << solver.error() << std::endl;
++i;
}
while (solver.info()!=
Success && i<100);
Note that such a step by step excution is slightly slower.
- See also
- class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
Definition at line 112 of file BiCGSTAB.h.