The Induced Dimension Reduction method (IDR(s)) is a short-recurrences Krylov method for sparse square problems.
More...
|
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) |
|
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 |
|
Preconditioner & | preconditioner () |
|
const Preconditioner & | preconditioner () 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 () |
|
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 () |
|
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 |
This class follows the sparse solver concept .
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 304 of file IDRS.h.