![]() |
ryujin 2.1.1 revision 9dcb748690310d6a540ebb8b066d1a0834fc7604
|
#include <source/time_integrator.h>
Public Types | |
using | HyperbolicSystem = typename Description::HyperbolicSystem |
using | ParabolicSystem = typename Description::ParabolicSystem |
using | HyperbolicSystemView = typename Description::HyperbolicSystem::template View< dim, Number > |
using | vector_type = MultiComponentVector< Number, problem_dimension > |
using | precomputed_type = MultiComponentVector< Number, n_precomputed_values > |
Public Member Functions | |
TimeIntegrator (const MPI_Comm &mpi_communicator, std::map< std::string, dealii::Timer > &computing_timer, 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 () |
Static Public Attributes | |
static constexpr unsigned int | problem_dimension |
static constexpr unsigned int | n_precomputed_values |
Functions for performing explicit time steps | |
Number | step (vector_type &U, Number t) |
auto & | time_stepping_scheme () const |
auto & | efficiency () const |
Number | step_ssprk_33 (vector_type &U, Number t) |
Number | step_erk_11 (vector_type &U, Number t) |
Number | step_erk_22 (vector_type &U, Number t) |
Number | step_erk_33 (vector_type &U, Number t) |
Number | step_erk_43 (vector_type &U, Number t) |
Number | step_erk_54 (vector_type &U, Number t) |
Number | step_strang_ssprk_33_cn (vector_type &U, Number t) |
Number | step_strang_erk_33_cn (vector_type &U, Number t) |
Number | step_strang_erk_43_cn (vector_type &U, Number t) |
The TimeIntegrator class implements IMEX timestepping strategies based on explicit and diagonally-implicit Runge Kutta schemes.
Definition at line 158 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::HyperbolicSystem = typename Description::HyperbolicSystem |
Definition at line 164 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::ParabolicSystem = typename Description::ParabolicSystem |
Definition at line 169 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::HyperbolicSystemView = typename Description::HyperbolicSystem::template View<dim, Number> |
Definition at line 174 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::vector_type = MultiComponentVector<Number, problem_dimension> |
Typedef for a MultiComponentVector storing the state U.
Definition at line 192 of file time_integrator.h.
using ryujin::TimeIntegrator< Description, dim, Number >::precomputed_type = MultiComponentVector<Number, n_precomputed_values> |
Typedef for a MultiComponentVector storing precomputed values.
Definition at line 197 of file time_integrator.h.
ryujin::TimeIntegrator< Description, dim, Number >::TimeIntegrator | ( | const MPI_Comm & | mpi_communicator, |
std::map< std::string, dealii::Timer > & | computing_timer, | ||
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 15 of file time_integrator.template.h.
References ryujin::bang_bang_control, ryujin::erk_33, ryujin::ScalarConservation::ParabolicSystem::is_identity, 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 62 of file time_integrator.template.h.
References ryujin::erk_11, ryujin::erk_22, ryujin::erk_33, ryujin::erk_43, ryujin::erk_54, ryujin::ScalarConservation::ParabolicSystem::is_identity, 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 | ( | vector_type & | U, |
Number | t | ||
) |
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.
Definition at line 183 of file time_integrator.template.h.
References ryujin::erk_11, ryujin::erk_22, ryujin::erk_33, ryujin::erk_43, and ryujin::ssprk_33.
|
inline |
The selected time-stepping scheme.
Definition at line 239 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 247 of file time_integrator.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.
If the parameter tau
is set to a nonzero value then the supplied value is used for time stepping instead of the computed maximal time step size.
Definition at line 244 of file time_integrator.template.h.
|
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.
Definition at line 272 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.
Definition at line 290 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.
Definition at line 318 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.
Definition at line 356 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.
Definition at line 404 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.
Definition at line 485 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) time step and an implicit Crank-Nicolson step (and store the result in U). The function returns the chosen time step size tau.
Definition at line 537 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(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.
Definition at line 605 of file time_integrator.template.h.
|
staticconstexpr |
Definition at line 180 of file time_integrator.h.
|
staticconstexpr |
Definition at line 186 of file time_integrator.h.