foodie_integrator_runge_kutta_lssp Module

module~~foodie_integrator_runge_kutta_lssp~~UsesGraph module~foodie_integrator_runge_kutta_lssp foodie_integrator_runge_kutta_lssp module~foodie_integrator_multistage_object foodie_integrator_multistage_object module~foodie_integrator_multistage_object->module~foodie_integrator_runge_kutta_lssp penf penf penf->module~foodie_integrator_runge_kutta_lssp 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_lssp module~foodie_integrator_object->module~foodie_integrator_multistage_object module~foodie_error_codes->module~foodie_integrator_runge_kutta_lssp module~foodie_integrand_object->module~foodie_integrator_runge_kutta_lssp 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 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:

Ns stages, Ns-1 order

$$ 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:

Computation of coefficients
  • allocate coefficients array, allocate(alpha(1:Ns));
  • initialize the first 2 elements with the coefficients of the two stages first order methods, alpha(1)=0, alpha(2)=1;
  • for i in [3,Ns] do:
  • alpha(i) = (2 / i) * alpha(i-1)
  • for j in [i-1,2,-1] do:
    • alpha(j) = (2 / (j-1)) * alpha(j-1)
  • for j in [2,i] do:
    • alpha(1) = 1 - alpha(j)

Ns stages, Ns order

$$ 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:

Computation of coefficients
  • allocate coefficients array, allocate(alpha(1:Ns));
  • initialize the first element with the coefficient of the one stage first order method, alpha(1)=1;
  • for i in [2,Ns] do:
  • alpha(i) = 1 / i!
  • for j in [i-1,2,-1] do:
    • alpha(j) = (1 / (j-1)) * alpha(j-1)
  • for j in [2,i] do:
    • alpha(1) = 1 - alpha(j)

Bibliography

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

Used By

module~~foodie_integrator_runge_kutta_lssp~~UsedByGraph module~foodie_integrator_runge_kutta_lssp foodie_integrator_runge_kutta_lssp module~foodie foodie module~foodie_integrator_runge_kutta_lssp->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_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

abstract interface

Abstract interfaces of integrator_runge_kutta_lssp methods.

  • private subroutine integrate_fast_interface(self, U, Dt, t)

    Integrate field with Linear SSP Runge-Kutta scheme, fast mode.

    Arguments

    Type IntentOptional AttributesName
    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 interface

Abstract interfaces of integrator_runge_kutta_lssp methods.

  • private subroutine integrate_interface(self, U, Dt, t)

    Integrate field with Linear SSP Runge-Kutta scheme.

    Arguments

    Type IntentOptional AttributesName
    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.


Derived Types

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

Components

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

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


Functions

private pure function class_name(self)

Return the class name of schemes.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_lssp), 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_lssp), 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_lssp), 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_lssp), 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_lssp), intent(inout) :: self

Integrator.

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

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

Arguments

Type IntentOptional AttributesName
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.

private elemental subroutine initialize_order_s(self)

Initialize integrator for s-th order formula.

Arguments

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

Integrator.

private elemental subroutine initialize_order_s_1(self)

Initialize integrator for (s-1)-th order formula.

Arguments

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

Integrator.

private subroutine integr_assign_integr(lhs, rhs)

Operator =.

Arguments

Type IntentOptional AttributesName
class(integrator_runge_kutta_lssp), 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 integrand field by Linear SSP Runge-Kutta methods.

Arguments

Type IntentOptional AttributesName
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.

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

Integrate integrand field by Linear SSP Runge-Kutta methods.

Arguments

Type IntentOptional AttributesName
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.

private subroutine integrate_order_s(self, U, Dt, t)

Integrate integrand field by s-th order formula.

Arguments

Type IntentOptional AttributesName
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.

private subroutine integrate_order_s_1(self, U, Dt, t)

Integrate integrand field by (s-1)-th order formula.

Arguments

Type IntentOptional AttributesName
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.

private subroutine integrate_order_s_1_fast(self, U, Dt, t)

Integrate integrand field by (s-1)-th order formula.

Arguments

Type IntentOptional AttributesName
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.

private subroutine integrate_order_s_fast(self, U, Dt, t)

Integrate integrand field by s-th order formula, fast mode.

Arguments

Type IntentOptional AttributesName
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.