rko10_class Derived Type

type, public, extends(rk_fixed_step_class) :: rko10_class

10th order Runge-Kutta Ono


Inherits

type~~rko10_class~~InheritsGraph type~rko10_class rko10_class type~rk_fixed_step_class rk_fixed_step_class type~rko10_class->type~rk_fixed_step_class type~rk_class rk_class type~rk_fixed_step_class->type~rk_class root_method root_method type~rk_class->root_method solver

Type-Bound Procedures

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 :: stop => rk_class_stop

user-callable method to stop the integration

  • private subroutine rk_class_stop(me)

    User-callable method to stop the integration.

    Arguments

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

procedure, public :: status => rk_class_status

get status code and message

  • private subroutine rk_class_status(me, istatus, message)

    Get the status of an integration.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_class), intent(in) :: me
    integer, intent(out), optional :: istatus

    status code (<0 means an error)

    character(len=:), intent(out), optional, allocatable :: message

    status message

procedure, public :: failed

  • private function failed(me)

    Returns true if there was an error. Can use rk_class_status to get more info.

    Arguments

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

    Return Value logical

procedure, public :: initialize => initialize_fixed_step

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

  • private subroutine initialize_fixed_step(me, n, f, report, g, stop_on_errors, max_number_of_steps, report_rate, solver)

    Initialize the rk_fixed_step_class.

    Arguments

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

    logical, intent(in), optional :: stop_on_errors

    stop the program for any errors (default is False)

    integer, intent(in), optional :: max_number_of_steps

    max number of steps allowed

    integer, intent(in), optional :: report_rate

    how often to call report function: 0 : no reporting (same as not associating report), 1 : report every point, 2 : report every other point, etc. The first and last point are always reported.

    class(root_method), intent(in), optional :: solver

    the root-finding method to use for even finding. if not present, then brent_solver is used.

procedure, public :: integrate => integrate_fixed_step

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

    Main integration routine for the rk_class.

    Arguments

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

    abs(time step)

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

    final time

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

    final state

procedure, public :: integrate_to_event => integrate_to_event_fixed_step

  • private subroutine integrate_to_event_fixed_step(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_fixed_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

    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(:) :: xf

    final state (at tf)

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

    g value at tf

procedure, public :: step => rko10

  • interface

    private module subroutine rko10(me, t, x, h, xf)

    Arguments

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

procedure, public :: properties => rko10_properties

  • interface

    private pure module function rko10_properties(me) result(p)

    Arguments

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

    Return Value type(rklib_properties)

    properties of the method