foodie_integrator_runge_kutta_ssp Module

module~~foodie_integrator_runge_kutta_ssp~~UsesGraph module~foodie_integrator_runge_kutta_ssp foodie_integrator_runge_kutta_ssp module~foodie_integrator_multistage_object foodie_integrator_multistage_object module~foodie_integrator_multistage_object->module~foodie_integrator_runge_kutta_ssp penf penf penf->module~foodie_integrator_runge_kutta_ssp penf->module~foodie_integrator_multistage_object module~foodie_integrator_object foodie_integrator_object penf->module~foodie_integrator_object module~foodie_error_codes foodie_error_codes penf->module~foodie_error_codes module~foodie_integrand_object foodie_integrand_object penf->module~foodie_integrand_object module~foodie_integrator_object->module~foodie_integrator_runge_kutta_ssp module~foodie_integrator_object->module~foodie_integrator_multistage_object module~foodie_error_codes->module~foodie_integrator_runge_kutta_ssp module~foodie_integrand_object->module~foodie_integrator_runge_kutta_ssp module~foodie_integrand_object->module~foodie_integrator_multistage_object iso_fortran_env iso_fortran_env iso_fortran_env->module~foodie_integrator_multistage_object iso_fortran_env->module~foodie_integrator_object
Help

FOODIE integrator: provide an explicit class of SSP Runge-Kutta schemes, from 1st to 4th order accurate.

The integrators provided have the Total Variation Diminishing (TVD) property or the Strong Stability Preserving (SSP) one. The schemes are explicit and defined through the Butcher's table syntax, see[1] .

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:

$$ U^{n+1} = U^n +\Delta t \sum_{s=1}^{Ns}\beta^s K^s $$

where Ns is the number of stages used and \(K^s\) is the \(s^{th}\) stage computed as:

$$ K^s = R\left( t^n+\gamma^s \Delta t, U^n +\Delta t \sum_{i=1}^{s-1}\alpha^{s,i} K^i \right) $$

The schemes are explicit thus the above summation is up to \(s-1\). The coefficients \(\beta\), \(\alpha\) and \(\gamma\) are given in the Butcher table form:

  gamma^1    | alpha^{1,1}       alpha^{1,2}       ...        alpha^{1,Ns}
  gamma^2    | alpha^{2,1}       alpha^{2,2}       ...        alpha^{2,Ns}
  .          | .                 .                 .          .
  .          | .                 .                  .         .
  .          | .                 .                   .        .
  gamma^{Ns} | alpha^{Ns,1}      alpha^{Ns,2}      ...        alpha^{Ns,Ns}
 ------------|-------------------------------------------------------------
             | beta^1            beta^2            ...        beta^{Ns}

Because only explicit schemes are considered the Butcher table reduces to diagonal matrix:

  gamma^1    | 0                 0                 ...        0
  gamma^2    | alpha^{2,1}       0                 ...        0
  .          | .                 .                 .          .
  .          | .                 .                  .         .
  .          | .                 .                   .        .
  gamma^{Ns} | alpha^{Ns,1}      alpha^{Ns,2}      ...        0
 ------------|-------------------------------------------------------------
             | beta^1            beta^2            ...        beta^{Ns}

Moreover the following relation always holds: \( \gamma^s = \sum_{i=1}^{Ns}\alpha^{s,i} \)

The different schemes are selected accordingly to the number of stages used. Currently the following schemes are available:

1 stage, Explicit Forward Euler, 1st order

This scheme is TVD and reverts to Explicit Forward Euler, it being 1st order.

  0 | 0
 ---|---
    | 1
2 stages, SSP, 2nd order

This scheme is an optmial SSP(2, 2) without low-storage algorithm, see [2].

  0 | 0     0
  1 | 1     0
 ---|-----------
    | 1/2   1/2
3 stages, SSP, 3rd order

This scheme is an optmial SSP(3, 3) without low-storage algorithm, see [2].

  0   | 0     0     0
  1   | 1     0     0
  1/2 | 1/4   1/4   0
 -----|-----------------
      | 1/6   1/6   1/3
5 stages, SSP, 4th order

This scheme is an optmial SSP(5, 4) without low-storage algorithm, see [2].

  0                | 0                  0                  0                  0                  0
  0.39175222700392 | 0.39175222700392   0                  0                  0                  0
  0.58607968896780 | 0.21766909633821   0.36841059262959   0                  0                  0
  0.47454236302687 | 0.08269208670950   0.13995850206999   0.25189177424738   0                  0
  0.93501063100924 | 0.06796628370320   0.11503469844438   0.20703489864929   0.54497475021237   0
 ------------------|---------------------------------------------------------------------------------------------
                   | 0.14681187618661   0.24848290924556   0.10425883036650   0.27443890091960   0.22600748319395

Bibliography

[1] Coefficients for the study of Runge-Kutta integration processes, Butcher, J.C., J. Austral. Math. Soc., Vol. 3, pages: 185–201, 1963.

[2] High Order Strong Stability Preserving Time Discretizations, Gottlieb, S., Ketcheson, D. I., Shu, C.W., Journal of Scientific Computing, vol. 38, N. 3, 2009, pp. 251-289.

Used By

