Eigen::IDRSTABL< MatrixType_, Preconditioner_ > Class Template Reference

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...

+ Inheritance diagram for Eigen::IDRSTABL< MatrixType_, Preconditioner_ >:

Public Types

typedef MatrixType_ MatrixType
 
typedef Preconditioner_ Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
- Public Types inherited from Eigen::IterativeSolverBase< IDRSTABL< MatrixType_, Preconditioner_ > >
typedef internal::traits< Derived >::MatrixType MatrixType
 
typedef internal::traits< Derived >::Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
typedef MatrixType::StorageIndex StorageIndex
 

Public Member Functions

template<typename Rhs , typename Dest >
void _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const
 
 IDRSTABL ()
 
template<typename MatrixDerived >
 IDRSTABL (const EigenBase< MatrixDerived > &A)
 
void setL (Index L)
 
void setS (Index S)
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< IDRSTABL< MatrixType_, Preconditioner_ > >
void _solve_impl (const Rhs &b, Dest &x) 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_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const
 
Derived & analyzePattern (const EigenBase< MatrixDerived > &A)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
Derived & compute (const EigenBase< MatrixDerived > &A)
 
Derived & derived ()
 
const Derived & derived () const
 
RealScalar error () const
 
Derived & factorize (const EigenBase< MatrixDerived > &A)
 
ComputationInfo info () const
 
Index iterations () const
 
 IterativeSolverBase ()
 
 IterativeSolverBase (const EigenBase< MatrixDerived > &A)
 
 IterativeSolverBase (IterativeSolverBase &&)=default
 
Index maxIterations () const
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
Derived & setMaxIterations (Index maxIters)
 
Derived & setTolerance (const RealScalar &tolerance)
 
const SolveWithGuess< Derived, Rhs, Guess > solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
RealScalar tolerance () const
 
 ~IterativeSolverBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< class >
Derived & derived ()
 
const Derived & derived () const
 
const Solve< Derived, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< Derived, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 SparseSolverBase (SparseSolverBase &&other)
 
 ~SparseSolverBase ()
 

Private Types

typedef IterativeSolverBase< IDRSTABLBase
 

Private Member Functions

const ActualMatrixTypematrix () const
 

Private Attributes

RealScalar m_error
 
ComputationInfo m_info
 
bool m_isInitialized
 
Index m_iterations
 
Index m_L
 
Index m_S
 

Additional Inherited Members

- Public Attributes inherited from Eigen::IterativeSolverBase< IDRSTABL< MatrixType_, Preconditioner_ > >
 ColsAtCompileTime
 
 MaxColsAtCompileTime
 
- Protected Types inherited from Eigen::IterativeSolverBase< IDRSTABL< MatrixType_, Preconditioner_ > >
typedef MatrixWrapper::ActualMatrixType ActualMatrixType
 
typedef SparseSolverBase< Derived > Base
 
typedef internal::generic_matrix_wrapper< MatrixTypeMatrixWrapper
 
- Protected Member Functions inherited from Eigen::IterativeSolverBase< IDRSTABL< MatrixType_, Preconditioner_ > >
void grab (const InputType &A)
 
void init ()
 
const ActualMatrixTypematrix () const
 
- Protected Attributes inherited from Eigen::IterativeSolverBase< IDRSTABL< MatrixType_, Preconditioner_ > >
bool m_analysisIsOk
 
RealScalar m_error
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
bool m_isInitialized
 
Index m_iterations
 
MatrixWrapper m_matrixWrapper
 
Index m_maxIterations
 
Preconditioner m_preconditioner
 
RealScalar m_tolerance
 
- Protected Attributes inherited from Eigen::SparseSolverBase< class >
bool m_isInitialized
 

Detailed Description

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

This class follows the sparse solver concept .

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

Definition at line 412 of file IDRSTABL.h.

Member Typedef Documentation

◆ Base

template<typename MatrixType_ , typename Preconditioner_ >
typedef IterativeSolverBase<IDRSTABL> Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::Base
private

