compute_initial_step Function

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

Compute the initial step size.

Type Bound

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


Calls

proc~~compute_initial_step~~CallsGraph proc~compute_initial_step rklib_module::rk_variable_step_class%compute_initial_step f f proc~compute_initial_step->f proc~hinit rklib_module::rk_variable_step_class%hinit proc~compute_initial_step->proc~hinit proc~hstart rklib_module::rk_variable_step_class%hstart proc~compute_initial_step->proc~hstart proc~hinit->f proc~order rklib_module::rk_variable_step_class%order proc~hinit->proc~order proc~hstart->f proc~hstart->proc~order properties properties proc~order->properties

Called by

proc~~compute_initial_step~~CalledByGraph proc~compute_initial_step rklib_module::rk_variable_step_class%compute_initial_step proc~integrate_to_event_variable_step rklib_module::rk_variable_step_class%integrate_to_event_variable_step proc~integrate_to_event_variable_step->proc~compute_initial_step proc~integrate_variable_step rklib_module::rk_variable_step_class%integrate_variable_step proc~integrate_variable_step->proc~compute_initial_step program~rklib_example rklib_example program~rklib_example->proc~integrate_variable_step

Source Code

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

        class(rk_variable_step_class),intent(inout) :: me
        real(wp),intent(in) :: h0 !! user-input initial step size (if zero, then one is computed)
        real(wp),intent(in) :: t0 !! initial time
        real(wp),intent(in) :: tf !! final time
        real(wp) :: dt !! step size to use

        real(wp),dimension(me%n) :: x0 !! initial state
        real(wp),dimension(me%n) :: etol !! tolerance vector
        real(wp),dimension(me%n) :: f0 !! initial derivative

        if (abs(h0)<=zero) then
            ! compute an appropriate initial step size:
            etol = me%rtol * me%stepsize_method%norm(x0) + me%atol
            call me%f(t0,x0,f0)  ! get initial dx/dt
            select case (me%hinit_method) ! value was checked in initialize_variable_step
            case(1); call me%hstart(t0,tf,x0,f0,etol,dt)
            case(2); dt = me%hinit(t0,x0,sign(1.0_wp,tf-t0),f0,&
                                me%stepsize_method%hmax,&
                                me%atol,me%rtol)
            end select
        else
            ! user-specified initial step size:
            dt = sign(h0,tf-t0)  ! (correct sign)
        end if

    end function compute_initial_step