![]() |
ryujin 2.1.1 revision 9391072059490dd712e0ea92785f21acd6605f00
|
This module contains classes and functions used for timestepping and running the program.
|
strong |
Controls the chosen invariant domain / CFL recovery strategy.
Definition at line 23 of file time_integrator.h.
|
strong |
Controls the chosen time-stepping scheme.
Enumerator | |
---|---|
ssprk_22 | The strong stability preserving Runge Kutta method of order 2, SSPRK(2,2;1/2), with the following butcher tableau \begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{2} & \tfrac{1}{2} & 0 \\ \hline 1 & 1 & 0 \end{array} \end{align*} |
ssprk_33 | The strong stability preserving Runge Kutta method of order 3, SSPRK(3,3;1/3), with the following butcher tableau \begin{align*} \begin{array}{c|ccc} 0 & 0 \\ 1 & 1 & 0 \\ \tfrac{1}{2} & \tfrac{1}{4} & \tfrac{1}{4} & 0\\ \hline 1 & \tfrac{1}{6} & \tfrac{1}{6} & \tfrac{2}{3} \end{array} \end{align*} |
erk_11 | The explicit Runge-Kutta method RK(1,1;1), aka a simple, forward Euler step. |
erk_22 | The explicit Runge-Kutta method RK(2,2;1) with the butcher tableau \begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{2} & \tfrac{1}{2} & 0 \\ \hline 1 & 0 & 1 \end{array} \end{align*} |
erk_33 | The explicit Runge-Kutta method RK(3,3;1) with the butcher tableau \begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{3} & \tfrac{1}{3} & 0 \\ \tfrac{2}{3} & 0 & \tfrac{2}{3} & 0 \\ \hline 1 & \tfrac{1}{4} & 0 & \tfrac{3}{4} \end{array} \end{align*} |
erk_43 | The explicit Runge-Kutta method RK(4,3;1) with the butcher tableau \begin{align*} \begin{array}{c|ccc} 0 & 0 \\ \tfrac{1}{4} & \tfrac{1}{4} & 0 \\ \tfrac{1}{2} & 0 & \tfrac{1}{2} & 0 \\ \tfrac{3}{4} & 0 & \tfrac{1}{4} & \tfrac{1}{2} & 0 \\ \hline 1 & 0 & \tfrac{2}{3} & -\tfrac{1}{3} & \tfrac{2}{3} \end{array} \end{align*} |
erk_54 | The explicit Runge-Kutta method RK(5,4;1) with the butcher tableau TODO |
strang_ssprk_33_cn | A Strang split using ssprk 33 for the hyperbolic subproblem and Crank-Nicolson for the parabolic subproblem |
strang_erk_33_cn | A Strang split using erk 33 for the hyperbolic subproblem and Crank-Nicolson for the parabolic subproblem |
strang_erk_43_cn | A Strang split using erk 43 for the hyperbolic subproblem and Crank-Nicolson for the parabolic subproblem |
Definition at line 45 of file time_integrator.h.