ryujin 2.1.1 revision 955e869188d49b3c97ca7b1cf4fd9ceb0e6f46ef
|
#include <source/time_integrator.h>
Public Types | |
Typedefs and constexpr constants | |
using | HyperbolicSystem = typename Description::HyperbolicSystem |
using | View = typename Description::template HyperbolicSystemView< dim, Number > |
using | ParabolicSystem = typename Description::ParabolicSystem |
using | StateVector = typename View::StateVector |
Public Member Functions | |
Constructor and setup | |
TimeIntegrator (const MPIEnsemble &mpi_ensemble, const OfflineData< dim, Number > &offline_data, const HyperbolicModule< Description, dim, Number > &hyperbolic_module, const ParabolicModule< Description, dim, Number > ¶bolic_module, const std::string &subsection="/TimeIntegrator") | |
void | prepare () |
Functions for performing explicit time steps | |
Number | step (StateVector &state_vector, Number t, Number t_final=std::numeric_limits< Number >::max()) |
const auto & | time_stepping_scheme () const |
const auto & | efficiency () const |
Number | step_ssprk_22 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_ssprk_33 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_erk_11 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_erk_22 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_erk_33 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_erk_43 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_erk_54 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_strang_ssprk_33_cn (StateVector &state_vector, Number t, Number tau_max) |
Number | step_strang_erk_33_cn (StateVector &state_vector, Number t, Number tau_max) |
Number | step_strang_erk_43_cn (StateVector &state_vector, Number t, Number tau_max) |
Number | step_imex_11 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_imex_22 (StateVector &state_vector, Number t, Number tau_max) |
Number | step_imex_33 (StateVector &state_vector, Number t, Number tau_max) |
The TimeIntegrator class implements IMEX timestepping strategies based on explicit and diagonally-implicit Runge Kutta schemes.
Definition at line 242 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem |
Definition at line 250 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::View = typename Description::template HyperbolicSystemView<dim, Number> |
Definition at line 252 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::ParabolicSystem = typename Description::ParabolicSystem |
Definition at line 255 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::StateVector = typename View::StateVector |
Definition at line 257 of file time_integrator.h.
ryujin::TimeIntegrator< Description, dim, Number >::TimeIntegrator | ( | const MPIEnsemble & | mpi_ensemble, |
const OfflineData< dim, Number > & | offline_data, | ||
const HyperbolicModule< Description, dim, Number > & | hyperbolic_module, | ||
const ParabolicModule< Description, dim, Number > & | parabolic_module, | ||
const std::string & | subsection = "/TimeIntegrator< Description, dim, Number >" |
||
) |
Constructor.
Definition at line 29 of file time_integrator.template.h.
References ryujin::bang_bang_control, ryujin::erk_33, and ryujin::strang_erk_33_cn.
void ryujin::TimeIntegrator< Description, dim, Number >::prepare |
Prepare time integration. A call to prepare() allocates temporary storage and is necessary before any of the following time-stepping functions can be called.
Definition at line 75 of file time_integrator.template.h.
References ryujin::erk_11, ryujin::erk_22, ryujin::erk_33, ryujin::erk_43, ryujin::erk_54, ryujin::imex_11, ryujin::imex_22, ryujin::imex_33, ryujin::ssprk_22, ryujin::ssprk_33, ryujin::strang_erk_33_cn, ryujin::strang_erk_43_cn, and ryujin::strang_ssprk_33_cn.
Number ryujin::TimeIntegrator< Description, dim, Number >::step | ( | StateVector & | state_vector, |
Number | t, | ||
Number | t_final = std::numeric_limits<Number>::max() |
||
) |
Given a reference to a previous state vector U performs an explicit time step (and store the result in U). The function returns the chosen time step size tau. The time step size tau is selected such that $t + tau <= t_final$.
Definition at line 207 of file time_integrator.template.h.
References ryujin::bang_bang_control, ryujin::erk_11, ryujin::erk_22, ryujin::erk_33, ryujin::erk_43, ryujin::erk_54, ryujin::imex_11, ryujin::imex_22, ryujin::imex_33, ryujin::none, ryujin::raise_exception, ryujin::ssprk_22, ryujin::ssprk_33, ryujin::strang_erk_33_cn, ryujin::strang_erk_43_cn, ryujin::strang_ssprk_33_cn, and ryujin::warn.
|
inline |
The selected time-stepping scheme.
Definition at line 308 of file time_integrator.h.
|
inline |
The eficiency of the selected time-stepping scheme expressed as the ratio of step size of the combined method to step size of an elementary forward Euler step. For example, SSPRK33 has an efficiency ratio of 1 whereas ERK33 has an efficiency ratio of 3.
Definition at line 316 of file time_integrator.h.
|
protected |
Given a reference to a previous state vector U performs an explicit second-order strong-stability preserving Runge-Kutta SSPRK(2,2;1/2) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 279 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs an explicit third-order strong-stability preserving Runge-Kutta SSPRK(3,3;1/3) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 302 of file time_integrator.template.h.
References ryujin::sadd().
|
protected |
Given a reference to a previous state vector U performs an explicit first-order Euler step ERK(1,1;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 332 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs an explicit second-order Runge-Kutta ERK(2,2;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 350 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs an explicit third-order Runge-Kutta ERK(3,3;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 373 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs an explicit 4 stage third-order Runge-Kutta ERK(4,3;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 407 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs an explicit 4 stage fourth-order Runge-Kutta ERK(5,4;1) time step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 446 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs a combined explicit implicit Strang split using a third-order Runge-Kutta ERK(3,3;1/3) time step and an implicit Crank-Nicolson step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 516 of file time_integrator.template.h.
References ryujin::sadd().
|
protected |
Given a reference to a previous state vector U performs a combined explicit implicit Strang split using a third-order Runge-Kutta ERK(3,3;1) time step and an implicit Crank-Nicolson step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 564 of file time_integrator.template.h.
References ryujin::sadd().
|
protected |
Given a reference to a previous state vector U performs a combined explicit implicit Strang split using a third-order Runge-Kutta ERK(4,3;1) time step and an implicit Crank-Nicolson step (and store the result in U). The function returns the chosen time step size tau, which is guaranteed to be less than or equal to the parameter tau_max
.
Definition at line 620 of file time_integrator.template.h.
References ryujin::sadd().
|
protected |
Given a reference to a previous state vector U, performs an implicit-explicit step IMEX(1,1;1) using a forward euler scheme for the hyperbolic subproblem and backward euler scheme for the parabolic subproblem.
Definition at line 683 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs an implicit-explicit IMEX(2,2;1) step using a two stage midpoint rule for the hyperbolic subproblem and a two stage midpoint rule for the parabolic subproblem. The function returns the chosen time step size tau.
Definition at line 704 of file time_integrator.template.h.
|
protected |
Given a reference to a previous state vector U performs an implicit-explicit IMEX(3,3;1) step using a three stage ERK tableau for the hyperbolic subproblem and a three stage DIRK tableau for the parabolic subproblem. The function returns the chosen time step size tau.
Definition at line 737 of file time_integrator.template.h.