integrate_fixed_step Subroutine

private subroutine integrate_fixed_step(me, t0, x0, h, tf, xf)

Main integration routine for the rk_class.

Type Bound

rk_fixed_step_class

Arguments

Type IntentOptional Attributes Name
class(rk_fixed_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

abs(time step)

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

final time

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

final state


Calls

proc~~integrate_fixed_step~~CallsGraph proc~integrate_fixed_step rklib_module::rk_fixed_step_class%integrate_fixed_step proc~begin_integration_rk_fixed_step_class rklib_module::rk_fixed_step_class%begin_integration_rk_fixed_step_class proc~integrate_fixed_step->proc~begin_integration_rk_fixed_step_class proc~export_point rklib_module::rk_class%export_point proc~integrate_fixed_step->proc~export_point proc~raise_exception rklib_module::rk_class%raise_exception proc~integrate_fixed_step->proc~raise_exception step step proc~integrate_fixed_step->step begin begin proc~begin_integration_rk_fixed_step_class->begin

Source Code

    subroutine integrate_fixed_step(me,t0,x0,h,tf,xf)

    implicit none

    class(rk_fixed_step_class),intent(inout) :: me
    real(wp),intent(in)               :: t0    !! initial time
    real(wp),dimension(:),intent(in)  :: x0    !! initial state
    real(wp),intent(in)               :: h     !! abs(time step)
    real(wp),intent(in)               :: tf    !! final time
    real(wp),dimension(:),intent(out) :: xf    !! final state

    real(wp) :: t  !! current time value
    real(wp) :: dt !! time step from `t` to `t2`
    real(wp) :: t2 !! time to step to from `t`
    real(wp),dimension(me%n) :: x !! state vector
    logical :: last !! if it is the last step

    if (.not. associated(me%f)) then
        call me%raise_exception(RKLIB_ERROR_F_NOT_ASSOCIATED)
        return
    end if
    if (abs(h)<=zero) then
        call me%raise_exception(RKLIB_ERROR_INVALID_H)
        return
    end if

    call me%begin_integration()

    call me%export_point(t0,x0,.true.)  !first point

    if (abs(h)>zero) then

        t = t0
        x = x0
        dt = sign(h,tf-t0)  !time step (correct sign)
        do
            t2 = t + dt
            last = ((dt>=zero .and. t2>=tf) .or. &  !adjust last time step
                    (dt<zero .and. t2<=tf))         !
            if (last) dt = tf-t                     !
            call me%step(t,x,dt,xf)
            if (me%stopped) return
            me%num_steps = me%num_steps + 1
            if (me%num_steps > me%max_number_of_steps) then
                call me%raise_exception(RKLIB_ERROR_TOO_MANY_STEPS)
                return
            end if
            if (last) exit
            call me%export_point(t2,xf)   !intermediate point
            x = xf
            t = t2
        end do

    else
        xf = x0
    end if

    call me%export_point(tf,xf,.true.)   !last point

    end subroutine integrate_fixed_step