stepsize_class Derived Type

type, public :: stepsize_class

Algorithms for adjusting the step size for variable-step Runge-Kutta integrators.


Inherited by

type~~stepsize_class~~InheritedByGraph type~stepsize_class stepsize_class type~rk_variable_step_class rk_variable_step_class type~rk_variable_step_class->type~stepsize_class stepsize_method 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
logical, private :: fixed_step_mode = .false.

if true, then the method runs in fixed step mode with not error estimation

real(kind=wp), private :: hmax = 1.0e+6_wp

maximum allowed step size

real(kind=wp), private :: hmin = 1.0e-6_wp

minimum allowed step size

real(kind=wp), private :: hfactor_reject = 0.5_wp

minimum allowed factor for decreasing step size after rejected step

real(kind=wp), private :: hfactor_accept = 2.0_wp

maximum allowed factor for increasing step size after accepted step

integer, private :: accept_mode = 2

method to determine if step is accepted [1,2]

integer, private :: max_attempts = 10000

maximum number of attempts to decrease step size before giving up

logical, private :: relative_err = .false.

to use tol*h in the hfactor equation

real(kind=wp), private :: safety_factor = 0.9_wp

for hfactor equation (>0)

integer, private :: p_exponent_offset = 1

p + this value in the exponent (0 or 1)

procedure(norm_func), private, nopass, pointer :: norm => maxval_func

routine for computing the norm of the state


Type-Bound Procedures

procedure, public :: initialize => stepsize_class_constructor

  • private pure subroutine stepsize_class_constructor(me, hmin, hmax, hfactor_reject, hfactor_accept, norm, accept_mode, relative_err, safety_factor, p_exponent_offset, max_attempts, fixed_step_mode)

    Constructor for a stepsize_class.

    Arguments

    Type IntentOptional Attributes Name
    class(stepsize_class), intent(inout) :: me
    real(kind=wp), intent(in), optional :: hmin

    minimum allowed step size (>0)

    real(kind=wp), intent(in), optional :: hmax

    maximum allowed step size (>0)

    real(kind=wp), intent(in), optional :: hfactor_reject

    minimum allowed factor for decreasing step size after rejected step (>0)

    real(kind=wp), intent(in), optional :: hfactor_accept

    maximum allowed factor for decreasing step size after accepted step (>0)

    procedure(norm_func), optional :: norm

    the user-specified function

    integer, intent(in), optional :: accept_mode

    method to determine if step is accepted [1,2]

    logical, intent(in), optional :: relative_err

    to use tol*h in the hfactor equation

    real(kind=wp), intent(in), optional :: safety_factor

    for hfactor equation (>0)

    integer, intent(in), optional :: p_exponent_offset

    p + this value in the exponent (0 or 1)

    integer, intent(in), optional :: max_attempts

    max step size change attempts after rejected step

    logical, intent(in), optional :: fixed_step_mode

    if true, then the method runs in fixed step mode with not error estimation. All the other inputs are ignored. Note that this requires a dt /= 0 input for the integrator.

procedure, public :: compute_stepsize

  • private subroutine compute_stepsize(me, n, h, tol, err, p, hnew, accept)

    Compute the new step size using the specific method.

    Arguments

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

    number of variables

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

    current step size (<>0)

    real(kind=wp), intent(in), dimension(n) :: tol

    abs error tolerance (>0)

    real(kind=wp), intent(in), dimension(n) :: err

    truncation error estimate (>0)

    integer, intent(in) :: p

    order of the method

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

    new step size (<>0)

    logical, intent(out) :: accept

    if the step is accepted

procedure, public :: destroy => destroy_stepsize_class

Source Code

    type,public :: stepsize_class

        !! Algorithms for adjusting the step size for variable-step
        !! Runge-Kutta integrators.

        private

        logical :: fixed_step_mode = .false. !! if true, then the method runs in
                                             !! fixed step mode with not error estimation

        real(wp) :: hmax           = 1.0e+6_wp  !! maximum allowed step size
        real(wp) :: hmin           = 1.0e-6_wp  !! minimum allowed step size
        real(wp) :: hfactor_reject = 0.5_wp     !! minimum allowed factor for decreasing step size after rejected step
        real(wp) :: hfactor_accept = 2.0_wp     !! maximum allowed factor for increasing step size after accepted step
        integer  :: accept_mode    = 2          !! method to determine if step is accepted [1,2]
        integer  :: max_attempts   = 10000      !! maximum number of attempts to decrease step size before giving up

        ! for the `hfactor` equation:
        logical  :: relative_err      = .false. !! to use `tol*h` in the `hfactor` equation
        real(wp) :: safety_factor     = 0.9_wp  !! for `hfactor` equation (>0)
        integer  :: p_exponent_offset = 1       !! `p` + this value in the exponent (0 or 1)

        procedure(norm_func),nopass,pointer :: norm => maxval_func
            !! routine for computing the norm of the state

        contains

        private

        procedure,public :: initialize => stepsize_class_constructor
        procedure,public :: compute_stepsize
        procedure,public :: destroy => destroy_stepsize_class

    end type stepsize_class