Definition at line 413 of file IDRSTABL.h.

◆ MatrixType

template<typename MatrixType_ , typename Preconditioner_ >
typedef MatrixType_ Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::MatrixType

Definition at line 423 of file IDRSTABL.h.

◆ Preconditioner

template<typename MatrixType_ , typename Preconditioner_ >
typedef Preconditioner_ Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::Preconditioner

Definition at line 426 of file IDRSTABL.h.

◆ RealScalar

template<typename MatrixType_ , typename Preconditioner_ >
typedef MatrixType::RealScalar Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::RealScalar

Definition at line 425 of file IDRSTABL.h.

◆ Scalar

template<typename MatrixType_ , typename Preconditioner_ >
typedef MatrixType::Scalar Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::Scalar

Definition at line 424 of file IDRSTABL.h.

Constructor & Destructor Documentation

◆ IDRSTABL() [1/2]

template<typename MatrixType_ , typename Preconditioner_ >
Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::IDRSTABL ( )
inline

Default constructor.

Definition at line 430 of file IDRSTABL.h.

430 : m_L(2), m_S(4) {}

◆ IDRSTABL() [2/2]

template<typename MatrixType_ , typename Preconditioner_ >
template<typename MatrixDerived >
Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::IDRSTABL ( const EigenBase< MatrixDerived > &  A)
inlineexplicit

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.

Definition at line 443 of file IDRSTABL.h.

443 : Base(A.derived()), m_L(2), m_S(4) {}
SparseMatrix< double > A(n, n)
IterativeSolverBase< IDRSTABL > Base
Definition: IDRSTABL.h:413
SparseMatrix< double, Options_, StorageIndex_ > & derived()

Member Function Documentation

◆ _solve_vector_with_guess_impl()

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:

  1. sets the tolerance and maxIterations
  2. Calls the function that has the core solver routine

Definition at line 452 of file IDRSTABL.h.

452  {
456 
458  }
RealScalar m_error
const ActualMatrixType & matrix() const
ComputationInfo m_info
Preconditioner m_preconditioner
NumericalIssue
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

◆ matrix()

template<typename MatrixType_ , typename Preconditioner_ >
const ActualMatrixType & Eigen::IterativeSolverBase< class >::matrix
private

◆ setL()

template<typename MatrixType_ , typename Preconditioner_ >
void Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::setL ( Index  L)
inline

Sets the parameter L, indicating the amount of minimize residual steps are used.

Definition at line 462 of file IDRSTABL.h.

462  {
463  eigen_assert(L >= 1 && "L needs to be positive");
464  m_L = L;
465  }
MatrixXd L
#define eigen_assert(x)

◆ setS()

template<typename MatrixType_ , typename Preconditioner_ >
void Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::setS ( Index  S)
inline

Sets the parameter S, indicating the dimension of the shadow residual space..

Definition at line 468 of file IDRSTABL.h.

468  {
469  eigen_assert(S >= 1 && "S needs to be positive");
470  m_S = S;
471  }

Member Data Documentation

◆ m_error

template<typename MatrixType_ , typename Preconditioner_ >
RealScalar Eigen::IterativeSolverBase< class >::m_error
private

◆ m_info

template<typename MatrixType_ , typename Preconditioner_ >
ComputationInfo Eigen::IterativeSolverBase< class >::m_info
private

◆ m_isInitialized

template<typename MatrixType_ , typename Preconditioner_ >
bool Eigen::IterativeSolverBase< class >::m_isInitialized
private

◆ m_iterations

template<typename MatrixType_ , typename Preconditioner_ >
Index Eigen::IterativeSolverBase< class >::m_iterations
private

◆ m_L

template<typename MatrixType_ , typename Preconditioner_ >
Index Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::m_L
private

Definition at line 419 of file IDRSTABL.h.

◆ m_S

template<typename MatrixType_ , typename Preconditioner_ >
Index Eigen::IDRSTABL< MatrixType_, Preconditioner_ >::m_S
private

Definition at line 420 of file IDRSTABL.h.


The documentation for this class was generated from the following file: