Public Types | Public Member Functions | Private Types | Private Attributes | List of all members
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 >:
Inheritance graph
[legend]

Public Types

typedef _MatrixType MatrixType
 
typedef _Preconditioner Preconditioner
 
typedef MatrixType::RealScalar RealScalar
 
typedef MatrixType::Scalar Scalar
 
- 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::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
 
 IDRS ()
 
template<typename MatrixDerived >
 IDRS (const EigenBase< MatrixDerived > &A)
 
void setAngle (RealScalar angle)
 
void setResidualUpdate (bool update)
 
void setS (Index S)
 
void setSmoothing (bool smoothing)
 
- Public Member Functions inherited from Eigen::IterativeSolverBase< IDRS< _MatrixType, _Preconditioner > >
void _solve_impl (const Rhs &b, Dest &x) const
 
void _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const
 
internal::enable_if< Rhs::ColsAtCompileTime!=1 &&DestDerived::ColsAtCompileTime!=1 >::type _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const
 
internal::enable_if< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 >::type _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const
 
IDRS< _MatrixType, _Preconditioner > & analyzePattern (const EigenBase< MatrixDerived > &A)
 
EIGEN_CONSTEXPR Index cols () const EIGEN_NOEXCEPT
 
IDRS< _MatrixType, _Preconditioner > & compute (const EigenBase< MatrixDerived > &A)
 
RealScalar error () const
 
IDRS< _MatrixType, _Preconditioner > & factorize (const EigenBase< MatrixDerived > &A)
 
ComputationInfo info () const
 
Index iterations () const
 
 IterativeSolverBase ()
 
 IterativeSolverBase (const EigenBase< MatrixDerived > &A)
 
Index maxIterations () const
 
Preconditionerpreconditioner ()
 
const Preconditionerpreconditioner () const
 
EIGEN_CONSTEXPR Index rows () const EIGEN_NOEXCEPT
 
IDRS< _MatrixType, _Preconditioner > & setMaxIterations (Index maxIters)
 
IDRS< _MatrixType, _Preconditioner > & setTolerance (const RealScalar &tolerance)
 
const SolveWithGuess< IDRS< _MatrixType, _Preconditioner >, Rhs, Guess > solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const
 
RealScalar tolerance () const
 
 ~IterativeSolverBase ()
 
- Public Member Functions inherited from Eigen::SparseSolverBase< IDRS< _MatrixType, _Preconditioner > >
void _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const
 
IDRS< _MatrixType, _Preconditioner > & derived ()
 
const IDRS< _MatrixType, _Preconditioner > & derived () const
 
const Solve< IDRS< _MatrixType, _Preconditioner >, Rhs > solve (const MatrixBase< Rhs > &b) const
 
const Solve< IDRS< _MatrixType, _Preconditioner >, Rhs > solve (const SparseMatrixBase< Rhs > &b) const
 
 SparseSolverBase ()
 
 ~SparseSolverBase ()
 

Private Types

typedef IterativeSolverBase< IDRSBase
 

Private Attributes

RealScalar m_angle
 
bool m_residual
 
Index m_S
 
bool m_smoothing
 

Additional Inherited Members

- Protected Types inherited from Eigen::IterativeSolverBase< IDRS< _MatrixType, _Preconditioner > >
typedef MatrixWrapper::ActualMatrixType ActualMatrixType
 
typedef SparseSolverBase< IDRS< _MatrixType, _Preconditioner > > Base
 
typedef internal::generic_matrix_wrapper< MatrixTypeMatrixWrapper
 
- Protected Member Functions inherited from Eigen::IterativeSolverBase< IDRS< _MatrixType, _Preconditioner > >
void grab (const InputType &A)
 
void init ()
 
const ActualMatrixTypematrix () const
 
- Protected Attributes inherited from Eigen::IterativeSolverBase< IDRS< _MatrixType, _Preconditioner > >
bool m_analysisIsOk
 
RealScalar m_error
 
bool m_factorizationIsOk
 
ComputationInfo m_info
 
Index m_iterations
 
MatrixWrapper m_matrixWrapper
 
Index m_maxIterations
 
Preconditioner m_preconditioner
 
RealScalar m_tolerance
 
- Protected Attributes inherited from Eigen::SparseSolverBase< IDRS< _MatrixType, _Preconditioner > >
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
_MatrixTypethe type of the sparse matrix A, can be a dense or a sparse matrix.
_Preconditionerthe 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.

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

Definition at line 275 of file IDRS.h.

Member Typedef Documentation

◆ Base

template<typename _MatrixType, typename _Preconditioner>
typedef IterativeSolverBase<IDRS> Eigen::IDRS< _MatrixType, _Preconditioner >::Base
private

Definition at line 341 of file IDRS.h.

◆ MatrixType

template<typename _MatrixType, typename _Preconditioner>
typedef _MatrixType Eigen::IDRS< _MatrixType, _Preconditioner >::MatrixType

Definition at line 335 of file IDRS.h.

◆ Preconditioner

template<typename _MatrixType, typename _Preconditioner>
typedef _Preconditioner Eigen::IDRS< _MatrixType, _Preconditioner >::Preconditioner

Definition at line 338 of file IDRS.h.

◆ RealScalar

template<typename _MatrixType, typename _Preconditioner>
typedef MatrixType::RealScalar Eigen::IDRS< _MatrixType, _Preconditioner >::RealScalar

Definition at line 337 of file IDRS.h.

◆ Scalar

template<typename _MatrixType, typename _Preconditioner>
typedef MatrixType::Scalar Eigen::IDRS< _MatrixType, _Preconditioner >::Scalar

Definition at line 336 of file IDRS.h.

Constructor & Destructor Documentation

◆ IDRS() [1/2]

template<typename _MatrixType, typename _Preconditioner>
Eigen::IDRS< _MatrixType, _Preconditioner >::IDRS ( )
inline

Default constructor.

Definition at line 354 of file IDRS.h.

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

Definition at line 367 of file IDRS.h.

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 tolerence and maxIterations
  2. Calls the function that has the core solver routine

Definition at line 377 of file IDRS.h.

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

Definition at line 419 of file IDRS.h.

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

Definition at line 427 of file IDRS.h.

◆ 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

Definition at line 388 of file IDRS.h.

◆ 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

Definition at line 404 of file IDRS.h.

Member Data Documentation

◆ m_angle

template<typename _MatrixType, typename _Preconditioner>
RealScalar Eigen::IDRS< _MatrixType, _Preconditioner >::m_angle
private

Definition at line 349 of file IDRS.h.

◆ m_residual

template<typename _MatrixType, typename _Preconditioner>
bool Eigen::IDRS< _MatrixType, _Preconditioner >::m_residual
private

Definition at line 350 of file IDRS.h.

◆ m_S

template<typename _MatrixType, typename _Preconditioner>
Index Eigen::IDRS< _MatrixType, _Preconditioner >::m_S
private

Definition at line 347 of file IDRS.h.

◆ m_smoothing

template<typename _MatrixType, typename _Preconditioner>
bool Eigen::IDRS< _MatrixType, _Preconditioner >::m_smoothing
private

Definition at line 348 of file IDRS.h.


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


gtsam
Author(s):
autogenerated on Tue Jul 4 2023 02:41:47