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

The Induced Dimension Reduction method (IDR(s)) is a short-recurrences Krylov method for sparse square problems. More...

#include <IDRS.h>

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

Public Types

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

Public Member Functions

 IDRS ()
 
template<typename MatrixDerived >
 IDRS (const EigenBase< MatrixDerived > &A)
 
template<typename Rhs , typename Dest >
void _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const
 
void setS (Index S)
 
void setSmoothing (bool smoothing)
 
void setAngle (RealScalar angle)
 
void setResidualUpdate (bool update)
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< IDRS< MatrixType_, Preconditioner_ > >
 IterativeSolverBase ()
 
 IterativeSolverBase (const EigenBase< MatrixDerived > &A)
 
 IterativeSolverBase (IterativeSolverBase &&)=default
 
 ~IterativeSolverBase ()
 
IDRS< MatrixType_, Preconditioner_ > & analyzePattern (const EigenBase< MatrixDerived > &A)
 
IDRS< MatrixType_, Preconditioner_ > & factorize (const EigenBase< MatrixDerived > &A)
 
IDRS< MatrixType_, Preconditioner_ > & compute (const EigenBase< MatrixDerived > &A)
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
RealScalar tolerance () const
 
IDRS< MatrixType_, Preconditioner_ > & setTolerance (const RealScalar &tolerance)
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
Index maxIterations () const
 
IDRS< MatrixType_, Preconditioner_ > & setMaxIterations (Index maxIters)
 
Index iterations () const
 
RealScalar error () const
 
const SolveWithGuess< IDRS< 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
 
IDRS< MatrixType_, Preconditioner_ > & derived ()
 
const IDRS< MatrixType_, Preconditioner_ > & derived () const
 
- Public Member Functions inherited from Eigen::SparseSolverBase< Derived >
 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
 

Private Types

typedef IterativeSolverBase< IDRSBase
 

Private Member Functions

const ActualMatrixTypematrix () const
 

Private Attributes

Index m_S
 
bool m_smoothing
 
RealScalar m_angle
 
bool m_residual
 
RealScalar m_error
 
ComputationInfo m_info
 
bool m_isInitialized
 
Index m_iterations
 

Additional Inherited Members

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

Detailed Description

template<typename MatrixType_, typename Preconditioner_>
class Eigen::IDRS< MatrixType_, Preconditioner_ >

The Induced Dimension Reduction method (IDR(s)) is a short-recurrences Krylov method for sparse square problems.

This class allows to solve for A.x = b sparse linear problems. The vectors x and b can be either dense or sparse. he Induced Dimension Reduction method, IDR(), is a robust and efficient short-recurrence Krylov subspace method for solving large nonsymmetric systems of linear equations.

For indefinite systems IDR(S) outperforms both BiCGStab and BiCGStab(L). Additionally, IDR(S) can handle matrices with complex eigenvalues more efficiently than BiCGStab.

Many problems that do not converge for BiCGSTAB converge for IDR(s) (for larger values of s). And if both methods converge the convergence for IDR(s) is typically much faster for difficult systems (for example indefinite problems).

IDR(s) is a limited memory finite termination method. In exact arithmetic it converges in at most N+N/s iterations, with N the system size. It uses a fixed number of 4+3s vector. In comparison, BiCGSTAB terminates in 2N iterations and uses 7 vectors. GMRES terminates in at most N iterations, and uses I+3 vectors, with I the number of iterations. Restarting GMRES limits the memory consumption, but destroys the finite termination property.

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

The tolerance corresponds to the relative residual error: |Ax-b|/|b|

Performance: when using sparse matrices, best performance is achied 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) can also be used in a matrix-free context, see the following example .

See also
class SimplicialCholesky, DiagonalPreconditioner, IdentityPreconditioner

Member Typedef Documentation

◆ Base

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

◆ MatrixType

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

◆ Preconditioner

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

◆ RealScalar

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

◆ Scalar

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

Constructor & Destructor Documentation

◆ IDRS() [1/2]

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

Default constructor.

326 : m_S(4), m_smoothing(false), m_angle(RealScalar(0.7)), m_residual(false) {}
bool m_residual
Definition: IDRS.h:322
MatrixType::RealScalar RealScalar
Definition: IDRS.h:309
bool m_smoothing
Definition: IDRS.h:320
Index m_S
Definition: IDRS.h:319
RealScalar m_angle
Definition: IDRS.h:321

◆ IDRS() [2/2]

template<typename MatrixType_ , typename Preconditioner_ >
template<typename MatrixDerived >
Eigen::IDRS< MatrixType_, Preconditioner_ >::IDRS ( 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.
340  : Base(A.derived()), m_S(4), m_smoothing(false), m_angle(RealScalar(0.7)), m_residual(false) {}
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
IterativeSolverBase< IDRS > Base
Definition: IDRS.h:313

Member Function Documentation

◆ _solve_vector_with_guess_impl()

template<typename MatrixType_ , typename Preconditioner_ >
template<typename Rhs , typename Dest >
void Eigen::IDRS< 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
348  {
351 
353  m_residual);
354 
356  }
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 idrs(const MatrixType &A, const Rhs &b, Dest &x, const Preconditioner &precond, Index &iter, typename Dest::RealScalar &relres, Index S, bool smoothing, typename Dest::RealScalar angle, bool replacement)
Definition: IDRS.h:58
list x
Definition: plotDoE.py:28

