FOODIE integrator: provide an explicit class of low storage Runge-Kutta schemes, from 1st to 4th order accurate.
The integrators provided have the low storage property allowing for an efficient use of the memory. Following Williamson approach [1], the LSRK(5,4)2N (solution 3) scheme of Carpenter et al. [2] is implemented.
Considering the following ODE system:
$$ U_t = R(t,U) $$
where \(U_t = \frac{dU}{dt}\), U is the vector of state variables being a function of the time-like independent variable t, R is the (vectorial) residual function, the class of schemes implemented are written in the form:
$$\begin{matrix} K_1 = U^n \\ K_2 = 0 \\ \left.\begin{matrix} K_2 = A_s K_2 + \Delta t R(t^n + C_s \Delta t, K_1) \\ K_1 = K_1 + B_s K_2 \end{matrix}\right\} s=1,2,...N_s\\ U^{n+1} = K_1 \end{matrix}$$
where Ns is the number of stages used and \(K_1, K_2\) are the 2 registers used for stages computation.
The value of \(\Delta t\) must be provided, it not being computed by the integrator.
The schemes are explicit thus \(A_1=C_1=0\). The coefficients \(A_s, B_s, C_s\) are given in the Williamson low storage table form.
The different schemes are selected accordingly to the number of stages used. Currently the following schemes are available:
This scheme is TVD and reverts to Explicit Forward Euler, it being 1st order. It is not a real low storage scheme, this being meaningless for a first order scheme. However it is added for safety reason.
Stage |
A |
B |
C |
---|---|---|---|
1 |
0 |
1 |
0 |
This scheme is a low storage RK(5, 4), based on the solution 3 proposed in [2].
Stage |
A |
B |
C |
---|---|---|---|
1 |
0 |
1432997174477/9575080441755 |
0 |
2 |
-567301805773 /1357537059087 |
5161836677717/13612068292357 |
1432997174477/9575080441755 |
3 |
-2404267990393/2016746695238 |
1720146321549/2090206949498 |
2526269341429/6820363962896 |
4 |
-3550918686646/2091501179385 |
3134564353537/4481467310338 |
2006345519317/3224310063776 |
5 |
-1275806237668/842570457699 |
2277821191437/14882151754819 |
2802321613138/2924317926251 |
This scheme is a low storage RK(6, 4), by [3].
Stage |
A |
B |
C |
---|---|---|---|
1 |
0 |
0.122000000000 |
0 |
2 |
-0.691750960670 |
0.477263056358 |
0.122000000000 |
3 |
-1.727127405211 |
0.381941220320 |
0.269115878630 |
4 |
-0.694890150986 |
0.447757195744 |
0.447717183551 |
5 |
-1.039942756197 |
0.498614246822 |
0.749979795490 |
6 |
-1.531977447611 |
0.186648570846 |
0.898555413085 |
This scheme is a low storage RK(7, 4), by [3].
Stage |
A |
B |
C |
---|---|---|---|
1 |
0 |
0.117322146869 |
0 |
2 |
-0.647900745934 |
0.503270262127 |
0.117322146869 |
3 |
-2.704760863204 |
0.233663281658 |
0.294523230758 |
4 |
-0.460080550118 |
0.283419634625 |
0.305658622131 |
5 |
-0.500581787785 |
0.540367414023 |
0.582864148403 |
6 |
-1.906532255913 |
0.371499414620 |
0.858664273599 |
7 |
-1.450000000000 |
0.136670099385 |
0.868664273599 |
This scheme is a low storage RK(12, 4), by [4].
Stage |
A |
B |
C |
---|---|---|---|
1 |
0 |
0.0650008435125904 |
0 |
2 |
-0.0923311242368072 |
0.0161459902249842 |
0.0650008435125904 |
3 |
-0.9441056581158819 |
0.5758627178358159 |
0.0796560563081853 |
4 |
-4.3271273247576394 |
0.1649758848361671 |
0.1620416710085376 |
5 |
-2.1557771329026072 |
0.3934619494248182 |
0.2248877362907778 |
6 |
-0.9770727190189062 |
0.0443509641602719 |
0.2952293985641261 |
7 |
-0.7581835342571139 |
0.2074504268408778 |
0.3318332506149405 |
8 |
-1.7977525470825499 |
0.6914247433015102 |
0.4094724050198658 |
9 |
-2.6915667972700770 |
0.3766646883450449 |
0.6356954475753369 |
10 |
-4.6466798960268143 |
0.0757190350155483 |
0.6806551557645497 |
11 |
-0.1539613783825189 |
0.2027862031054088 |
0.7143773712418350 |
12 |
-0.5943293901830616 |
0.2167029365631842 |
0.9032588871651854 |
This scheme is a low storage RK(13, 4), by [4].
Stage |
A |
B |
C |
---|---|---|---|
1 |
0 |
0.0271990297818803 |
0 |
2 |
-0.6160178650170565 |
0.1772488819905108 |
0.0271990297818803 |
3 |
-0.4449487060774118 |
0.0378528418949694 |
0.0952594339119365 |
4 |
-1.0952033345276178 |
0.6086431830142991 |
0.1266450286591127 |
5 |
-1.2256030785959187 |
0.2154313974316100 |
0.1825883045699772 |
6 |
-0.2740182222332805 |
0.2066152563885843 |
0.3737511439063931 |
7 |
-0.0411952089052647 |
0.0415864076069797 |
0.5301279418422206 |
8 |
-0.1797084899153560 |
0.0219891884310925 |
0.5704177433952291 |
9 |
-1.1771530652064288 |
0.9893081222650993 |
0.5885784947099155 |
10 |
-0.4078831463120878 |
0.0063199019859826 |
0.6160769826246714 |
11 |
-0.8295636426191777 |
0.3749640721105318 |
0.6223252334314046 |
12 |
-4.7895970584252288 |
1.6080235151003195 |
0.6897593128753419 |
13 |
-0.6606671432964504 |
0.0961209123818189 |
0.9126827615920843 |
This scheme is a low storage RK(14, 4), by [4].
Stage |
A |
B |
C |
---|---|---|---|
1 |
0 |
0.0367762454319673 |
0 |
2 |
-0.7188012108672410 |
0.3136296607553959 |
0.0367762454319673 |
3 |
-0.7785331173421570 |
0.1531848691869027 |
0.1249685262725025 |
4 |
-0.0053282796654044 |
0.0030097086818182 |
0.2446177702277698 |
5 |
-0.8552979934029281 |
0.3326293790646110 |
0.2476149531070420 |
6 |
-3.9564138245774565 |
0.2440251405350864 |
0.2969311120382472 |
7 |
-1.5780575380587385 |
0.3718879239592277 |
0.3978149645802642 |
8 |
-2.0837094552574054 |
0.6204126221582444 |
0.5270854589440328 |
9 |
-0.7483334182761610 |
0.1524043173028741 |
0.6981269994175695 |
10 |
-0.7032861106563359 |
0.0760894927419266 |
0.8190890835352128 |
11 |
0.0013917096117681 |
0.0077604214040978 |
0.8527059887098624 |
12 |
-0.0932075369637460 |
0.0024647284755382 |
0.8604711817462826 |
13 |
-0.9514200470875948 |
0.0780348340049386 |
0.8627060376969976 |
14 |
-7.1151571693922548 |
5.5059777270269628 |
0.8734213127600976 |
[1] Low-Storage Runge-Kutta Schemes, J. H. Williamson, Journal of Computational Physics, vol. 35, 1980, pp. 48–56.
[2] Fourth-Order 2N-Storage Runge-Kutta Schemes, Mark H. Carpenter, Christopher A. Kennedy, NASA Technical Memorandum 109112, June 1994.
[3] High-accuracy large-step explicit Runge–Kutta (HALE-RK) schemes for computational aeroacoustics, Vasanth Allampalli and Ray Hixon and M. Nallasamy and Scott D. Sawyer, Journal of Computational Physics, vol. 228, 2009, pp. 3837–3850.
[4] Efficient low-storage Runge–Kutta schemes with optimized stability regions, Jens Niegemann and Richard Diehl and Kurt Busch, Journal of Computational Physics, vol. 231, 2012, pp. 364–372.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
character(len=99), | private, | parameter | :: | class_name_ | = | 'runge_kutta_ls' | Name of the class of schemes. |
logical, | private, | parameter | :: | has_fast_mode_ | = | .true. | Flag to check if integrator provides fast mode integrate. |
integer(kind=I_P), | private, | parameter | :: | registers | = | 2 | Registers used (2N schemes). |
character(len=99), | private, | parameter | :: | supported_schemes_(1:7) | = | [trim(class_name_)//'_stages_1_order_1 ', trim(class_name_)//'_stages_5_order_4 ', trim(class_name_)//'_stages_6_order_4 ', trim(class_name_)//'_stages_7_order_4 ', trim(class_name_)//'_stages_12_order_4', trim(class_name_)//'_stages_13_order_4', trim(class_name_)//'_stages_14_order_4'] | List of supported schemes. |
FOODIE integrator: provide an explicit class of low storage Runge-Kutta schemes, from 1st to 4th order accurate.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=R_P), | private, | allocatable | :: | A(:) | Low storage A coefficients. |
||
real(kind=R_P), | private, | allocatable | :: | B(:) | Low storage B coefficients. |
||
real(kind=R_P), | private, | allocatable | :: | C(:) | Low storage C coefficients. |
||
class(integrand_object), | public, | allocatable | :: | buffer | Buffer used for fast integration. |
||
character(len=:), | public, | allocatable | :: | description_ | Informative description of the integrator. |
||
integer(kind=I_P), | public | :: | error | = | 0 | Error status code. |
|
character(len=:), | public, | allocatable | :: | error_message | Error message, hopefully meaningful. |
||
integer(kind=I_P), | public | :: | registers | Number of registers used for stages. |
|||
class(integrand_object), | public, | allocatable | :: | stage(:) | Stages. |
||
integer(kind=I_P), | public | :: | stages | Number of stages. |
procedure, public, pass(self) :: allocate_integrand_members | Allocate integrand members. |
procedure, public, pass(lhs) :: assign_abstract | Assign ony members of abstract integrator_object type. |
procedure, public, pass(lhs) :: assign_multistage | Assign members of integrator_multistage_object and parents. |
generic, public :: assignment(=) => integr_assign_integr | Overload |
procedure, public, pass(self) :: check_error | Check for error occurrencies. |
procedure, public, pass(self) :: class_name | Return the class name of schemes. |
procedure, public, pass(self) :: description | Return informative integrator description. |
procedure, public, pass(self) :: destroy | Destroy the integrator. |
procedure, public, pass(self) :: destroy_abstract | Destroy only members of abstract integrator_object type. |
procedure, public, pass(self) :: destroy_multistage | Destroy the integrator. |
procedure, public, pass(self) :: has_fast_mode | Return .true. if the integrator class has fast mode integrate. |
procedure, public, pass(self) :: initialize | Initialize (create) the integrator. |
procedure, public, pass(lhs) :: integr_assign_integr | Operator |
procedure, public, pass(self) :: integrate | Integrate integrand field. |
procedure, public, pass(self) :: integrate_fast | Integrate integrand field, fast mode. |
procedure, public, pass(self) :: is_multistage | Return .true. for multistage integrator. |
procedure, public, pass(self) :: is_multistep | Return .true. for multistep integrator. |
procedure, public, pass(self) :: is_supported | Return .true. if the integrator class support the given scheme. |
procedure, public, nopass :: registers_number | Return the number of registers used. |
procedure, public, pass(self) :: stages_number | Return number of stages used. |
procedure, public, pass(self) :: steps_number | Return number of steps used. |
procedure, public, pass(self) :: supported_schemes | Return the list of supported schemes. |
procedure, public, pass(self) :: trigger_error | Trigger an error. |
Return the class name of schemes.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(in) | :: | self | Integrator. |
Class name.
Return .true. if the integrator class has fast mode integrate.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(in) | :: | self | Integrator. |
Inquire result.
Return .true. if the integrator class support the given scheme.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(in) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
Inquire result.
Return the number of registers used.
Number of registers used.
Return the list of supported schemes.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(in) | :: | self | Integrator. |
Queried scheme.
Destroy the integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(inout) | :: | self | Integrator. |
Create the actual RK integrator: initialize the Butcher' low storage table coefficients.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(inout) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
||
class(integrand_object), | intent(in), | optional | :: | U | Integrand molding prototype. |
|
logical, | intent(in), | optional | :: | stop_on_fail | Stop execution if initialization fail. |
Operator =
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(inout) | :: | lhs | Left hand side. |
||
class(integrator_object), | intent(in) | :: | rhs | Right hand side. |
Integrate field with explicit low storage Runge-Kutta scheme.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(inout) | :: | self | Integrator. |
||
class(integrand_object), | intent(inout) | :: | U | Field to be integrated. |
||
real(kind=R_P), | intent(in) | :: | Dt | Time step. |
||
real(kind=R_P), | intent(in) | :: | t | Time. |
||
real(kind=R_P), | intent(out), | optional | :: | new_Dt | New adapted time step. |
Integrate field with explicit low storage Runge-Kutta scheme, fast mode.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_ls), | intent(inout) | :: | self | Integrator. |
||
class(integrand_object), | intent(inout) | :: | U | Field to be integrated. |
||
real(kind=R_P), | intent(in) | :: | Dt | Time step. |
||
real(kind=R_P), | intent(in) | :: | t | Time. |
||
real(kind=R_P), | intent(out), | optional | :: | new_Dt | New adapted time step. |