Afivo  0.3
Data Types | Functions/Subroutines | Variables
m_af_advance Module Reference

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
 
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...
 

Detailed Description

Module with methods to perform time integration.

Function/Subroutine Documentation

◆ af_advance()

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.

Todo:
check whether time in sub-steps is accurate with some test
Parameters
[in]dtCurrent time step
[out]dt_limTime step limit
[in,out]timeCurrent time
[in]i_ccIndex of cell-centered variables
[in]time_integratorOne of the pre-defined time integrators (e.g. af_heuns_method)
forward_eulerForward Euler method provided by the user
implicit_solverImplicit solver to be used with IMEX schemes

Definition at line 118 of file m_af_advance.f90.

Variable Documentation

◆ af_num_integrators

integer, parameter, public m_af_advance::af_num_integrators = 8

Definition at line 9 of file m_af_advance.f90.

◆ af_forward_euler

integer, parameter, public m_af_advance::af_forward_euler = 1

Forward Euler method.

Definition at line 11 of file m_af_advance.f90.

◆ af_heuns_method

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 14 of file m_af_advance.f90.

◆ af_midpoint_method

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 16 of file m_af_advance.f90.

◆ af_ssprk33_method

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 19 of file m_af_advance.f90.

◆ af_ssprk43_method

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 22 of file m_af_advance.f90.

◆ af_imex_euler

integer, parameter, public m_af_advance::af_imex_euler = 6

1st order IMEX method, forward Euler and then implicit solve

Definition at line 24 of file m_af_advance.f90.

◆ af_imex_trapezoidal

integer, parameter, public m_af_advance::af_imex_trapezoidal = 7

Trapezoidal IMEX method (2nd order)

Definition at line 26 of file m_af_advance.f90.

◆ af_rk4_method

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.

Definition at line 29 of file m_af_advance.f90.

◆ af_integrator_names

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 31 of file m_af_advance.f90.

◆ af_advance_num_steps

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 37 of file m_af_advance.f90.

◆ req_copies

integer, dimension(af_num_integrators), parameter m_af_advance::req_copies = af_advance_num_steps
private

How many variable copies are required for the time integrators.

Definition at line 41 of file m_af_advance.f90.

◆ req_implicit

logical, dimension(af_num_integrators), parameter m_af_advance::req_implicit = [.false., .false., .false., .false., .false., .true., .true., .false.]
private

Whether an implicit solver is required for the scheme.

Definition at line 44 of file m_af_advance.f90.