FOODIE integrator: provide a predictor-corrector class of Adams-Bashforth-Moutlon multi-step schemes, from 1st to 4rd order accurate.
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 Adams-Bashforth-Moulton class scheme implemented is:
$$ U^{n+N_s^p} = U^{n+N_s^p-1} +\Delta t \left[ \sum_{s=1}^{N_s^p}{ b_s^p \cdot R(t^{n+s-1}, U^{n+s-1}) } \right] $$
$$ U^{n+N_s^c} = U^{n+N_s^c-1} +\Delta t \left[ \sum_{s=0}^{N_s^c}{ b_s^c \cdot R(t^{n+s}, U^{n+s}) } \right] $$
where \(N_s^p\) is the number of previous steps considered for the predictor and \(N_s^c\) is the number of previous steps considered for the corrector.
The value of \(\Delta t\) must be provided, it not being computed by the integrator.
The coefficients \(b_s^{p,c}\) define the actual scheme, that is selected accordingly to the number of steps used. The predictor and corrector schemes should have the same formal order of accuracy, thus \(N_s^p=N_s^c+1\) should hold.
Currently, the following schemes are available:
This scheme is TVD and reverts to Explicit/Implicit Farward/Backward Euler, it being 1st order. The b coefficient is: $$b^p = \left[b_1\right] = \left[1\right]$$ $$b^c = \left[b_0\right] = \left[1\right]$$ The scheme is: $$ U^{n+1,p} = U^{n} + \Delta t R(t^{n},U^{n}) $$ $$ U^{n+1,c} = U^{n} + \Delta t R(t^{n+1},U^{n+1,p}) $$
This scheme is 2nd order. The b coefficients are: $$b^p = \left[ {\begin{array}{*{20}{c}} b_1 & b_2 \end{array}} \right] = \left[ {\begin{array}{*{20}{c}} -\frac{1}{2} & \frac{3}{2} \end{array}} \right]$$ $$b^c = \left[ {\begin{array}{*{20}{c}} b_0 & b_1 \end{array}} \right] = \left[ {\begin{array}{*{20}{c}} \frac{1}{2} & \frac{1}{2} \end{array}} \right]$$ The scheme is: $$ U^{n+2,p} = U^{n+1} +\Delta t \left[ \frac{3}{2} R(t^{n+1}, U^{n+1})-\frac{1}{2} R(t^{n}, U^{n}) \right] $$ $$ U^{n+2,c} = U^{n+1} +\Delta t \left[ \frac{1}{2} R(t^{n+2,p}, U^{n+2})+\frac{1}{2} R(t^{n+1}, U^{n+1}) \right] $$
This scheme is 3rd order. The b coefficients are: $$b^p = \left[ {\begin{array}{*{20}{c}} b_1 & b_2 & b_3 \end{array}} \right] = \left[ {\begin{array}{*{20}{c}} \frac{5}{12} & -\frac{4}{3} & \frac{23}{12} \end{array}} \right]$$ $$b^c = \left[ {\begin{array}{*{20}{c}} b_0 & b_1 & b_2 \end{array}} \right] = \left[ {\begin{array}{*{20}{c}} -\frac{1}{12} & \frac{2}{3} & \frac{5}{12} \end{array}} \right]$$ The scheme is: $$ U^{n+3,p} = U^{n+2} +\Delta t \left[ \frac{23}{12}R(t^{n+2}, U^{n+2}) - \frac{4}{3}R(t^{n+1}, U^{n+1}) +\frac{5}{12} R(t^{n}, U^{n}) \right] $$ $$ U^{n+3,c} = U^{n+2} +\Delta t \left[ \frac{5}{12}R(t^{n+3}, U^{n+3,p}) + \frac{2}{3}R(t^{n+2}, U^{n+2}) -\frac{1}{12} R(t^{n+1}, U^{n+1}) \right] $$
This scheme is 4th order. The b coefficients are: $$b^p = \left[ {\begin{array}{*{20}{c}} b_1 & b_2 & b_3 & b_4 \end{array}} \right] = \left[ {\begin{array}{*{20}{c}} -\frac{9}{24} & \frac{37}{24} & -\frac{59}{24} & \frac{55}{24} \end{array}} \right]$$ $$b^c = \left[ {\begin{array}{*{20}{c}} b_0 & b_1 & b_2 & b_3 \end{array}} \right] = \left[ {\begin{array}{*{20}{c}} \frac{1}{24} & -\frac{5}{24} & \frac{19}{24} & \frac{9}{24} \end{array}} \right]$$ The scheme is: $$ U^{n+4,p} = U^{n+3} +\Delta t \left[ \frac{55}{24}R(t^{n+3}, U^{n+3}) - \frac{59}{24}R(t^{n+2}, U^{n+2}) +\frac{37}{24} R(t^{n+1}, U^{n+1}) - \frac{9}{24} R(t^{n}, U^{n}) \right] $$ $$ U^{n+4,c} = U^{n+3} +\Delta t \left[ \frac{9}{24}R(t^{n+4}, U^{n+4,p}) + \frac{19}{24}R(t^{n+3}, U^{n+3}) -\frac{5}{24} R(t^{n+2}, U^{n+2}) + \frac{1}{24} R(t^{n+1}, U^{n+1}) \right] $$
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
character(len=99), | private, | parameter | :: | class_name_ | = | 'adams_bashforth_moulton' | 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:16) | = | [trim(class_name_)//'_1 ', trim(class_name_)//'_2 ', trim(class_name_)//'_3 ', trim(class_name_)//'_4 ', trim(class_name_)//'_5 ', trim(class_name_)//'_6 ', trim(class_name_)//'_7 ', trim(class_name_)//'_8 ', trim(class_name_)//'_9 ', trim(class_name_)//'_10', trim(class_name_)//'_11', trim(class_name_)//'_12', trim(class_name_)//'_13', trim(class_name_)//'_14', trim(class_name_)//'_15', trim(class_name_)//'_16'] | List of supported schemes. |
FOODIE integrator: provide an explicit class of Adams-Bashforth-Moulton multi-step schemes, from 1st to 4rd order accurate.
Type | Visibility | Attributes | Name | Initial | |||
---|---|---|---|---|---|---|---|
real(kind=R_P), | public, | allocatable | :: | Dt(:) | Previous time steps. |
||
logical, | public | :: | autoupdate | Perform cyclic autoupdate of previous time steps buffers. |
|||
class(integrand_object), | public, | allocatable | :: | buffer | Buffer used for fast integration. |
||
type(integrator_adams_moulton), | private | :: | corrector | Corrector solver. |
|||
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 | :: | iterations | Implicit iterations. |
|||
type(integrator_adams_bashforth), | private | :: | predictor | Predictor solver. |
|||
class(integrand_object), | public, | allocatable | :: | previous(:) | Previous steps. |
||
integer(kind=I_P), | public | :: | registers | Number of registers used for steps. |
|||
integer(kind=I_P), | public | :: | steps | Number of time steps. |
|||
real(kind=R_P), | public, | allocatable | :: | t(:) | Previous times. |
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_multistep | Assign members of integrator_multistep_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_multistep | 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. |
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) :: scheme_number | Return the scheme number in the list of supported schemes. |
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. |
procedure, public, nopass :: update_previous | Cyclic update previous time steps. |
Return the class name of schemes.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(in) | :: | self | Integrator. |
Class name.
Return .true. if the integrator class has fast mode integrate.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(in) | :: | self | Integrator. |
Inquire result.
Return .true. if the integrator class support the given scheme.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(in) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
Inquire result.
Return the scheme number in the list of supported schemes.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(in) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
Scheme number in the list of supported schemes.
Return the list of supported schemes.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(in) | :: | self | Integrator. |
Queried scheme.
Destroy the integrator.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(inout) | :: | self | Integrator. |
Create the actual Adams-Bashforth-Moulton integrator: initialize the b coefficients.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | intent(inout) | :: | self | Integrator. |
||
character(len=*), | intent(in) | :: | scheme | Selected scheme. |
||
integer(kind=I_P), | intent(in), | optional | :: | iterations | Implicit iterations. |
|
logical, | intent(in), | optional | :: | autoupdate | Enable cyclic autoupdate of previous time steps. |
|
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_adams_bashforth_moulton), | intent(inout) | :: | lhs | Left hand side. |
||
class(integrator_object), | intent(in) | :: | rhs | Right hand side. |
Integrate field with Adams-Bashforth-Moulton class scheme.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | 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. |
Integrate field with Adams-Bashforth-Moulton class scheme, fast mode.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
class(integrator_adams_bashforth_moulton), | 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. |