contd8 Function

private function contd8(me, ii, x) result(y)

this function can be used for continuous output in connection with the output-subroutine for dop853. it provides an approximation to the ii-th component of the solution at x.

Arguments

Type IntentOptional Attributes Name
class(dop853_class), intent(in) :: me
integer, intent(in) :: ii
real(kind=wp), intent(in) :: x

Return Value real(kind=wp)


Contents

Source Code


Source Code

    function contd8(me,ii,x) result(y)

    implicit none

    class(dop853_class),intent(in) :: me
    integer,intent(in)             :: ii
    real(wp),intent(in)            :: x
    real(wp)                       :: y

    real(wp) :: conpar, s, s1
    integer :: i,j,nd,ierr

    ! compute place of ii-th component
    i = 0
    do j = 1, me%nrdens
        if ( me%icomp(j)==ii ) i = j
    end do
    if ( i==0 ) then
        !always report this message, since it is an invalid use of the code.
        if (me%iprint==0) then
            ierr = error_unit
        else
            ierr = me%iprint
        end if
        write (ierr,*) &
            ' Error in contd8: no dense output available for component:', ii
        y = 0.0_wp
    else
        nd = me%nrdens
        s = (x-me%xold)/me%hout
        s1 = 1.0_wp - s
        conpar = me%cont(i+nd*4) + &
                 s*(me%cont(i+nd*5)+ &
                 s1*(me%cont(i+nd*6)+s*me%cont(i+nd*7)))
        y = me%cont(i) + &
            s*(me%cont(i+nd)+ &
            s1*(me%cont(i+nd*2)+&
            s*(me%cont(i+nd*3)+s1*conpar)))
    end if

    end function contd8