rk_class Derived Type

type, public :: rk_class

main integration class


Inherits

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

Inherited by

type~~rk_class~~InheritedByGraph type~rk_class rk_class type~rk_fixed_step_class rk_fixed_step_class type~rk_fixed_step_class->type~rk_class type~rk_variable_step_class rk_variable_step_class type~rk_variable_step_class->type~rk_class type~dverk65_class dverk65_class type~dverk65_class->type~rk_variable_step_class type~dverk78_class dverk78_class type~dverk78_class->type~rk_variable_step_class type~euler_class euler_class type~euler_class->type~rk_fixed_step_class type~heun_class heun_class type~heun_class->type~rk_fixed_step_class type~midpoint_class midpoint_class type~midpoint_class->type~rk_fixed_step_class type~rk3_class rk3_class type~rk3_class->type~rk_fixed_step_class type~rk4_class rk4_class type~rk4_class->type~rk_fixed_step_class type~rk5_class rk5_class type~rk5_class->type~rk_fixed_step_class type~rk7_class rk7_class type~rk7_class->type~rk_fixed_step_class type~rk8_10_class rk8_10_class type~rk8_10_class->type~rk_fixed_step_class type~rk8_12_class rk8_12_class type~rk8_12_class->type~rk_fixed_step_class type~rk_variable_step_fsal_class rk_variable_step_fsal_class type~rk_variable_step_fsal_class->type~rk_variable_step_class type~rkb109_class rkb109_class type~rkb109_class->type~rk_variable_step_class type~rkb6_class rkb6_class type~rkb6_class->type~rk_fixed_step_class type~rkbs54_class rkbs54_class type~rkbs54_class->type~rk_variable_step_class type~rkc108_class rkc108_class type~rkc108_class->type~rk_variable_step_class type~rkc5_class rkc5_class type~rkc5_class->type~rk_fixed_step_class type~rkc65_class rkc65_class type~rkc65_class->type~rk_variable_step_class type~rkck54_class rkck54_class type~rkck54_class->type~rk_variable_step_class type~rkcv8_class rkcv8_class type~rkcv8_class->type~rk_fixed_step_class type~rkdp65_class rkdp65_class type~rkdp65_class->type~rk_variable_step_class type~rkdp85_class rkdp85_class type~rkdp85_class->type~rk_variable_step_class type~rkdp87_class rkdp87_class type~rkdp87_class->type~rk_variable_step_class type~rkev87_class rkev87_class type~rkev87_class->type~rk_variable_step_class type~rkf108_class rkf108_class type~rkf108_class->type~rk_variable_step_class type~rkf1210_class rkf1210_class type~rkf1210_class->type~rk_variable_step_class type~rkf1412_class rkf1412_class type~rkf1412_class->type~rk_variable_step_class type~rkf45_class rkf45_class type~rkf45_class->type~rk_variable_step_class type~rkf78_class rkf78_class type~rkf78_class->type~rk_variable_step_class type~rkf89_class rkf89_class type~rkf89_class->type~rk_variable_step_class type~rkh10_class rkh10_class type~rkh10_class->type~rk_fixed_step_class type~rkk87_class rkk87_class type~rkk87_class->type~rk_variable_step_class type~rkl5_class rkl5_class type~rkl5_class->type~rk_fixed_step_class type~rklk5a_class rklk5a_class type~rklk5a_class->type~rk_fixed_step_class type~rklk5b_class rklk5b_class type~rklk5b_class->type~rk_fixed_step_class type~rkls44_class rkls44_class type~rkls44_class->type~rk_fixed_step_class type~rkls54_class rkls54_class type~rkls54_class->type~rk_fixed_step_class type~rko10_class rko10_class type~rko10_class->type~rk_fixed_step_class type~rko129_class rko129_class type~rko129_class->type~rk_variable_step_class type~rkr4_class rkr4_class type~rkr4_class->type~rk_fixed_step_class type~rks1110a_class rks1110a_class type~rks1110a_class->type~rk_variable_step_class type~rks4_class rks4_class type~rks4_class->type~rk_fixed_step_class type~rks5_class rks5_class type~rks5_class->type~rk_fixed_step_class type~rks98_class rks98_class type~rks98_class->type~rk_variable_step_class type~rkss54_class rkss54_class type~rkss54_class->type~rk_variable_step_class type~rkss76_class rkss76_class type~rkss76_class->type~rk_variable_step_class type~rkssp22_class rkssp22_class type~rkssp22_class->type~rk_fixed_step_class type~rkssp33_class rkssp33_class type~rkssp33_class->type~rk_fixed_step_class type~rkssp43_class rkssp43_class type~rkssp43_class->type~rk_variable_step_class type~rkssp53_class rkssp53_class type~rkssp53_class->type~rk_fixed_step_class type~rkssp54_class rkssp54_class type~rkssp54_class->type~rk_fixed_step_class type~rkt98a_class rkt98a_class type~rkt98a_class->type~rk_variable_step_class type~rktmy7_class rktmy7_class type~rktmy7_class->type~rk_variable_step_class type~rktmy7s_class rktmy7s_class type~rktmy7s_class->type~rk_variable_step_class type~rktp64_class rktp64_class type~rktp64_class->type~rk_variable_step_class type~rktp75_class rktp75_class type~rktp75_class->type~rk_variable_step_class type~rktp86_class rktp86_class type~rktp86_class->type~rk_variable_step_class type~rkv65_class rkv65_class type~rkv65_class->type~rk_variable_step_class type~rkv76e_class rkv76e_class type~rkv76e_class->type~rk_variable_step_class type~rkv76r_class rkv76r_class type~rkv76r_class->type~rk_variable_step_class type~rkv78_class rkv78_class type~rkv78_class->type~rk_variable_step_class type~rkv87e_class rkv87e_class type~rkv87e_class->type~rk_variable_step_class type~rkv87r_class rkv87r_class type~rkv87r_class->type~rk_variable_step_class type~rkv89_class rkv89_class type~rkv89_class->type~rk_variable_step_class type~rkv98e_class rkv98e_class type~rkv98e_class->type~rk_variable_step_class type~rkv98r_class rkv98r_class type~rkv98r_class->type~rk_variable_step_class type~rkz10_class rkz10_class type~rkz10_class->type~rk_fixed_step_class type~rkbs32_class rkbs32_class type~rkbs32_class->type~rk_variable_step_fsal_class type~rkdp54_class rkdp54_class type~rkdp54_class->type~rk_variable_step_fsal_class type~rkpp54_class rkpp54_class type~rkpp54_class->type~rk_variable_step_fsal_class type~rkpp54b_class rkpp54b_class type~rkpp54b_class->type~rk_variable_step_fsal_class type~rks54_class rks54_class type~rks54_class->type~rk_variable_step_fsal_class type~rkt54_class rkt54_class type~rkt54_class->type~rk_variable_step_fsal_class type~rktf65_class rktf65_class type~rktf65_class->type~rk_variable_step_fsal_class type~rkv65e_class rkv65e_class type~rkv65e_class->type~rk_variable_step_fsal_class type~rkv65r_class rkv65r_class type~rkv65r_class->type~rk_variable_step_fsal_class