References b, Eigen::internal::idrs(), Eigen::IDRS< MatrixType_, Preconditioner_ >::m_angle, Eigen::IDRS< MatrixType_, Preconditioner_ >::m_error, Eigen::IDRS< MatrixType_, Preconditioner_ >::m_info, Eigen::IDRS< MatrixType_, Preconditioner_ >::m_iterations, Eigen::IterativeSolverBase< Derived >::m_preconditioner, Eigen::IDRS< MatrixType_, Preconditioner_ >::m_residual, Eigen::IDRS< MatrixType_, Preconditioner_ >::m_S, Eigen::IDRS< MatrixType_, Preconditioner_ >::m_smoothing, Eigen::IterativeSolverBase< Derived >::m_tolerance, Eigen::IDRS< MatrixType_, Preconditioner_ >::matrix(), Eigen::IterativeSolverBase< Derived >::maxIterations(), Eigen::NoConvergence, Eigen::NumericalIssue, ret, Eigen::Success, and plotDoE::x.

◆ matrix()

template<typename MatrixType_ , typename Preconditioner_ >
const ActualMatrixType& Eigen::IterativeSolverBase< Derived >::matrix
inlineprivate

◆ setAngle()

template<typename MatrixType_ , typename Preconditioner_ >
void Eigen::IDRS< MatrixType_, Preconditioner_ >::setAngle ( RealScalar  angle)
inline

The angle must be a real scalar. In IDR(s), a value for the iteration parameter omega must be chosen in every s+1th step. The most natural choice is to select a value to minimize the norm of the next residual. This corresponds to the parameter omega = 0. In practice, this may lead to values of omega that are so small that the other iteration parameters cannot be computed with sufficient accuracy. In such cases it is better to increase the value of omega sufficiently such that a compromise is reached between accurate computations and reduction of the residual norm. The parameter angle =0.7 (”maintaining the convergence strategy”) results in such a compromise.

385 { m_angle = angle; }
double angle(const double &t)
Angular position as a function of time t.
Definition: jeffery_orbit.cc:98

References Jeffery_Solution::angle(), and Eigen::IDRS< MatrixType_, Preconditioner_ >::m_angle.

◆ setResidualUpdate()

template<typename MatrixType_ , typename Preconditioner_ >
void Eigen::IDRS< MatrixType_, Preconditioner_ >::setResidualUpdate ( bool  update)
inline

The parameter replace is a logical that determines whether a residual replacement strategy is employed to increase the accuracy of the solution.

390 { m_residual = update; }

References Eigen::IDRS< MatrixType_, Preconditioner_ >::m_residual.

◆ setS()

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

Sets the parameter S, indicating the dimension of the shadow space. Default is 4

359  {
360  if (S < 1) {
361  S = 4;
362  }
363 
364  m_S = S;
365  }
@ S
Definition: quadtree.h:62

References Eigen::IDRS< MatrixType_, Preconditioner_ >::m_S, and oomph::QuadTreeNames::S.

◆ setSmoothing()

template<typename MatrixType_ , typename Preconditioner_ >
void Eigen::IDRS< MatrixType_, Preconditioner_ >::setSmoothing ( bool  smoothing)
inline

Switches off and on smoothing. Residual smoothing results in monotonically decreasing residual norms at the expense of two extra vectors of storage and a few extra vector operations. Although monotonic decrease of the residual norms is a desirable property, the rate of convergence of the unsmoothed process and the smoothed process is basically the same. Default is off

373 { m_smoothing = smoothing; }

References Eigen::IDRS< MatrixType_, Preconditioner_ >::m_smoothing.

Member Data Documentation

◆ m_angle

template<typename MatrixType_ , typename Preconditioner_ >
RealScalar Eigen::IDRS< MatrixType_, Preconditioner_ >::m_angle
private

◆ m_error

template<typename MatrixType_ , typename Preconditioner_ >
RealScalar Eigen::IterativeSolverBase< Derived >::m_error
mutableprivate

◆ m_info

template<typename MatrixType_ , typename Preconditioner_ >
ComputationInfo Eigen::IterativeSolverBase< Derived >::m_info
mutableprivate

◆ m_isInitialized

template<typename MatrixType_ , typename Preconditioner_ >
bool Eigen::SparseSolverBase< Derived >::m_isInitialized
mutableprivate

◆ m_iterations

template<typename MatrixType_ , typename Preconditioner_ >
Index Eigen::IterativeSolverBase< Derived >::m_iterations
mutableprivate

◆ m_residual

template<typename MatrixType_ , typename Preconditioner_ >
bool Eigen::IDRS< MatrixType_, Preconditioner_ >::m_residual
private

◆ m_S

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

◆ m_smoothing

template<typename MatrixType_ , typename Preconditioner_ >
bool Eigen::IDRS< MatrixType_, Preconditioner_ >::m_smoothing
private

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