Eigen-unsupported  3.4.90 (git rev 67eeba6e720c5745abc77ae6c92ce0a44aa7b7ae)
Eigen::IDRSTABL< MatrixType_, Preconditioner_ > Class Template Reference

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
+ Inheritance diagram for Eigen::IDRSTABL< MatrixType_, Preconditioner_ >:

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)
 

Constructor & Destructor Documentation

◆ IDRSTABL() [1/2]

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

Default constructor.

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

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

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

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


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