FOODIE integrator: provide an explicit class of Linear SSP Runge-Kutta schemes, from 1st to s-th order accurate.
The integrators provided have the Total Variation Diminishing (TVD) property or the Strong Stability Preserving (SSP) one. The schemes are explicit.
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 following 2 forms:
$$ U^0 = U^n $$ $$ U^k = U^{k-1} + \frac{1}{2} \Delta t R(U^{k-1}) \quad k=1,..,Ns $$ $$ U^{Ns} = U^{Ns} + \frac{1}{2} \Delta t R(U^{Ns}) $$ $$ U^{n+1} = \sum_{k=1}^{Ns}\alpha_k U^k $$
where Ns is the number of stages used and \(U^k\) is the \(k^{th}\) stage. The minimum number of stages is 2. The coefficients \(\alpha\) are computed by means of the following recursive algorith:
allocate(alpha(1:Ns))
;alpha(1)=0, alpha(2)=1
;i
in [3,Ns] do:alpha(i) = (2 / i) * alpha(i-1)
j
in [i-1,2,-1] do:alpha(j) = (2 / (j-1)) * alpha(j-1)
j
in [2,i] do:alpha(1) = 1 - alpha(j)
$$ U^0 = U^n $$ $$ U^k = U^{k-1} + \Delta t R(U^{k-1}) \quad k=1,..,Ns $$ $$ U^{Ns} = U^{Ns} +\Delta t R(U^{Ns}) $$ $$ U^{n+1} = \sum_{k=1}^{Ns}\alpha_k U^k $$
where Ns is the number of stages used and \(U^k\) is the \(k^{th}\) stage. The minimum number of stages is 1. The coefficients \(\alpha\) are computed by means of the following recursive algorith:
allocate(alpha(1:Ns))
;alpha(1)=1
;i
in [2,Ns] do:alpha(i) = 1 / i!
j
in [i-1,2,-1] do:alpha(j) = (1 / (j-1)) * alpha(j-1)
j
in [2,i] do:alpha(1) = 1 - alpha(j)
The value of \(\Delta t\) must be provided, it not being computed by the integrator.
[2] Strong Stability Preserving Runge-Kutta and Multistep Time Discretizations, S. Gottlieb, D. Ketcheson, C.W. Shu, 2011, 978-981-4289-26-9, doi:10.1142/7498, World Scientific Publishing Co. Pte. Ltd.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
character(len=99), | private, | parameter | :: | class_name_ | = | 'runge_kutta_lssp' | Name of the class of schemes. |
logical, | private, | parameter | :: | has_fast_mode_ | = | .true. | Flag to check if integrator provides fast mode integrate. |
character(len=99), | private, | parameter | :: | supported_schemes_(1:2) | = | [trim(class_name_)//'_stages_s_order_s_1', trim(class_name_)//'_stages_s_order_s '] | List of supported schemes. |
Abstract interfaces of integrator_runge_kutta_lssp methods.
Integrate field with Linear SSP Runge-Kutta scheme, fast mode.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(inout) | :: | self | Integrator. |
||
class(integrand_object), | intent(inout) | :: | U | Field to be integrated. |
||
real(kind=R_P), | intent(in) | :: | Dt | Time steps. |
||
real(kind=R_P), | intent(in) | :: | t | Times. |
Abstract interfaces of integrator_runge_kutta_lssp methods.
Integrate field with Linear SSP Runge-Kutta scheme.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(inout) | :: | self | Integrator. |
||
class(integrand_object), | intent(inout) | :: | U | Field to be integrated. |
||
real(kind=R_P), | intent(in) | :: | Dt | Time steps. |
||
real(kind=R_P), | intent(in) | :: | t | Times. |
FOODIE integrator: provide an explicit class of Linear SSP Runge-Kutta schemes, from 1st to s-th order accurate.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=R_P), | private, | allocatable | :: | alpha(:) | \(\alpha\) 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. |
||
procedure(integrate_interface), | private, | pointer | :: | integrate_ | => | integrate_order_s_1 | Integrate integrand field. |
procedure(integrate_fast_interface), | private, | pointer | :: | integrate_fast_ | => | integrate_order_s_1_fast | Integrate integrand field, fast. |
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, private, pass(self) :: initialize_order_s | Integrate integrator for s-th order formula. |
procedure, private, pass(self) :: initialize_order_s_1 | Integrate integrator for (s-1)-th order formula. |
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, private, pass(self) :: integrate_order_s | Integrate integrand field by s-th order formula. |
procedure, private, pass(self) :: integrate_order_s_1 | Integrate integrand field by (s-1)-th order formula. |
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, 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_lssp), | 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_lssp), | 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_lssp), | intent(in) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
Inquire result.
Return the list of supported schemes.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(in) | :: | self | Integrator. |
Queried scheme.
Destroy the integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(inout) | :: | self | Integrator. |
Create the actual RK integrator: initialize the Butcher' table coefficients.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(inout) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
||
class(integrand_object), | intent(in), | optional | :: | U | Integrand molding prototype. |
|
integer(kind=I_P), | intent(in), | optional | :: | stages | Stages number. |
|
logical, | intent(in), | optional | :: | stop_on_fail | Stop execution if initialization fail. |
Initialize integrator for s-th order formula.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(inout) | :: | self | Integrator. |
Initialize integrator for (s-1)-th order formula.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(inout) | :: | self | Integrator. |
Operator =
.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | intent(inout) | :: | lhs | Left hand side. |
||
class(integrator_object), | intent(in) | :: | rhs | Right hand side. |
Integrate integrand field by Linear SSP Runge-Kutta methods.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | 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 integrand field by Linear SSP Runge-Kutta methods.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | 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 integrand field by s-th order formula.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | 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. |
Integrate integrand field by (s-1)-th order formula.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | 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. |
Integrate integrand field by (s-1)-th order formula.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | 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. |
Integrate integrand field by s-th order formula, fast mode.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_runge_kutta_lssp), | 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. |