A bi conjugate gradient stabilized solver for sparse square problems. More...
#include <BiCGSTAB.h>
Public Types | |
typedef MatrixType::Index | Index |
typedef _MatrixType | MatrixType |
typedef _Preconditioner | Preconditioner |
typedef MatrixType::RealScalar | RealScalar |
typedef MatrixType::Scalar | Scalar |
Public Member Functions | |
template<typename Rhs , typename Dest > | |
void | _solve (const Rhs &b, Dest &x) const |
void | _solve_sparse (const Rhs &b, SparseMatrix< DestScalar, DestOptions, DestIndex > &dest) const |
template<typename Rhs , typename Dest > | |
void | _solveWithGuess (const Rhs &b, Dest &x) const |
BiCGSTAB< _MatrixType, _Preconditioner > & | analyzePattern (const MatrixType &A) |
BiCGSTAB () | |
BiCGSTAB (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 |
int | maxIterations () const |
Preconditioner & | preconditioner () |
const Preconditioner & | preconditioner () const |
Index | rows () const |
BiCGSTAB< _MatrixType, _Preconditioner > & | setMaxIterations (int maxIters) |
BiCGSTAB< _MatrixType, _Preconditioner > & | setTolerance (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 () | |
Protected Member Functions | |
void | init () |
Protected Attributes | |
bool | m_analysisIsOk |
RealScalar | m_error |
bool | m_factorizationIsOk |
ComputationInfo | m_info |
bool | m_isInitialized |
int | m_iterations |
int | m_maxIterations |
Preconditioner | m_preconditioner |
RealScalar | m_tolerance |
const MatrixType * | mp_matrix |
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(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. Here is a step by step execution example starting with a random guess and printing the evolution of the estimated error: *
x = VectorXd::Random(n); 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.
typedef MatrixType::Index Index |
Reimplemented from IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >.
typedef _MatrixType MatrixType |
Reimplemented from IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >.
typedef _Preconditioner Preconditioner |
Reimplemented from IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >.
typedef MatrixType::RealScalar RealScalar |
Reimplemented from IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >.
typedef MatrixType::Scalar Scalar |
Reimplemented from IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >.
BiCGSTAB | ( | ) | [inline] |
Default constructor.
Referenced by BiCGSTAB< _MatrixType, _Preconditioner >::solveWithGuess().
BiCGSTAB | ( | const MatrixType & | A | ) | [inline] |
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 | ( | ) | [inline] |
void _solve | ( | const Rhs & | b, |
Dest & | x | ||
) | const [inline] |
void _solve_sparse | ( | const Rhs & | b, |
SparseMatrix< DestScalar, DestOptions, DestIndex > & | dest | ||
) | const [inline, inherited] |
void _solveWithGuess | ( | const Rhs & | b, |
Dest & | x | ||
) | const [inline] |
References Eigen::internal::bicgstab(), IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >::m_error, IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >::m_info, IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >::m_isInitialized, IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >::m_iterations, IterativeSolverBase< Derived >::m_preconditioner, IterativeSolverBase< Derived >::m_tolerance, IterativeSolverBase< Derived >::maxIterations(), IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >::mp_matrix, Eigen::NoConvergence, Eigen::NumericalIssue, and Eigen::Success.
Referenced by BiCGSTAB< _MatrixType, _Preconditioner >::_solve().
BiCGSTAB< _MatrixType, _Preconditioner > & analyzePattern | ( | const MatrixType & | 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 MatrixType & | 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 MatrixType & | 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] |
void init | ( | ) | [inline, protected, 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 | ( | 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] |
const internal::sparse_solve_retval<IterativeSolverBase, Rhs> solve | ( | const SparseMatrixBase< Rhs > & | b | ) | const [inline, inherited] |
const internal::solve_retval_with_guess<BiCGSTAB, Rhs, Guess> solveWithGuess | ( | const MatrixBase< Rhs > & | b, |
const Guess & | x0 | ||
) | const [inline] |
References BiCGSTAB< _MatrixType, _Preconditioner >::BiCGSTAB(), IterativeSolverBase< BiCGSTAB< _MatrixType, _Preconditioner > >::m_isInitialized, and IterativeSolverBase< Derived >::rows().
RealScalar tolerance | ( | ) | const [inline, inherited] |
bool m_analysisIsOk [mutable, protected, inherited] |
RealScalar m_error [mutable, protected, inherited] |
Referenced by BiCGSTAB< _MatrixType, _Preconditioner >::_solveWithGuess().
bool m_factorizationIsOk [mutable, protected, inherited] |
ComputationInfo m_info [mutable, protected, inherited] |
Referenced by BiCGSTAB< _MatrixType, _Preconditioner >::_solveWithGuess().
bool m_isInitialized [mutable, protected, inherited] |
int m_iterations [mutable, protected, inherited] |
Referenced by BiCGSTAB< _MatrixType, _Preconditioner >::_solveWithGuess().
int m_maxIterations [protected, inherited] |
Preconditioner m_preconditioner [protected, inherited] |
RealScalar m_tolerance [protected, inherited] |
const MatrixType* mp_matrix [protected, inherited] |
Referenced by BiCGSTAB< _MatrixType, _Preconditioner >::_solveWithGuess().