Dune Core Modules (2.8.0)

implements the Flexible Generalized Minimal Residual (FGMRes) method (right preconditioned) More...

#include <dune/istl/solvers.hh>

Public Types

typedef X domain_type
 Type of the domain of the operator to be inverted.
 
typedef X range_type
 Type of the range of the operator to be inverted.
 
typedef X::field_type field_type
 The field type of the operator.
 
typedef FieldTraits< field_type >::real_type real_type
 The real type of the field type (is the same if using real numbers, but differs for std::complex)
 
typedef Simd::Scalar< real_typescalar_real_type
 scalar type underlying the field_type
 

Public Member Functions

void apply (X &x, Y &b, double reduction, InverseOperatorResult &res) override
 Apply inverse operator. More...
 
virtual void apply (X &x, X &b, InverseOperatorResult &res)
 Apply inverse operator. More...
 
virtual void apply (X &x, X &b, double reduction, InverseOperatorResult &res)
 Apply inverse operator. More...
 
virtual SolverCategory::Category category () const
 Category of the solver (see SolverCategory::Category)
 

Protected Member Functions

void printHeader (std::ostream &s) const
 helper function for printing header of solver output
 
void printOutput (std::ostream &s, const CountType &iter, const DataType &norm, const DataType &norm_old) const
 helper function for printing solver output
 
void printOutput (std::ostream &s, const CountType &iter, const DataType &norm) const
 helper function for printing solver output
 

Detailed Description

template<class X, class Y = X, class F = Y>
class Dune::RestartedFlexibleGMResSolver< X, Y, F >

implements the Flexible Generalized Minimal Residual (FGMRes) method (right preconditioned)

FGMRes solves the right-preconditioned unsymmetric linear system Ax = b using the Flexible Generalized Minimal Residual method. It is flexible because the preconditioner can change in every iteration, which allows to use Krylov solvers without fixed number of iterations as preconditioners. Needs more memory than GMRes.

Template Parameters
Xtrial vector, vector type of the solution
Ytest vector, vector type of the RHS
Fvector type for orthonormal basis of Krylov space

Member Function Documentation

◆ apply() [1/3]

virtual void Dune::RestartedGMResSolver< X, X , X >::apply ( X &  x,
X &  b,
double  reduction,
InverseOperatorResult res 
)
inlinevirtualinherited

Apply inverse operator.

apply inverse operator, with given convergence criteria.

Warning
Right hand side b may be overwritten!
Parameters
xThe left hand side to store the result in.
bThe right hand side
reductionThe minimum defect reduction to achieve.
resObject to store the statistics about applying the operator.
Exceptions
SolverAbortWhen the solver detects a problem and cannot continue
Note
Currently, the RestartedGMResSolver aborts when it detects a breakdown.

Reimplemented from Dune::IterativeSolver< X, X >.

◆ apply() [2/3]

virtual void Dune::RestartedGMResSolver< X, X , X >::apply ( X &  x,
X &  b,
InverseOperatorResult res 
)
inlinevirtualinherited

Apply inverse operator.

Apply inverse operator,.

Warning
Note: right hand side b may be overwritten!
Parameters
xThe left hand side to store the result in.
bThe right hand side
resObject to store the statistics about applying the operator.
Exceptions
SolverAbortWhen the solver detects a problem and cannot continue
Note
Currently, the RestartedGMResSolver aborts when it detects a breakdown.

Implements Dune::InverseOperator< X, X >.

◆ apply() [3/3]

template<class X , class Y = X, class F = Y>
void Dune::RestartedFlexibleGMResSolver< X, Y, F >::apply ( X &  x,
Y &  b,
double  reduction,
InverseOperatorResult res 
)
inlineoverride

Apply inverse operator.

apply inverse operator, with given convergence criteria.

Warning
Right hand side b may be overwritten!
Parameters
xThe left hand side to store the result in.
bThe right hand side
reductionThe minimum defect reduction to achieve.
resObject to store the statistics about applying the operator.
Exceptions
SolverAbortWhen the solver detects a problem and cannot continue
Note
Currently, the RestartedFlexibleGMResSolver aborts when it detects a breakdown.

References Dune::Simd::allTrue(), Dune::InverseOperatorResult::converged, and DUNE_THROW.


The documentation for this class was generated from the following file:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Dec 22, 23:30, 2024)