Discontinuous Galerkin Library
#include "dg/algorithm.h"

Abstract timeloop independent of stepper and ODE. More...

Inheritance diagram for dg::aTimeloop< ContainerType >:
[legend]

Public Types

using value_type = dg::get_value_type< ContainerType >
 
using container_type = ContainerType
 

Public Member Functions

void integrate (value_type t0, const ContainerType &u0, value_type t1, ContainerType &u1)
 Integrate a differential equation between given bounds. More...
 
void integrate (value_type &t0, const ContainerType &u0, value_type t1, ContainerType &u1, enum to mode)
 Build your own timeloop. More...
 
value_type get_dt () const
 The current timestep. More...
 
virtual aTimeloopclone () const =0
 Abstract copy method that returns a copy of *this on the heap. More...
 
virtual ~aTimeloop ()
 

Protected Member Functions

 aTimeloop ()
 empty More...
 
 aTimeloop (const aTimeloop &)
 empty More...
 
aTimeloopoperator= (const aTimeloop &)
 return *this More...
 

Detailed Description

template<class ContainerType>
struct dg::aTimeloop< ContainerType >

Abstract timeloop independent of stepper and ODE.

This class enables to write abstract time-loops that are independent of the used stepper (e.g. dg::RungeKutta, dg::ExplicitMultistep, ...) and the differential equation in use. The recommended way to implement it is using a std::function that erases the latter tpyes and emulates the step function of the stepper type.

Template Parameters
ContainerTypeAny class for which a specialization of TensorTraits exists and which fulfills the requirements of the there defined data and execution policies derived from AnyVectorTag and AnyPolicyTag. Among others
  • dg::HVec (serial), dg::DVec (cuda / omp), dg::MHVec (mpi + serial) or dg::MDVec (mpi + cuda / omp)
  • std::vector<dg::DVec> (vector of shared device vectors), std::array<double, 4> (array of 4 doubles) or std::map < std::string, dg::DVec> ( a map of named vectors)
  • double (scalar) and other primitive types ...
If there are several ContainerTypes in the argument list, then TensorTraits must exist for all of them
See also
See The dg dispatch system for a detailed explanation of our type dispatch system

Member Typedef Documentation

◆ container_type

template<class ContainerType >
using dg::aTimeloop< ContainerType >::container_type = ContainerType

◆ value_type

template<class ContainerType >
using dg::aTimeloop< ContainerType >::value_type = dg::get_value_type<ContainerType>

Constructor & Destructor Documentation

◆ ~aTimeloop()

template<class ContainerType >
virtual dg::aTimeloop< ContainerType >::~aTimeloop ( )
inlinevirtual

◆ aTimeloop() [1/2]

template<class ContainerType >
dg::aTimeloop< ContainerType >::aTimeloop ( )
inlineprotected

empty

◆ aTimeloop() [2/2]

template<class ContainerType >
dg::aTimeloop< ContainerType >::aTimeloop ( const aTimeloop< ContainerType > &  )
inlineprotected

empty

Member Function Documentation

◆ clone()

template<class ContainerType >
virtual aTimeloop * dg::aTimeloop< ContainerType >::clone ( ) const
pure virtual

Abstract copy method that returns a copy of *this on the heap.

Returns
a copy of *this on the heap
See also
dg::ClonePtr

Implemented in dg::AdaptiveTimeloop< ContainerType >, dg::MultistepTimeloop< ContainerType >, and dg::SinglestepTimeloop< ContainerType >.

◆ get_dt()

template<class ContainerType >
value_type dg::aTimeloop< ContainerType >::get_dt ( ) const
inline

The current timestep.

Returns
The dt value at the end of the latest call to integrate. If integrate fails for some reason then return the timestep at which the failure happens. Undefined if integrate has not been called at least once.

◆ integrate() [1/2]

template<class ContainerType >
void dg::aTimeloop< ContainerType >::integrate ( value_type t0,
const ContainerType &  u0,
value_type  t1,
ContainerType &  u1,
enum to  mode 
)
inline

Build your own timeloop.

Integrate an ode from t = t0 until t >= t1 using a set of discrete steps with or without forcing a certain timestep to match t1 exactly. For example, if t0=0 and t1=1 and the timestep is dt=0.6 then the timestepper without forcing stops at t1=1.2 This behaviour is useful especially for adaptive timesteppers because it allows to implement timeloops with minimal interference with the controller

double time = t_begin;
double deltaT = (t_end - t_begin) / (double)maxout;
for( unsigned u=1; u<=maxout; u++)
{
timeloop.integrate( time, u0, t_begin + u*deltaT, u0,
// Takes as many steps as an uninterupted integrate
// from t_begin to t_end
// Warning: the following version does not(!) finish at t_end:
//timeloop.integrate( time, u0, time+deltaT, u0,
// u<maxout ? dg::to::at_least : dg::to::exact);
}
@ exact
match the ending exactly
@ at_least
integrate to the end or further
Note
This function is re-entrant, i.e. if you integrate from t0 to t1 and then from t1 to t2 the timestepper does not need to re-initialize on the second call
Parameters
t0(read-write) initial time on entry; on exit it is the value to where the ode is actually integrated, corresponding to u1
u0initial value at t0
t1(read-only) end time
u1(write only) contains the result corresponding to exactly or at least t1 on output (may alias u0)
modeeither integrate exactly to t1 or at least to t1. In dg::at_least mode the timestep is bound only by t1-t0
Attention
The function may throw dg::Error (or anything derived from std::exception).

◆ integrate() [2/2]

template<class ContainerType >
void dg::aTimeloop< ContainerType >::integrate ( value_type  t0,
const ContainerType &  u0,
value_type  t1,
ContainerType &  u1 
)
inline

Integrate a differential equation between given bounds.

Integrate an ode from t = t0 until t == t1 using a set of discrete steps and forcing the last timestep to land on t1 exactly.

Note
This function is re-entrant, i.e. if you integrate from t0 to t1 and then from t1 to t2 the timestepper does not need to re-initialize on the second call
Parameters
t0initial time
u0initial value at t0
t1end time
u1(write only) contains the result corresponding to t1 on output. May alias u0.
Note
May not work for a multistep integrator if the interval does not evenly multiply the (fixed) timestep
Attention
The function may throw dg::Error (or anything derived from std::exception).

◆ operator=()

template<class ContainerType >
aTimeloop & dg::aTimeloop< ContainerType >::operator= ( const aTimeloop< ContainerType > &  )
inlineprotected

return *this


The documentation for this struct was generated from the following file: