Dune::LoopSolver< X > Class Template Reference
[Iterative Solvers]
#include <solvers.hh>

Detailed Description
template<class X>
class Dune::LoopSolver< X >
Preconditioned loop solver.
Implements a preconditioned loop. Unsing this class every Preconditioner can be turned into a solver. The solver will apply one preconditioner step in each iteration loop.
Public Types | |
typedef X | domain_type |
The domain type of the operator that we do the inverse for. | |
typedef X | range_type |
The range type of the operator that we do the inverse for. | |
typedef X::field_type | field_type |
The field type of the operator that we do the inverse for. | |
Public Member Functions | |
template<class L, class P> | |
LoopSolver (L &op, P &prec, double reduction, int maxit, int verbose) | |
Set up Loop solver. | |
template<class L, class S, class P> | |
LoopSolver (L &op, S &sp, P &prec, double reduction, int maxit, int verbose) | |
Set up loop solver. | |
virtual void | apply (X &x, X &b, InverseOperatorResult &res) |
virtual void | apply (X &x, X &b, double reduction, InverseOperatorResult &res) |
Protected Member Functions | |
void | printHeader (std::ostream &s) const |
helper function for printing header of solver output | |
void | printOutput (std::ostream &s, const int iter, const DataType &norm, const DataType &norm_old) const |
helper function for printing solver output | |
void | printOutput (std::ostream &s, const int iter, const DataType &norm) const |
helper function for printing solver output |
Constructor & Destructor Documentation
Dune::LoopSolver< X >::LoopSolver | ( | L & | op, | |
P & | prec, | |||
double | reduction, | |||
int | maxit, | |||
int | verbose | |||
) | [inline] |
Set up Loop solver.
- Parameters:
-
op The operator we solve. prec The preconditioner to apply in each iteration of the loop. Has to inherit from Preconditioner. reduction The relative defect reduction to achieve when applying the operator. maxit The maximum number of iteration steps allowed when applying the operator. verbose The verbosity level.
- 0 : print nothing
- 1 : print initial and final defect and statistics
- 2 : print line for each iteration
References Dune::SolverCategory::sequential.
Dune::LoopSolver< X >::LoopSolver | ( | L & | op, | |
S & | sp, | |||
P & | prec, | |||
double | reduction, | |||
int | maxit, | |||
int | verbose | |||
) | [inline] |
Set up loop solver.
- Parameters:
-
op The operator we solve. sp The scalar product to use, e. g. SeqScalarproduct. prec The preconditioner to apply in each iteration of the loop. Has to inherit from Preconditioner. reduction The relative defect reduction to achieve when applying the operator. maxit The maximum number of iteration steps allowed when applying the operator. verbose The verbosity level.
- 0 : print nothing
- 1 : print initial and final defect and statistics
- 2 : print line for each iteration
Member Function Documentation
virtual void Dune::LoopSolver< X >::apply | ( | X & | x, | |
X & | b, | |||
InverseOperatorResult & | res | |||
) | [inline, virtual] |
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.
Implements Dune::InverseOperator< X, X >.
References Dune::Preconditioner< X, Y >::apply(), Dune::LinearOperator< X, Y >::applyscaleadd(), Dune::InverseOperatorResult::clear(), Dune::InverseOperatorResult::conv_rate, Dune::InverseOperatorResult::converged, Dune::InverseOperatorResult::elapsed, Dune::InverseOperatorResult::iterations, Dune::Preconditioner< X, Y >::post(), Dune::Preconditioner< X, Y >::pre(), Dune::InverseOperator< X, X >::printHeader(), Dune::InverseOperator< X, X >::printOutput(), and Dune::InverseOperatorResult::reduction.
virtual void Dune::LoopSolver< X >::apply | ( | X & | x, | |
X & | b, | |||
double | reduction, | |||
InverseOperatorResult & | res | |||
) | [inline, virtual] |
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.
Implements Dune::InverseOperator< X, X >.
The documentation for this class was generated from the following file: