The IDR(s)STAB(l) is a combination of IDR(s) and BiCGSTAB(l). It is a short-recurrences Krylov method for sparse square problems. It can outperform both IDR(s) and BiCGSTAB(l). IDR(s)STAB(l) generally closely follows the optimal GMRES convergence in terms of the number of Matrix-Vector products. However, without the increasing cost per iteration of GMRES. IDR(s)STAB(l) is suitable for both indefinite systems and systems with complex eigenvalues.
More...
|
| IDRSTABL () |
|
template<typename MatrixDerived > |
| IDRSTABL (const EigenBase< MatrixDerived > &A) |
|
template<typename Rhs , typename Dest > |
void | _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const |
|
void | setL (Index L) |
|
void | setS (Index S) |
|
| IterativeSolverBase () |
|
| IterativeSolverBase (const EigenBase< MatrixDerived > &A) |
|
| IterativeSolverBase (IterativeSolverBase &&)=default |
|
| ~IterativeSolverBase () |
|
IDRSTABL< MatrixType_, Preconditioner_ > & | analyzePattern (const EigenBase< MatrixDerived > &A) |
|
IDRSTABL< MatrixType_, Preconditioner_ > & | factorize (const EigenBase< MatrixDerived > &A) |
|
IDRSTABL< MatrixType_, Preconditioner_ > & | compute (const EigenBase< MatrixDerived > &A) |
|
EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
|
EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
|
RealScalar | tolerance () const |
|
IDRSTABL< MatrixType_, Preconditioner_ > & | setTolerance (const RealScalar &tolerance) |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () const |
|
Index | maxIterations () const |
|
IDRSTABL< MatrixType_, Preconditioner_ > & | setMaxIterations (Index maxIters) |
|
Index | iterations () const |
|
RealScalar | error () const |
|
const SolveWithGuess< IDRSTABL< MatrixType_, Preconditioner_ >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
|
ComputationInfo | info () const |
|
void | _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const |
|
std::enable_if_t< Rhs::ColsAtCompileTime !=1 &&DestDerived::ColsAtCompileTime !=1 > | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const |
|
std::enable_if_t< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 > | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const |
|
void | _solve_impl (const Rhs &b, Dest &x) const |
|
IDRSTABL< MatrixType_, Preconditioner_ > & | derived () |
|
const IDRSTABL< MatrixType_, Preconditioner_ > & | derived () const |
|
| SparseSolverBase () |
|
| SparseSolverBase (SparseSolverBase &&other) |
|
| ~SparseSolverBase () |
|
Derived & | derived () |
|
const Derived & | derived () const |
|
template<typename Rhs > |
const Solve< Derived, Rhs > | solve (const MatrixBase< Rhs > &b) const |
|
template<typename Rhs > |
const Solve< Derived, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
|
template<typename Rhs , typename Dest > |
void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
|
template<typename MatrixType_, typename Preconditioner_>
class Eigen::IDRSTABL< MatrixType_, Preconditioner_ >
The IDR(s)STAB(l) is a combination of IDR(s) and BiCGSTAB(l). It is a short-recurrences Krylov method for sparse square problems. It can outperform both IDR(s) and BiCGSTAB(l). IDR(s)STAB(l) generally closely follows the optimal GMRES convergence in terms of the number of Matrix-Vector products. However, without the increasing cost per iteration of GMRES. IDR(s)STAB(l) is suitable for both indefinite systems and systems with complex eigenvalues.
This class allows solving for A.x = b sparse linear problems. 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 |
\implsparsesolverconcept
The maximum 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 maximum number of iterations and NumTraits<Scalar>::epsilon() for the tolerance.
The tolerance is the maximum relative residual error: |Ax-b|/|b| for which the linear system is considered solved.
Performance: When using sparse matrices, best performance is achieved 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.
By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method.
IDR(s)STAB(l) can also be used in a matrix-free context, see the following example .
- See also
- class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner
template<typename MatrixType_ , typename Preconditioner_ >
template<typename Rhs , typename Dest >
void Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::_solve_vector_with_guess_impl |
( |
const Rhs & |
b, |
|
|
Dest & |
x |
|
) |
| const |
|
inline |
Loops over the number of columns of b and does the following:
- sets the tolerance and maxIterations
- Calls the function that has the core solver routine
Scalar * b
Definition: benchVecAdd.cpp:17
ComputationInfo m_info
Definition: IterativeSolverBase.h:389
RealScalar m_error
Definition: IterativeSolverBase.h:387
Index m_iterations
Definition: IterativeSolverBase.h:388
const ActualMatrixType & matrix() const
Definition: IterativeSolverBase.h:374
Index maxIterations() const
Definition: IterativeSolverBase.h:251
Preconditioner m_preconditioner
Definition: IterativeSolverBase.h:382
RealScalar m_tolerance
Definition: IterativeSolverBase.h:385
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
@ NoConvergence
Definition: Constants.h:444
Eigen::DenseIndex ret
Definition: level1_cplx_impl.h:43
bool idrstabl(const MatrixType &mat, const Rhs &rhs, Dest &x, const Preconditioner &precond, Index &iters, typename Dest::RealScalar &tol_error, Index L, Index S)
Definition: IDRSTABL.h:46
list x
Definition: plotDoE.py:28
References b, Eigen::internal::idrstabl(), Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::m_error, Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::m_info, Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::m_iterations, Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::m_L, Eigen::IterativeSolverBase< Derived >::m_preconditioner, Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::m_S, Eigen::IterativeSolverBase< Derived >::m_tolerance, Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::matrix(), Eigen::IterativeSolverBase< Derived >::maxIterations(), Eigen::NoConvergence, Eigen::NumericalIssue, ret, Eigen::Success, and plotDoE::x.