stepsize_class_constructor Subroutine

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.

Type Bound

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.


Source Code

    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)

    implicit none

    class(stepsize_class),intent(inout)       :: me
    real(wp),intent(in),optional              :: hmin             !! minimum allowed step size (>0)
    real(wp),intent(in),optional              :: hmax             !! maximum allowed step size (>0)
    real(wp),intent(in),optional              :: hfactor_reject   !! minimum allowed factor for
                                                                  !! decreasing step size after
                                                                  !! rejected step (>0)
    real(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 \( ||x|| \)
                                                                  !! function
    integer,intent(in),optional               :: accept_mode      !! method to determine if step
                                                                  !! is accepted [1,2]
    integer,intent(in),optional               :: max_attempts     !! max step size change attempts
                                                                  !! after rejected step
    logical,intent(in),optional   :: relative_err       !! to use `tol*h` in the `hfactor` equation
    real(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)
    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.

    if (present(hmin))                me%hmin                = abs(hmin)
    if (present(hmax))                me%hmax                = abs(hmax)
    if (present(hfactor_reject))      me%hfactor_reject      = abs(hfactor_reject)
    if (present(hfactor_accept))      me%hfactor_accept      = abs(hfactor_accept)
    if (present(norm))                me%norm                => norm
    if (present(accept_mode))         me%accept_mode         = accept_mode
    if (present(max_attempts))        me%max_attempts        = max_attempts

    !if (present(compute_h_factor)) me%compute_h_factor => compute_h_factor
    if (present(relative_err     )) me%relative_err      = relative_err
    if (present(safety_factor    )) me%safety_factor     = abs(safety_factor    )
    if (present(p_exponent_offset)) me%p_exponent_offset = abs(p_exponent_offset)

    if (present(fixed_step_mode)) me%fixed_step_mode = fixed_step_mode

    end subroutine stepsize_class_constructor