DUNE PDELab (git)

Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC > Class Template Reference

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...
 
apply (T time, T dt, TrlV &xold, TrlV &xnew)
 do one step; More...
 
template<typename Limiter >
apply (T time, T dt, TrlV &xold, TrlV &xnew, Limiter &limiter)
 do one step; More...
 
template<typename F >
apply (T time, T dt, TrlV &xold, F &f, TrlV &xnew)
 do one step; More...
 
template<typename F , typename Limiter >
apply (T time, T dt, TrlV &xold, F &f, TrlV &xnew, Limiter &limiter)
 do one step; More...
 

Detailed Description

template<class T, class IGOS, class LS, class TrlV, class TstV = TrlV, class TC = SimpleTimeController<T>>
class Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >

Do one step of an explicit time-stepping scheme.

Template Parameters
Ttype to represent time values
IGOSassembler for instationary problems
LSbackend to solve linear system
TrlVvector type to represent coefficients of solutions
TstVvector type to represent residuals
TCtime controller class

Constructor & Destructor Documentation

◆ ExplicitOneStepMethod() [1/4]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::ExplicitOneStepMethod ( const TimeSteppingParameterInterface< T > &  method_,
IGOS &  igos_,
LS &  ls_ 
)
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]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::ExplicitOneStepMethod ( const TimeSteppingParameterInterface< T > &  method_,
IGOS &  igos_,
LS &  ls_,
double  ls_reduct_ 
)
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
See also
ExplicitOneStepMethod::setReduction

References DUNE_THROW, and Dune::PDELab::TimeSteppingParameterInterface< R >::implicit().

◆ ExplicitOneStepMethod() [3/4]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::ExplicitOneStepMethod ( const TimeSteppingParameterInterface< T > &  method_,
IGOS &  igos_,
LS &  ls_,
TC &  tc_ 
)
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]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::ExplicitOneStepMethod ( const TimeSteppingParameterInterface< T > &  method_,
IGOS &  igos_,
LS &  ls_,
TC &  tc_,
double  ls_reduct_ 
)
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
See also
ExplicitOneStepMethod::setReduction

References DUNE_THROW, and Dune::PDELab::TimeSteppingParameterInterface< R >::implicit().

Member Function Documentation

◆ apply() [1/4]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
template<typename F >
T Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::apply ( time,
dt,
TrlV &  xold,
F &  f,
TrlV &  xnew 
)
inline

do one step;

Parameters
[in]timestart of time step
[in]dtsuggested time step size
[in]xoldvalue at begin of time step
[in]functionto interpolate boundary condition from
[in,out]xnewvalue 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]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
template<typename F , typename Limiter >
T Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::apply ( time,
dt,
TrlV &  xold,
F &  f,
TrlV &  xnew,
Limiter &  limiter 
)
inline

do one step;

Parameters
[in]timestart of time step
[in]dtsuggested time step size
[in]xoldvalue at begin of time step
[in]functionto interpolate boundary condition from
[in,out]xnewvalue at end of time step; contains initial guess for first substep on entry
[in]limiter
Returns
time step size

◆ apply() [3/4]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
T Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::apply ( time,
dt,
TrlV &  xold,
TrlV &  xnew 
)
inline

do one step;

Parameters
[in]timestart of time step
[in]dtsuggested time step size
[in]xoldvalue at begin of time step
[in,out]xnewvalue 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]

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
template<typename Limiter >
T Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::apply ( time,
dt,
TrlV &  xold,
TrlV &  xnew,
Limiter &  limiter 
)
inline

do one step;

Parameters
[in]timestart of time step
[in]dtsuggested time step size
[in]xoldvalue at begin of time step
[in,out]xnewvalue 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::Hybrid::min, Dune::PDELab::TimeSteppingParameterInterface< R >::name(), Dune::PDELab::TimeSteppingParameterInterface< R >::s(), and Dune::PDELab::TimeControllerInterface< R >::suggestTimestep().

◆ setMethod()

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
void Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::setMethod ( const TimeSteppingParameterInterface< T > &  method_)
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()

template<class T , class IGOS , class LS , class TrlV , class TstV = TrlV, class TC = SimpleTimeController<T>>
void Dune::PDELab::ExplicitOneStepMethod< T, IGOS, LS, TrlV, TstV, TC >::setReduction ( const double &  ls_reduction_)
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:
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)