Afivo
0.3
|
Module with methods to perform time integration. More...
Data Types | |
interface | subr_feuler |
Interface for a generic forward Euler scheme for time integration. More... | |
interface | subr_implicit |
Interface for a implicit solver for integrating stiff terms over time. More... | |
Functions/Subroutines | |
subroutine, public | af_advance (tree, dt, dt_lim, time, i_cc, time_integrator, forward_euler, implicit_solver) |
Advance solution over dt using time_integrator. More... | |
Variables | |
integer, parameter, public | af_num_integrators = 8 |
[time_integration_schemes] More... | |
integer, parameter, public | af_forward_euler = 1 |
Forward Euler method. More... | |
integer, parameter, public | af_heuns_method = 2 |
Heun's method (AKA modified Euler's method, explicit trapezoidal rule), CFL coefficient of 1. See e.g. https://en.wikipedia.org/wiki/Heun%27s_method. More... | |
integer, parameter, public | af_midpoint_method = 3 |
Midpoint method, see e.g. https://en.wikipedia.org/wiki/Midpoint_method. More... | |
integer, parameter, public | af_ssprk33_method = 4 |
Optimal 3-stage third-order SSPRK method (Shu & Osher), CFL coefficient of 1, see e.g. https://doi.org/10.1137/S0036142902419284. More... | |
integer, parameter, public | af_ssprk43_method = 5 |
Optimal 4-stage third-order SSPRK method (Ruuth & Spiteri), CFL coefficient of 2, see e.g. https://doi.org/10.1137/S0036142902419284. More... | |
integer, parameter, public | af_imex_euler = 6 |
1st order IMEX method, forward Euler and then implicit solve More... | |
integer, parameter, public | af_imex_trapezoidal = 7 |
Trapezoidal IMEX method (2nd order) More... | |
integer, parameter, public | af_rk4_method = 8 |
Classic 4th order Runge Kutta method, see e.g. https://en.wikipedia.org/wiki/Runge%E2%80%93Kutta_methods. More... | |
character(len=af_nlen), dimension(af_num_integrators), public | af_integrator_names = [character(len=af_nlen) :: "forward_euler", "heuns_method", "midpoint_method", "ssprk33", "ssprk43", "imex_euler", "imex_trapezoidal", "rk4"] |
integer, dimension(af_num_integrators), parameter, public | af_advance_num_steps = [1, 2, 2, 3, 4, 1, 2, 4] |
How many steps the time integrators take. More... | |
integer, dimension(af_num_integrators), parameter | req_copies = af_advance_num_steps |
How many variable copies are required for the time integrators. More... | |
logical, dimension(af_num_integrators), parameter | req_implicit = [.false., .false., .false., .false., .false., .true., .true., .false.] |
Whether an implicit solver is required for the scheme. More... | |
Module with methods to perform time integration.
subroutine, public m_af_advance::af_advance | ( | type(af_t), intent(inout) | tree, |
real(dp), intent(in) | dt, | ||
real(dp), intent(out) | dt_lim, | ||
real(dp), intent(inout) | time, | ||
integer, dimension(:), intent(in) | i_cc, | ||
integer, intent(in) | time_integrator, | ||
procedure(subr_feuler) | forward_euler, | ||
procedure(subr_implicit), optional | implicit_solver | ||
) |
Advance solution over dt using time_integrator.
The user should supply a forward Euler method as documented in subr_feuler. The indices of the cell-centered variables that will be operated on should also be provided, so that higher-order schemes can be constructed automatically from the forward Euler method.
[in] | dt | Current time step |
[out] | dt_lim | Time step limit |
[in,out] | time | Current time |
[in] | i_cc | Index of cell-centered variables |
[in] | time_integrator | One of the pre-defined time integrators (e.g. af_heuns_method) |
forward_euler | Forward Euler method provided by the user | |
implicit_solver | Implicit solver to be used with IMEX schemes |
Definition at line 121 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_num_integrators = 8 |
[time_integration_schemes]
Definition at line 9 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_forward_euler = 1 |
Forward Euler method.
Definition at line 13 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_heuns_method = 2 |
Heun's method (AKA modified Euler's method, explicit trapezoidal rule), CFL coefficient of 1. See e.g. https://en.wikipedia.org/wiki/Heun%27s_method.
Definition at line 16 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_midpoint_method = 3 |
Midpoint method, see e.g. https://en.wikipedia.org/wiki/Midpoint_method.
Definition at line 18 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_ssprk33_method = 4 |
Optimal 3-stage third-order SSPRK method (Shu & Osher), CFL coefficient of 1, see e.g. https://doi.org/10.1137/S0036142902419284.
Definition at line 21 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_ssprk43_method = 5 |
Optimal 4-stage third-order SSPRK method (Ruuth & Spiteri), CFL coefficient of 2, see e.g. https://doi.org/10.1137/S0036142902419284.
Definition at line 24 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_imex_euler = 6 |
1st order IMEX method, forward Euler and then implicit solve
Definition at line 26 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_imex_trapezoidal = 7 |
Trapezoidal IMEX method (2nd order)
Definition at line 28 of file m_af_advance.f90.
integer, parameter, public m_af_advance::af_rk4_method = 8 |
Classic 4th order Runge Kutta method, see e.g. https://en.wikipedia.org/wiki/Runge%E2%80%93Kutta_methods.
[time_integration_schemes]
Definition at line 31 of file m_af_advance.f90.
character(len=af_nlen), dimension(af_num_integrators), public m_af_advance::af_integrator_names = [character(len=af_nlen) :: "forward_euler", "heuns_method", "midpoint_method", "ssprk33", "ssprk43", "imex_euler", "imex_trapezoidal", "rk4"] |
Definition at line 34 of file m_af_advance.f90.
integer, dimension(af_num_integrators), parameter, public m_af_advance::af_advance_num_steps = [1, 2, 2, 3, 4, 1, 2, 4] |
How many steps the time integrators take.
Definition at line 40 of file m_af_advance.f90.
|
private |
How many variable copies are required for the time integrators.
Definition at line 44 of file m_af_advance.f90.
|
private |
Whether an implicit solver is required for the scheme.
Definition at line 47 of file m_af_advance.f90.