module~~foodie_integrator_runge_kutta_ssp~~UsedByGraph module~foodie_integrator_runge_kutta_ssp foodie_integrator_runge_kutta_ssp module~foodie foodie module~foodie_integrator_runge_kutta_ssp->module~foodie program~integrate_burgers integrate_burgers module~foodie->program~integrate_burgers program~integrate_lorenz integrate_lorenz module~foodie->program~integrate_lorenz module~type_euler_1d type_euler_1D module~foodie->module~type_euler_1d module~foodie_test_integrand_tester_object foodie_test_integrand_tester_object module~foodie->module~foodie_test_integrand_tester_object module~foodie_test_lcce_test foodie_test_lcce_test module~foodie->module~foodie_test_lcce_test module~type_euler_1d_openmp type_euler_1D_openmp module~foodie->module~type_euler_1d_openmp program~integrate_euler_1d_caf~2 integrate_euler_1D_caf module~foodie->program~integrate_euler_1d_caf~2 module~foodie_test_oscillation_test foodie_test_oscillation_test module~foodie->module~foodie_test_oscillation_test program~integrate_euler_1d integrate_euler_1D module~foodie->program~integrate_euler_1d module~foodie_test_object foodie_test_object module~foodie->module~foodie_test_object module~foodie_test_integrand_oscillation foodie_test_integrand_oscillation module~foodie->module~foodie_test_integrand_oscillation module~foodie_test_integrand_lcce foodie_test_integrand_lcce module~foodie->module~foodie_test_integrand_lcce module~type_euler_1d_caf type_euler_1D_caf module~foodie->module~type_euler_1d_caf module~type_burgers type_burgers module~foodie->module~type_burgers module~foodie_test_integrand_ladvection foodie_test_integrand_ladvection module~foodie->module~foodie_test_integrand_ladvection program~integrate_euler_1d_openmp integrate_euler_1D_openmp module~foodie->program~integrate_euler_1d_openmp module~type_lorenz type_lorenz module~foodie->module~type_lorenz module~type_euler_1d->program~integrate_euler_1d module~foodie_test_integrand_tester_object->module~foodie_test_object module~foodie_test_integrand_tester_object->module~foodie_test_integrand_oscillation module~foodie_test_integrand_tester_object->module~foodie_test_integrand_lcce module~foodie_test_integrand_tester_object->module~foodie_test_integrand_ladvection program~foodie_test_lcce foodie_test_lcce module~foodie_test_lcce_test->program~foodie_test_lcce module~type_euler_1d_openmp->program~integrate_euler_1d_openmp program~foodie_test_oscillation foodie_test_oscillation module~foodie_test_oscillation_test->program~foodie_test_oscillation program~foodie_tester foodie_tester module~foodie_test_object->program~foodie_tester module~foodie_test_integrand_oscillation->module~foodie_test_oscillation_test module~foodie_test_integrand_oscillation->module~foodie_test_object module~foodie_test_integrand_lcce->module~foodie_test_lcce_test module~foodie_test_integrand_lcce->module~foodie_test_object module~type_euler_1d_caf->program~integrate_euler_1d_caf~2 module~type_burgers->program~integrate_burgers module~foodie_test_integrand_ladvection->module~foodie_test_object module~type_lorenz->program~integrate_lorenz
Help


Variables

TypeVisibility AttributesNameInitial
character(len=99), private, parameter:: class_name_ ='runge_kutta_ssp'

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:4) =[trim(class_name_)//'_stages_1_order_1', trim(class_name_)//'_stages_2_order_2', trim(class_name_)//'_stages_3_order_3', trim(class_name_)//'_stages_5_order_4']

List of supported schemes.


Derived Types

FOODIE integrator: provide an explicit class of SSP Runge-Kutta schemes, from 1st to 4th order accurate.

Components

TypeVisibility AttributesNameInitial
real(kind=R_P), private, allocatable:: alph(:,:)

\(\alpha\) Butcher's coefficients.

real(kind=R_P), private, allocatable:: beta(:)

\(\beta\) Butcher's 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.

real(kind=R_P), private, allocatable:: gamm(:)

\(\gamma\) Butcher's coefficients.

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.

Type-Bound Procedures

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


Functions

private pure function class_name(self)

Return the class name of schemes.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), intent(in) :: self

Integrator.

Return Value character(len=99)

Class name.

private elemental function has_fast_mode(self)

Return .true. if the integrator class has fast mode integrate.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), intent(in) :: self

Integrator.

Return Value logical

Inquire result.

private elemental function is_supported(self, scheme)

Return .true. if the integrator class support the given scheme.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), intent(in) :: self

Integrator.

character(len=*), intent(in) :: scheme

Selected scheme.

Return Value logical

Inquire result.

private pure function supported_schemes(self) result(schemes)

Return the list of supported schemes.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), intent(in) :: self

Integrator.

Return Value character(len=99), allocatable, (:)

Queried scheme.


Subroutines

private elemental subroutine destroy(self)

Destroy the integrator.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), intent(inout) :: self

Integrator.

private subroutine initialize(self, scheme, U, stop_on_fail)

Create the actual RK integrator: initialize the Butcher' table coefficients.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), 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.

private subroutine integr_assign_integr(lhs, rhs)

Operator =.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), intent(inout) :: lhs

Left hand side.

class(integrator_object), intent(in) :: rhs

Right hand side.

private subroutine integrate(self, U, Dt, t, new_Dt)

Integrate field with explicit SSP Runge-Kutta scheme.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), 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.

private subroutine integrate_fast(self, U, Dt, t, new_Dt)

Integrate field with explicit SSP Runge-Kutta scheme.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_ssp), 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.