dune-istl
2.4.1
|
implements the Generalized Minimal Residual (GMRes) method More...
#include <dune/istl/solvers.hh>
Public Types | |
typedef X | domain_type |
The domain type of the operator to be inverted. More... | |
typedef Y | range_type |
The range type of the operator to be inverted. More... | |
typedef X::field_type | field_type |
The field type of the operator to be inverted. More... | |
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) More... | |
typedef F | basis_type |
The field type of the basis vectors. More... | |
Public Member Functions | |
template<class L , class P > | |
RestartedGMResSolver (L &op, P &prec, real_type reduction, int restart, int maxit, int verbose, bool recalc_defect) | |
template<class L , class P > | |
RestartedGMResSolver (L &op, P &prec, real_type reduction, int restart, int maxit, int verbose) | |
Set up solver. More... | |
template<class L , class S , class P > | |
RestartedGMResSolver (L &op, S &sp, P &prec, real_type reduction, int restart, int maxit, int verbose, bool recalc_defect) | |
template<class L , class S , class P > | |
RestartedGMResSolver (L &op, S &sp, P &prec, real_type reduction, int restart, int maxit, int verbose) | |
Set up solver. More... | |
virtual void | apply (X &x, Y &b, InverseOperatorResult &res) |
Apply inverse operator,. More... | |
virtual void | apply (X &x, Y &b, double reduction, InverseOperatorResult &res) |
Apply inverse operator. More... | |
Protected Types | |
enum | { iterationSpacing = 5, normSpacing = 16 } |
Protected Member Functions | |
void | printHeader (std::ostream &s) const |
helper function for printing header of solver output More... | |
template<typename CountType , typename DataType > | |
void | printOutput (std::ostream &s, const CountType &iter, const DataType &norm, const DataType &norm_old) const |
helper function for printing solver output More... | |
template<typename CountType , typename DataType > | |
void | printOutput (std::ostream &s, const CountType &iter, const DataType &norm) const |
helper function for printing solver output More... | |
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).
X | trial vector, vector type of the solution |
Y | test vector, vector type of the RHS |
F | vector type for orthonormal basis of Krylov space |
typedef F Dune::RestartedGMResSolver< X, Y, F >::basis_type |
The field type of the basis vectors.
typedef X Dune::RestartedGMResSolver< X, Y, F >::domain_type |
The domain type of the operator to be inverted.
typedef X::field_type Dune::RestartedGMResSolver< X, Y, F >::field_type |
The field type of the operator to be inverted.
typedef Y Dune::RestartedGMResSolver< X, Y, F >::range_type |
The range type of the operator to be inverted.
typedef FieldTraits<field_type>::real_type Dune::RestartedGMResSolver< X, Y, F >::real_type |
The real type of the field type (is the same if using real numbers, but differs for std::complex)
|
inline |
References Dune::SolverCategory::sequential.
|
inline |
Set up solver.
restart | number of GMRes cycles before restart |
References Dune::SolverCategory::sequential.
|
inline |
|
inline |
Set up solver.
restart | number of GMRes cycles before restart |
|
inlinevirtual |
Apply inverse operator,.
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. |
Implements Dune::InverseOperator< X, Y >.
|
inlinevirtual |
Apply inverse operator.
apply inverse operator, with given convergence criteria.
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. |
Implements Dune::InverseOperator< X, Y >.
References Dune::InverseOperatorResult::clear(), Dune::InverseOperatorResult::conv_rate, Dune::InverseOperatorResult::converged, Dune::InverseOperatorResult::elapsed, Dune::InverseOperatorResult::iterations, Dune::InverseOperator< X, Y >::printHeader(), Dune::InverseOperator< X, Y >::printOutput(), and Dune::InverseOperatorResult::reduction.
|
inlineprotectedinherited |
helper function for printing header of solver output
Referenced by Dune::RestartedGMResSolver< X, Y, F >::apply().
|
inlineprotectedinherited |
helper function for printing solver output
Referenced by Dune::RestartedGMResSolver< X, Y, F >::apply().
|
inlineprotectedinherited |
helper function for printing solver output