Dune Core Modules (2.4.1)

Generalized preconditioned conjugate gradient solver. More...

#include <dune/istl/solvers.hh>

Public Types

typedef X domain_type
 The domain type of the operator to be inverted.
 
typedef X range_type
 The range type of the operator to be inverted.
 
typedef X::field_type field_type
 The field type of the operator to be inverted.
 
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)
 

Public Member Functions

template<class L , class P >
 GeneralizedPCGSolver (L &op, P &prec, real_type reduction, int maxit, int verbose, int restart=10)
 Set up nonlinear preconditioned conjugate gradient solver. More...
 
template<class L , class P , class S >
 GeneralizedPCGSolver (L &op, S &sp, P &prec, real_type reduction, int maxit, int verbose, int restart=10)
 Set up nonlinear preconditioned conjugate gradient solver. 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 with given reduction factor. More...
 

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 Dune::GeneralizedPCGSolver< X >

Generalized preconditioned conjugate gradient solver.

A preconditioned conjugate gradient that allows the preconditioner to change between iterations.

One example for such preconditioner is AMG when used without a direct coarse solver. In this case the number of iterations performed on the coarsest level might change between applications.

In contrast to CGSolver the search directions are stored and the orthogonalization is done explicitly.

Constructor & Destructor Documentation

◆ GeneralizedPCGSolver() [1/2]

template<class X >
template<class L , class P >
Dune::GeneralizedPCGSolver< X >::GeneralizedPCGSolver ( L &  op,
P &  prec,
real_type  reduction,
int  maxit,
int  verbose,
int  restart = 10 
)
inline

Set up nonlinear preconditioned conjugate gradient solver.

Parameters
restartWhen to restart the construction of the Krylov search space.

References Dune::SolverCategory::sequential.

◆ GeneralizedPCGSolver() [2/2]

template<class X >
template<class L , class P , class S >
Dune::GeneralizedPCGSolver< X >::GeneralizedPCGSolver ( L &  op,
S &  sp,
P &  prec,
real_type  reduction,
int  maxit,
int  verbose,
int  restart = 10 
)
inline

Set up nonlinear preconditioned conjugate gradient solver.

Parameters
restartWhen to restart the construction of the Krylov search space.

Member Function Documentation

◆ apply() [1/2]

template<class X >
virtual void Dune::GeneralizedPCGSolver< X >::apply ( X &  x,
X &  b,
double  reduction,
InverseOperatorResult res 
)
inlinevirtual

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.

Implements Dune::InverseOperator< X, X >.

◆ apply() [2/2]


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)