Components

Type Visibility Attributes Name Initial
integer, private :: istatus = 0

status code

logical, private :: stopped = .false.

if user has stopped the integration in f or report.

integer, private :: num_steps = 0

number of accepted steps taken

integer, private :: max_number_of_steps = huge(1)

maximum number of steps to take

integer, private :: report_rate = 1

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.

logical, private :: stop_on_errors = .false.

if true, then errors will stop the program

integer, private :: n = 0

user specified number of variables

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

user-specified derivative function

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

user-specified report function

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

event function (stop when this is zero)

type(root_method), private :: solver = root_method_brent

the root solver method to use for even finding

real(kind=wp), private, dimension(:,:), allocatable :: funcs

matrix to store the function evalutaions in the step function. this will be size (n x number_of_registers)


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, private :: init => initialize_rk_class

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

    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

    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, private :: begin => begin_integration_rk_class

procedure, private :: raise_exception

  • private subroutine raise_exception(me, error_code)

    Raise an exception.

    Arguments

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

    the error to raise

procedure, private :: clear_exception

  • private subroutine clear_exception(me)

    Clear any exception.

    Arguments

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

procedure, private :: export_point

  • private subroutine export_point(me, t, x, first_or_last)

    Wrapper for exporting points during integration.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_class), intent(inout) :: me
    real(kind=wp), intent(in) :: t
    real(kind=wp), intent(in), dimension(:) :: x
    logical, intent(in), optional :: first_or_last

    if this is the first or last point (always reported)

procedure(begin_func), private, deferred :: begin_integration

  • subroutine begin_func(me) Prototype

    routine called before integration begins to set up internal variables.

    Arguments

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

procedure(properties_func), public, deferred :: properties

  • pure function properties_func(me) result(p) Prototype

    Returns the properties of the method.

    Arguments

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

    Return Value type(rklib_properties)

    properties of the method

Source Code

    type,abstract,public :: rk_class

        !! main integration class

        private

        integer :: istatus = 0 !! status code
        logical :: stopped = .false. !! if user has stopped the integration in `f` or `report`.
        integer :: num_steps = 0 !! number of accepted steps taken
        integer :: max_number_of_steps  = huge(1) !! maximum number of steps to take
        integer :: report_rate = 1 !! 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.

        logical :: stop_on_errors = .false. !! if true, then errors will stop the program
        integer :: n = 0  !! user specified number of variables
        procedure(deriv_func),pointer  :: f      => null()  !! user-specified derivative function
        procedure(report_func),pointer :: report => null()  !! user-specified report function
        procedure(event_func),pointer  :: g      => null()  !! event function (stop when this is zero)
        type(root_method) :: solver = root_method_brent !! the root solver method to use for even finding

        real(wp),dimension(:,:),allocatable :: funcs !! matrix to store the function
                                                     !! evalutaions in the step function.
                                                     !! this will be size (`n` x `number_of_registers`)

        contains

        private

        procedure,public :: destroy !! destructor
        procedure,public :: stop => rk_class_stop !! user-callable method to stop the integration
        procedure,public :: status => rk_class_status !! get status code and message
        procedure,public :: failed

        procedure :: init => initialize_rk_class
        procedure :: begin => begin_integration_rk_class
        procedure :: raise_exception
        procedure :: clear_exception
        procedure :: export_point
        procedure(begin_func),deferred :: begin_integration
        procedure(properties_func),deferred,public :: properties

    end type rk_class