rk8_10_class Derived Type

type, public, extends(rk_class) :: rk8_10_class

8th order Runge-Kutta method.


Inherits

type~~rk8_10_class~~InheritsGraph type~rk8_10_class rk8_10_class type~rk_class rk_class type~rk8_10_class->type~rk_class

Components

Type Visibility Attributes Name Initial
integer, public :: n = 0

user specified number of variables

procedure(deriv_func), public, pointer :: f => null()

user-specified derivative function

procedure(report_func), public, pointer :: report => null()

user-specified report function

procedure(event_func), public, pointer :: g => null()

event function (stop when this is zero)


Type-Bound Procedures

procedure, public :: initialize

initialize the class (set n,f, and report)

  • private subroutine initialize(me, n, f, report, g)

    Initialize the rk_class.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_class), intent(inout) :: me
    integer, intent(in) :: n

    number of variables

    procedure(deriv_func) :: f

    derivative function

    procedure(report_func), optional :: report

    for reporting the steps

    procedure(event_func), optional :: g

    for stopping at an event

procedure, public, non_overridable :: integrate

main integration routine

  • private subroutine integrate(me, t0, x0, h, tf, xf)

    Main integration routine for the rk_class.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_class), intent(inout) :: me
    real(kind=wp), intent(in) :: t0

    initial time

    real(kind=wp), intent(in), dimension(me%n) :: x0

    initial state

    real(kind=wp), intent(in) :: h

    abs(time step)

    real(kind=wp), intent(in) :: tf

    final time

    real(kind=wp), intent(out), dimension(me%n) :: xf

    final state

procedure, public, non_overridable :: integrate_to_event

integration with event finding

  • private subroutine integrate_to_event(me, t0, x0, h, tmax, tol, tf, xf, gf)

    Event-finding integration routine for the rk_class. Integrates until g(t,x)=0, or until t=tf (whichever happens first).

    Read more…

    Arguments

    Type IntentOptional Attributes Name
    class(rk_class), intent(inout) :: me
    real(kind=wp), intent(in) :: t0

    initial time

    real(kind=wp), intent(in), dimension(me%n) :: x0

    initial state

    real(kind=wp), intent(in) :: h

    abs(time step)

    real(kind=wp), intent(in) :: tmax

    max final time if event not located

    real(kind=wp), intent(in) :: tol

    function tolerance for root finding

    real(kind=wp), intent(out) :: tf

    actual final time reached

    real(kind=wp), intent(out), dimension(me%n) :: xf

    final state (at tf)

    real(kind=wp), intent(out) :: gf

    g value at tf

procedure, public :: destroy

destructor

  • private subroutine destroy(me)

    Destructor for rk_class.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_class), intent(out) :: me

procedure, public :: step => rk8_10

  • private subroutine rk8_10(me, t, x, h, xf)

    Take one Runge Kutta 8 integration step: t -> t+h (x -> xf) This is Formula (8-10) from Reference [1].

    Read more…

    Arguments

    Type IntentOptional Attributes Name
    class(rk8_10_class), intent(inout) :: me
    real(kind=wp), intent(in) :: t

    initial time

    real(kind=wp), intent(in), dimension(me%n) :: x

    initial state

    real(kind=wp), intent(in) :: h

    time step

    real(kind=wp), intent(out), dimension(me%n) :: xf

    state at time t+h

Source Code

    type,extends(rk_class),public :: rk8_10_class
        !! 8th order Runge-Kutta method.
        contains
        procedure :: step => rk8_10
    end type rk8_10_class