DUNE PDELab (2.8)
Do one step of an explicit time-stepping scheme. More...
#include <dune/pdelab/instationary/explicitonestep.hh>
Public Member Functions | |
ExplicitOneStepMethod (const TimeSteppingParameterInterface< T > &method_, IGOS &igos_, LS &ls_) | |
construct a new one step scheme More... | |
ExplicitOneStepMethod (const TimeSteppingParameterInterface< T > &method_, IGOS &igos_, LS &ls_, double ls_reduct_) | |
Construct a new Explicit One Step Method object. More... | |
ExplicitOneStepMethod (const TimeSteppingParameterInterface< T > &method_, IGOS &igos_, LS &ls_, TC &tc_) | |
construct a new one step scheme More... | |
ExplicitOneStepMethod (const TimeSteppingParameterInterface< T > &method_, IGOS &igos_, LS &ls_, TC &tc_, double ls_reduct_) | |
Construct a new Explicit One Step Method object. More... | |
void | setVerbosityLevel (int level) |
change verbosity level; 0 means completely quiet | |
void | setStepNumber (int newstep) |
change number of current step | |
void | setReduction (const double &ls_reduction_) |
Set the Reduction for linear system soltion. More... | |
void | setMethod (const TimeSteppingParameterInterface< T > &method_) |
redefine the method to be used; can be done before every step More... | |
T | apply (T time, T dt, TrlV &xold, TrlV &xnew) |
do one step; More... | |
template<typename Limiter > | |
T | apply (T time, T dt, TrlV &xold, TrlV &xnew, Limiter &limiter) |
do one step; More... | |
template<typename F > | |
T | apply (T time, T dt, TrlV &xold, F &f, TrlV &xnew) |
do one step; More... | |
template<typename F , typename Limiter > | |
T | apply (T time, T dt, TrlV &xold, F &f, TrlV &xnew, Limiter &limiter) |
do one step; More... | |
Detailed Description
class Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >
Do one step of an explicit time-stepping scheme.
- Template Parameters
-
T type to represent time values IGOS assembler for instationary problems LS backend to solve linear system TrlV vector type to represent coefficients of solutions TstV vector type to represent residuals TC time controller class
Constructor & Destructor Documentation
◆ ExplicitOneStepMethod() [1/4]
|
inline |
construct a new one step scheme
- Parameters
-
method_ Parameter object. igos_ Assembler object (instationary grid operator space). pdesolver_ solver object (typically Newton).
The contructed method object stores references to the object it is constructed with, so these objects should be valid for as long as the constructed object is used (or until setMethod() is called, see there). Use SimpleTimeController that does not control the time step.
References DUNE_THROW, and Dune::PDELab::TimeSteppingParameterInterface< R >::implicit().
◆ ExplicitOneStepMethod() [2/4]
|
inline |
Construct a new Explicit One Step Method object.
- Note
- The contructed method object stores references to the object it is constructed with, so these objects should be valid for as long as the constructed object is used (or until setMethod() is called).
- Use SimpleTimeController that does not control the time step.
- Parameters
-
method_ A time stepping strategy igos_ An instationary local operator ls_ A linear solver backend ls_reduct_ Reduction for the linear solver
References DUNE_THROW, and Dune::PDELab::TimeSteppingParameterInterface< R >::implicit().
◆ ExplicitOneStepMethod() [3/4]
|
inline |
construct a new one step scheme
- Parameters
-
method_ Parameter object. igos_ Assembler object (instationary grid operator space). pdesolver_ solver object (typically Newton). tc_ a time controller object
The contructed method object stores references to the object it is constructed with, so these objects should be valid for as long as the constructed object is used (or until setMethod() is called, see there).
References DUNE_THROW, and Dune::PDELab::TimeSteppingParameterInterface< R >::implicit().
◆ ExplicitOneStepMethod() [4/4]
|
inline |
Construct a new Explicit One Step Method object.
- Note
- The contructed method object stores references to the object it is constructed with, so these objects should be valid for as long as the constructed object is used (or until setMethod() is called).
- Use SimpleTimeController that does not control the time step.
- Parameters
-
method_ A time stepping strategy igos_ An instationary local operator ls_ A linear solver backend ls_reduct_ Reduction for the linear solver
References DUNE_THROW, and Dune::PDELab::TimeSteppingParameterInterface< R >::implicit().
Member Function Documentation
◆ apply() [1/4]
|
inline |
do one step;
- Parameters
-
[in] time start of time step [in] dt suggested time step size [in] xold value at begin of time step [in] function to interpolate boundary condition from [in,out] xnew value at end of time step; contains initial guess for first substep on entry
- Returns
- time step size
References Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::apply().
◆ apply() [2/4]
|
inline |
do one step;
- Parameters
-
[in] time start of time step [in] dt suggested time step size [in] xold value at begin of time step [in] function to interpolate boundary condition from [in,out] xnew value at end of time step; contains initial guess for first substep on entry [in] limiter
- Returns
- time step size
◆ apply() [3/4]
|
inline |
do one step;
- Parameters
-
[in] time start of time step [in] dt suggested time step size [in] xold value at begin of time step [in,out] xnew value at end of time step; contains initial guess for first substep on entry
- Returns
- time step size
References Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::apply().
Referenced by Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::apply().
◆ apply() [4/4]
|
inline |
do one step;
- Parameters
-
[in] time start of time step [in] dt suggested time step size [in] xold value at begin of time step [in,out] xnew value at end of time step; contains initial guess for first substep on entry [in] limiter
- Returns
- time step size
References Dune::PDELab::TimeSteppingParameterInterface< R >::d(), Dune::Simd::Overloads::min(), Dune::PDELab::TimeSteppingParameterInterface< R >::name(), Dune::PDELab::TimeSteppingParameterInterface< R >::s(), and Dune::PDELab::TimeControllerInterface< R >::suggestTimestep().
◆ setMethod()
|
inline |
redefine the method to be used; can be done before every step
- Parameters
-
method_ Parameter object.
The OneStepMethod object stores a reference to the method_ object. The old method object is no longer referenced after this member function returns.
References DUNE_THROW, and Dune::PDELab::TimeSteppingParameterInterface< R >::implicit().
◆ setReduction()
|
inline |
Set the Reduction for linear system soltion.
This reduction is applied to the solution of the linear systems resulting from the shu-osher RK form for each stage. In cases when the mass matrix is known to be diagonal of block diagonal and an appropiate linear solver is provided (e.g. BlockJacobi, Gauss-Seidel), it's sufficient to set this to 0.99.
- Parameters
-
ls_reduction_ Reduction for the linear system solution
The documentation for this class was generated from the following file:
- dune/pdelab/instationary/explicitonestep.hh