DUNE PDELab (git)
Dune::CompleteFCGSolver< X > Class Template Reference
Complete flexible conjugate gradient 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_type > | scalar_real_type |
scalar type underlying the field_type | |
Public Member Functions | |
void | apply (X &x, X &b, 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... | |
SolverCategory::Category | category () const override |
Category of the solver (see SolverCategory::Category) | |
Detailed Description
template<class X>
class Dune::CompleteFCGSolver< X >
class Dune::CompleteFCGSolver< X >
Complete flexible conjugate gradient method.
This solver is a simple modification of the RestartedFCGSolver and, if possible, uses mmax old directions. It uses noticeably more memory, but provides more stability for preconditioner changes.
Member Function Documentation
◆ apply() [1/2]
|
inlineoverridevirtualinherited |
Apply inverse operator with given reduction factor.
apply inverse operator, with given convergence criteria.
- Warning
- Right hand side b may be overwritten!
- Parameters
-
x The left hand side to store the result in. b The right hand side reduction The minimum defect reduction to achieve. res Object to store the statistics about applying the operator.
- Exceptions
-
SolverAbort When the solver detects a problem and cannot continue
Implements Dune::InverseOperator< X, X >.
Reimplemented in Dune::RestartedGMResSolver< X, X >.
◆ apply() [2/2]
template<class X >
|
inlineoverridevirtual |
Apply inverse operator.
Apply inverse operator,.
- Warning
- Note: right hand side b may be overwritten!
- Parameters
-
x The left hand side to store the result in. b The right hand side res Object to store the statistics about applying the operator.
- Exceptions
-
SolverAbort When the solver detects a problem and cannot continue
- Note
- Currently, the RestartedFCGSolver aborts when a NaN or infinite defect is detected. However, -ffinite-math-only (implied by -ffast-math) can inhibit a result from becoming NaN that really should be NaN. E.g. numeric_limits<double>::quiet_NaN()*0.0==0.0 with gcc-5.3 -ffast-math.
Reimplemented from Dune::RestartedFCGSolver< X >.
References Dune::RestartedFCGSolver< X >::apply().
The documentation for this class was generated from the following file:
- dune/istl/solvers.hh
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 23, 23:29, 2024)