rkf1210_class Derived Type

type, public, extends(rk_variable_step_class) :: rkf1210_class

Runga-Kutta Feagin 12(10) method.


Inherits

type~~rkf1210_class~~InheritsGraph type~rkf1210_class rkf1210_class type~rk_variable_step_class rk_variable_step_class type~rkf1210_class->type~rk_variable_step_class type~stepsize_class stepsize_class type~rk_variable_step_class->type~stepsize_class stepsize_method

Type-Bound Procedures

procedure, public :: initialize

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

  • private subroutine initialize(me, n, f, rtol, atol, stepsize_method, hinit_method, report, g)

    Initialize the rk_variable_step_class.

    Arguments

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

    number of equations

    procedure(deriv_func) :: f

    derivative function

    real(kind=wp), intent(in), dimension(:) :: rtol

    relative tolerance (if size=1, then same tol used for all equations)

    real(kind=wp), intent(in), dimension(:) :: atol

    absolute tolerance (if size=1, then same tol used for all equations)

    class(stepsize_class), intent(in) :: stepsize_method

    method for varying the step size

    integer, intent(in), optional :: hinit_method

    which method to use for automatic initial step size computation. 1 = use hstart, 2 = use hinit.

    procedure(report_func), optional :: report

    for reporting the steps

    procedure(event_func), optional :: g

    for stopping at an event

procedure, public :: destroy

destructor

procedure, public, non_overridable :: integrate

main integration routine

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

    Main integration routine for the rk_variable_step_class.

    Arguments

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

    initial time

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

    initial state

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

    initial abs(time step)

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

    final time

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

    final state

    integer, intent(out), optional :: ierr

    0 = no errors, <0 = error. if not present, an error will stop program.

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, ierr)

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

    Read more…

    Arguments

    Type IntentOptional Attributes Name
    class(rk_variable_step_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

    integer, intent(out), optional :: ierr

    0 = no errors, <0 = error. if not present, an error will stop program.

procedure, public :: step => rkf1210

  • private subroutine rkf1210(me, t, x, h, xf, terr)

    Feagin's RK12(10) method -- a 12th-order method with an embedded 10th-order method.

    Read more…

    Arguments

    Type IntentOptional Attributes Name
    class(rkf1210_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

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

    truncation error estimate

procedure, public :: order => rkf1210_order

  • private pure function rkf1210_order(me) result(p)

    Returns the order of the rkf1210 method.

    Arguments

    Type IntentOptional Attributes Name
    class(rkf1210_class), intent(in) :: me

    Return Value integer

    order of the method

Source Code

    type,extends(rk_variable_step_class),public :: rkf1210_class
        !! Runga-Kutta Feagin 12(10) method.
        contains
        procedure :: step  => rkf1210
        procedure :: order => rkf1210_order
    end type rkf1210_class