![]() |
Eigen
3.2.8
|
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.
| _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); 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, and solve again x = solver.solve(b);
By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method.
Inheritance diagram for BiCGSTAB< _MatrixType, _Preconditioner >:Public Member Functions | |
| 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 |
| BiCGSTAB | ( | ) | [inline] |
Default constructor.
Referenced by BiCGSTAB< _MatrixType, _Preconditioner >::solveWithGuess().
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().
| BiCGSTAB< _MatrixType, _Preconditioner > & analyzePattern | ( | const EigenBase< InputDerived > & | A | ) | [inline, inherited] |
Initializes the iterative solver for the sparcity pattern of the matrix A for further solving Ax=b problems.
Currently, this function mostly call analyzePattern on the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products.
| BiCGSTAB< _MatrixType, _Preconditioner > & compute | ( | const EigenBase< InputDerived > & | A | ) | [inline, inherited] |
Initializes the iterative solver with the matrix A for further solving Ax=b problems.
Currently, this function mostly initialized/compute the preconditioner. In the future we might, for instance, implement column reodering for faster matrix vector products.
| RealScalar error | ( | ) | const [inline, inherited] |
| BiCGSTAB< _MatrixType, _Preconditioner > & factorize | ( | const EigenBase< InputDerived > & | A | ) | [inline, inherited] |
Initializes the iterative solver with the numerical values of the matrix A for further solving Ax=b problems.
Currently, this function mostly call factorize on the preconditioner.
| ComputationInfo info | ( | ) | const [inline, inherited] |
| int iterations | ( | ) | const [inline, inherited] |
| int maxIterations | ( | ) | const [inline, inherited] |
| Preconditioner& preconditioner | ( | ) | [inline, inherited] |
| const Preconditioner& preconditioner | ( | ) | const [inline, inherited] |
| BiCGSTAB< _MatrixType, _Preconditioner > & setMaxIterations | ( | int | maxIters | ) | [inline, inherited] |
Sets the max number of iterations
| BiCGSTAB< _MatrixType, _Preconditioner > & setTolerance | ( | const RealScalar & | tolerance | ) | [inline, inherited] |
Sets the tolerance threshold used by the stopping criteria
| const internal::solve_retval<BiCGSTAB< _MatrixType, _Preconditioner > , Rhs> solve | ( | const MatrixBase< Rhs > & | b | ) | const [inline, inherited] |
using the current decomposition of A.| const internal::sparse_solve_retval<IterativeSolverBase, Rhs> solve | ( | const SparseMatrixBase< Rhs > & | b | ) | const [inline, inherited] |
using the current decomposition of A.| const internal::solve_retval_with_guess<BiCGSTAB, Rhs, Guess> solveWithGuess | ( | const MatrixBase< Rhs > & | b, |
| const Guess & | x0 | ||
| ) | const [inline] |
using the current decomposition of A x0 as an initial solution.References BiCGSTAB< _MatrixType, _Preconditioner >::BiCGSTAB().
| RealScalar tolerance | ( | ) | const [inline, inherited] |