rk_variable_step_fsal_class Derived Type

type, public, extends(rk_variable_step_class) :: rk_variable_step_fsal_class

a variable step method with the "first same as last" (FSAL) property. Cache the last f and x vectors to use for the next step.

The assumption is that the nature of the function has not changed since the last step. If it has, the user would need to manually call destroy_fsal_cache so that the previous point was not reused.


Inherits

type~~rk_variable_step_fsal_class~~InheritsGraph type~rk_variable_step_fsal_class rk_variable_step_fsal_class type~rk_variable_step_class rk_variable_step_class type~rk_variable_step_fsal_class->type~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_fsal_class~~InheritedByGraph type~rk_variable_step_fsal_class rk_variable_step_fsal_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
real(kind=wp), private, allocatable :: t_saved

cached t

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

cached x

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

cached f


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, 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, public :: destroy_fsal_cache

procedure, public :: check_fsal_cache

  • private subroutine check_fsal_cache(me, t, x, f)

    Check the FSAL cache.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_variable_step_fsal_class), intent(inout) :: me
    real(kind=wp), intent(in) :: t
    real(kind=wp), intent(in), dimension(:) :: x
    real(kind=wp), intent(out), dimension(:) :: f

procedure, public :: set_fsal_cache

  • private subroutine set_fsal_cache(me, t, x, f)

    Compute the function and add it to the FSAL cache.

    Arguments

    Type IntentOptional Attributes Name
    class(rk_variable_step_fsal_class), intent(inout) :: me
    real(kind=wp), intent(in) :: t
    real(kind=wp), intent(in), dimension(:) :: x
    real(kind=wp), intent(out), dimension(:) :: f

Source Code

    type,extends(rk_variable_step_class),abstract,public :: rk_variable_step_fsal_class
        !! a variable step method with the "first same as last" (FSAL) property.
        !! Cache the last `f` and `x` vectors to use for the next step.
        !!
        !! The assumption is that the nature of the
        !! function has not changed since the last step.
        !! If it has, the user would need to manually call [[destroy_fsal_cache]]
        !! so that the previous point was not reused.
        private
        real(wp),allocatable :: t_saved !! cached `t`
        real(wp),dimension(:),allocatable :: x_saved  !! cached `x`
        real(wp),dimension(:),allocatable :: f_saved !! cached `f`
        contains
        private
        procedure,public :: destroy_fsal_cache
        procedure,public :: check_fsal_cache
        procedure,public :: set_fsal_cache
    end type rk_variable_step_fsal_class