rk_variable_step_class Derived Type

type, public, extends(rk_class) :: rk_variable_step_class

Main integration class for variable step size Runge-Kutta methods


Inherits

type~~rk_variable_step_class~~InheritsGraph type~rk_variable_step_class rk_variable_step_class type~rk_class rk_class type~rk_variable_step_class->type~rk_class type~stepsize_class stepsize_class type~rk_variable_step_class->type~stepsize_class stepsize_method root_method root_method type~rk_class->root_method solver

Inherited by

type~~rk_variable_step_class~~InheritedByGraph type~rk_variable_step_class rk_variable_step_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~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~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~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~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~rkk87_class rkk87_class type~rkk87_class->type~rk_variable_step_class type~rko129_class rko129_class type~rko129_class->type~rk_variable_step_class type~rks1110a_class rks1110a_class type~rks1110a_class->type~rk_variable_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~rkssp43_class rkssp43_class type~rkssp43_class->type~rk_variable_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~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
type(stepsize_class), private :: stepsize_method

the method for varying the step size

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

relative tolerance (size(n))

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

absolute tolerance (size(n))

integer, private :: hinit_method = 1

if automatically computing the inital step size, which method to use. 1 = hstart, 2 = hinit.

integer, private :: num_rejected_steps = 0

number of rejected steps

real(kind=wp), private :: last_accepted_step_size = zero

the last accepted step size dt from the integration (positive or negative)


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(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

procedure(step_func_variable), private, deferred :: step

the step routine for the rk method

  • subroutine step_func_variable(me, t, x, h, xf, xerr) Prototype

    rk step function for the variable-step methods.

    Arguments

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

    initial time

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

    initial state vector

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

    time step

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

    final state vector

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

    truncation error estimate

procedure, public :: initialize => initialize_variable_step

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

  • private subroutine initialize_variable_step(me, n, f, rtol, atol, stepsize_method, hinit_method, report, g, stop_on_errors, max_number_of_steps, report_rate, solver)

    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), optional, dimension(:) :: rtol

    relative tolerance (if size=1, then same tol used for all equations). If not present, a default of 100*eps is used

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

    absolute tolerance (if size=1, then same tol used for all equations). If not present, a default of 100*eps is used

    type(stepsize_class), intent(in), optional :: stepsize_method

    method for varying the step size

    integer, intent(in), optional :: hinit_method

    which method (1 or 2) 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

    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_variable_step

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

    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

procedure, public :: integrate_to_event => integrate_to_event_variable_step

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

    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

procedure, public :: info => info_variable_step

  • private subroutine info_variable_step(me, num_steps, num_rejected_steps, last_accepted_step_size)

    Return some info about the integration.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_variable_step_class), intent(in) :: me
    integer, intent(out), optional :: num_steps

    number of steps taken

    integer, intent(out), optional :: num_rejected_steps

    number of rejected steps

    real(kind=wp), intent(out), optional :: last_accepted_step_size

    the last accepted step size dt from the integration (positive or negative)

procedure, private :: hstart

for automatically computing the initial step size [this is from DDEABM]

  • private subroutine hstart(me, a, b, y, yprime, etol, h)

    Computes a starting step size to be used in solving initial value problems in ordinary differential equations.

    Read more…

    Arguments

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

    the initial point of integration.

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

    a value of the independent variable used to define the direction of integration. a reasonable choice is to set b to the first point at which a solution is desired. you can also use b, if necessary, to restrict the length of the first integration step because the algorithm will not compute a starting step length which is bigger than abs(b-a), unless b has been chosen too close to a. (it is presumed that hstart has been called with b different from a on the machine being used. also see the discussion about the parameter small.)

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

    the vector of initial values of the neq solution components at the initial point a.

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

    the vector of derivatives of the neq solution components at the initial point a. (defined by the differential equations in subroutine me%f)

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

    the vector of error tolerances corresponding to the neq solution components. it is assumed that all elements are positive. following the first integration step, the tolerances are expected to be used by the integrator in an error test which roughly requires that abs(local error) <= etol for each vector component.

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

    appropriate starting step size to be attempted by the differential equation method.

procedure, private :: hinit

for automatically computing the initial step size [this is from DOP853]

  • private function hinit(me, x, y, posneg, f0, hmax, atol, rtol)

    computation of an initial step size guess

    Read more…

    Arguments

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

    dimension(n)

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

    posneg = sign(1.0_wp,xend-x)

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

    dimension(n)

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

    Return Value real(kind=wp)

procedure, private :: begin_integration => begin_integration_rk_variable_step_class

procedure, private :: compute_initial_step

  • private function compute_initial_step(me, t0, tf, x0, h0) result(dt)

    Compute the initial step size.

    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) :: tf

    final time

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

    initial state

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

    user-input initial step size (if zero, then one is computed)

    Return Value real(kind=wp)

    step size to use

procedure, private :: order

returns p, the order of the method

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

    Returns the order of the RK method

    Arguments

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

    Return Value integer

    order of the method

Source Code

    type,extends(rk_class),abstract,public :: rk_variable_step_class

        !! Main integration class for variable step size Runge-Kutta methods

        private

        type(stepsize_class) :: stepsize_method  !! the method for varying the step size

        real(wp),dimension(:),allocatable :: rtol  !! relative tolerance (`size(n)`)
        real(wp),dimension(:),allocatable :: atol  !! absolute tolerance (`size(n)`)

        integer :: hinit_method = 1 !! if automatically computing the inital step size, which
                                    !! method to use. 1 = `hstart`, 2 = `hinit`.

        integer :: num_rejected_steps = 0 !! number of rejected steps
        real(wp) :: last_accepted_step_size = zero !! the last accepted step size `dt` from the integration
                                                   !! (positive or negative)

        contains

        private

        procedure(step_func_variable),deferred :: step !! the step routine for the rk method
        procedure,public :: initialize => initialize_variable_step  !! initialize the class (set n,f, and report)
        procedure,public :: integrate => integrate_variable_step
        procedure,public :: integrate_to_event => integrate_to_event_variable_step
        procedure,public :: info => info_variable_step

        procedure :: hstart  !! for automatically computing the initial step size [this is from DDEABM]
        procedure :: hinit   !! for automatically computing the initial step size [this is from DOP853]
        procedure :: begin_integration => begin_integration_rk_variable_step_class
        procedure :: compute_initial_step
        procedure :: order !! returns `p`, the order of the method

    end type rk_variable_step_class