DUNE-FEM (unstable)

implements the Generalized Minimal Residual (GMRes) method 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

 RestartedGMResSolver (const LinearOperator< X, Y > &op, Preconditioner< X, Y > &prec, scalar_real_type reduction, int restart, int maxit, int verbose)
 Set up RestartedGMResSolver solver. More...
 
 RestartedGMResSolver (const LinearOperator< X, Y > &op, const ScalarProduct< X > &sp, Preconditioner< X, Y > &prec, scalar_real_type reduction, int restart, int maxit, int verbose)
 Set up RestartedGMResSolver solver. More...
 
 RestartedGMResSolver (std::shared_ptr< const LinearOperator< X, Y > > op, std::shared_ptr< Preconditioner< X, X > > prec, const ParameterTree &configuration)
 Constructor. More...
 
 RestartedGMResSolver (std::shared_ptr< const LinearOperator< X, Y > > op, std::shared_ptr< const ScalarProduct< X > > sp, std::shared_ptr< Preconditioner< X, Y > > prec, scalar_real_type reduction, int restart, int maxit, int verbose)
 Set up RestartedGMResSolver solver. More...
 
void apply (X &x, Y &b, InverseOperatorResult &res) override
 Apply inverse operator. More...
 
void apply (X &x, Y &b, double reduction, InverseOperatorResult &res) override
 Apply inverse operator. More...
 
void apply (X &x, X &b, double reduction, InverseOperatorResult &res) override
 Apply inverse operator with given reduction factor. More...
 
virtual void apply (X &x, X &b, InverseOperatorResult &res)=0
 Apply inverse operator,. More...
 
SolverCategory::Category category () const override
 Category of the solver (see SolverCategory::Category)
 

Protected Types

using fAlloc = ReboundAllocatorType< X, field_type >
 field_type Allocator retrieved from domain type
 
using rAlloc = ReboundAllocatorType< X, real_type >
 real_type Allocator retrieved from domain type
 

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::RestartedGMResSolver< X, Y, F >

implements the Generalized Minimal Residual (GMRes) method

GMRes solves the unsymmetric linear system Ax = b using the Generalized Minimal Residual method as described the SIAM Templates book (http://www.netlib.org/templates/templates.pdf).

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

Constructor & Destructor Documentation

◆ RestartedGMResSolver() [1/4]

template<class X , class Y = X, class F = Y>
Dune::RestartedGMResSolver< X, Y, F >::RestartedGMResSolver ( const LinearOperator< X, Y > &  op,
Preconditioner< X, Y > &  prec,
scalar_real_type  reduction,
int  restart,
int  maxit,
int  verbose 
)
inline

Set up RestartedGMResSolver solver.

Parameters
restartnumber of GMRes cycles before restart

◆ RestartedGMResSolver() [2/4]

template<class X , class Y = X, class F = Y>
Dune::RestartedGMResSolver< X, Y, F >::RestartedGMResSolver ( const LinearOperator< X, Y > &  op,
const ScalarProduct< X > &  sp,
Preconditioner< X, Y > &  prec,
scalar_real_type  reduction,
int  restart,
int  maxit,
int  verbose 
)
inline

Set up RestartedGMResSolver solver.

Parameters
restartnumber of GMRes cycles before restart

◆ RestartedGMResSolver() [3/4]

template<class X , class Y = X, class F = Y>
Dune::RestartedGMResSolver< X, Y, F >::RestartedGMResSolver ( std::shared_ptr< const LinearOperator< X, Y > >  op,
std::shared_ptr< Preconditioner< X, X > >  prec,
const ParameterTree configuration 
)
inline

Constructor.

Additional parameter:

ParameterTree Key Meaning
restart number of GMRes cycles before restart

See ISTL_Factory for the ParameterTree layout and examples.

◆ RestartedGMResSolver() [4/4]

template<class X , class Y = X, class F = Y>
Dune::RestartedGMResSolver< X, Y, F >::RestartedGMResSolver ( std::shared_ptr< const LinearOperator< X, Y > >  op,
std::shared_ptr< const ScalarProduct< X > >  sp,
std::shared_ptr< Preconditioner< X, Y > >  prec,
scalar_real_type  reduction,
int  restart,
int  maxit,
int  verbose 
)
inline

Set up RestartedGMResSolver solver.

Parameters
restartnumber of GMRes cycles before restart

Member Function Documentation

◆ apply() [1/4]

void Dune::IterativeSolver< X, X >::apply ( X &  x,
X &  b,
double  reduction,
InverseOperatorResult res 
)
inlineoverridevirtualinherited

Apply inverse operator with given reduction factor.

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

Implements Dune::InverseOperator< X, X >.

Reimplemented in Dune::RestartedGMResSolver< X, X >.

◆ apply() [2/4]

virtual void Dune::InverseOperator< X, X >::apply ( X &  x,
X &  b,
InverseOperatorResult res 
)
pure virtualinherited

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

Implemented in Dune::LoopSolver< X >, Dune::GradientSolver< X >, Dune::CGSolver< X >, Dune::BiCGSTABSolver< X >, Dune::MINRESSolver< X >, Dune::RestartedGMResSolver< X, X >, Dune::GeneralizedPCGSolver< X >, Dune::RestartedFCGSolver< X >, and Dune::CompleteFCGSolver< X >.

◆ apply() [3/4]

template<class X , class Y = X, class F = Y>
void Dune::RestartedGMResSolver< 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 RestartedGMResSolver aborts when it detects a breakdown.

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

◆ apply() [4/4]

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

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.

References Dune::RestartedGMResSolver< X, Y, F >::apply(), and Dune::Simd::max().

Referenced by Dune::RestartedGMResSolver< X, Y, F >::apply().


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 (Nov 12, 23:30, 2024)