dintrv Subroutine

private pure subroutine dintrv(xt, lxt, xx, ilo, ileft, mflag, extrap)

Computes the largest integer ileft in 1 ileft lxt such that xt(ileft) x where xt(*) is a subdivision of the x interval. precisely,

         if            x < xt(1)   then ileft=1,   mflag=-1
         if   xt(i) <= x < xt(i+1) then ileft=i,   mflag=0
         if xt(lxt) <= x           then ileft=lxt, mflag=-2

that is, when multiplicities are present in the break point to the left of x, the largest index is taken for ileft.

History

  • interv written by carl de boor [5]
  • dintrv author: amos, d. e., (snla) : date written 800901
  • revision date 820801
  • Jacob Williams, 2/24/2015 : updated to free-form Fortran.
  • Jacob Williams, 2/17/2016 : additional refactoring (eliminated GOTOs).
  • Jacob Williams, 3/4/2017 : added extrapolation option.

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: xt

a knot or break point vector of length lxt

integer(kind=ip), intent(in) :: lxt

length of the xt vector

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

argument

integer(kind=ip), intent(inout) :: ilo

an initialization parameter which must be set to 1 the first time the spline array xt is processed by dintrv. ilo contains information for efficient processing after the initial call and ilo must not be changed by the user. distinct splines require distinct ilo parameters.

integer(kind=ip), intent(out) :: ileft

largest integer satisfying xt(ileft) x

integer(kind=ip), intent(out) :: mflag

signals when x lies out of bounds

logical, intent(in), optional :: extrap

if extrapolation is allowed (if not present, default is False)


Calls

proc~~dintrv~~CallsGraph proc~dintrv bspline_sub_module::dintrv proc~get_temp_x_for_extrap bspline_sub_module::get_temp_x_for_extrap proc~dintrv->proc~get_temp_x_for_extrap

Called by

proc~~dintrv~~CalledByGraph proc~dintrv bspline_sub_module::dintrv proc~db2val bspline_sub_module::db2val proc~db2val->proc~dintrv proc~dbvalu bspline_sub_module::dbvalu proc~db2val->proc~dbvalu proc~db3val bspline_sub_module::db3val proc~db3val->proc~dintrv proc~db3val->proc~dbvalu proc~db4val bspline_sub_module::db4val proc~db4val->proc~dintrv proc~db4val->proc~dbvalu proc~db5val bspline_sub_module::db5val proc~db5val->proc~dintrv proc~db5val->proc~dbvalu proc~db6val bspline_sub_module::db6val proc~db6val->proc~dintrv proc~db6val->proc~dbvalu proc~dbfqad bspline_sub_module::dbfqad proc~dbfqad->proc~dintrv proc~dbsgq8 bspline_sub_module::dbsgq8 proc~dbfqad->proc~dbsgq8 proc~dbsqad bspline_sub_module::dbsqad proc~dbsqad->proc~dintrv proc~dbsqad->proc~dbvalu proc~dbvalu->proc~dintrv proc~db1fqad bspline_sub_module::db1fqad proc~db1fqad->proc~dbfqad proc~db1sqad bspline_sub_module::db1sqad proc~db1sqad->proc~dbsqad proc~db1val_alt bspline_sub_module::db1val_alt proc~db1val_alt->proc~dbvalu proc~db1val_default bspline_sub_module::db1val_default proc~db1val_default->proc~dbvalu proc~dbsgq8->proc~dbvalu proc~evaluate_2d bspline_oo_module::bspline_2d%evaluate_2d proc~evaluate_2d->proc~db2val proc~evaluate_3d bspline_oo_module::bspline_3d%evaluate_3d proc~evaluate_3d->proc~db3val proc~evaluate_4d bspline_oo_module::bspline_4d%evaluate_4d proc~evaluate_4d->proc~db4val proc~evaluate_5d bspline_oo_module::bspline_5d%evaluate_5d proc~evaluate_5d->proc~db5val proc~evaluate_6d bspline_oo_module::bspline_6d%evaluate_6d proc~evaluate_6d->proc~db6val interface~db1val bspline_sub_module::db1val interface~db1val->proc~db1val_alt interface~db1val->proc~db1val_default proc~fintegral_1d bspline_oo_module::bspline_1d%fintegral_1d proc~fintegral_1d->proc~db1fqad proc~integral_1d bspline_oo_module::bspline_1d%integral_1d proc~integral_1d->proc~db1sqad proc~evaluate_1d bspline_oo_module::bspline_1d%evaluate_1d proc~evaluate_1d->interface~db1val

Source Code

    pure subroutine dintrv(xt,lxt,xx,ilo,ileft,mflag,extrap)

    implicit none

    integer(ip),intent(in)             :: lxt    !! length of the `xt` vector
    real(wp),dimension(:),intent(in)   :: xt     !! a knot or break point vector of length `lxt`
    real(wp),intent(in)                :: xx     !! argument
    integer(ip),intent(inout)          :: ilo    !! an initialization parameter which must be set
                                                 !! to 1 the first time the spline array `xt` is
                                                 !! processed by dintrv. `ilo` contains information for
                                                 !! efficient processing after the initial call and `ilo`
                                                 !! must not be changed by the user.  distinct splines
                                                 !! require distinct `ilo` parameters.
    integer(ip),intent(out)            :: ileft  !! largest integer satisfying `xt(ileft)` \( \le \) `x`
    integer(ip),intent(out)            :: mflag  !! signals when `x` lies out of bounds
    logical,intent(in),optional        :: extrap !! if extrapolation is allowed
                                                 !! (if not present, default is False)

    integer(ip) :: ihi, istep, middle
    real(wp) :: x

    x = get_temp_x_for_extrap(xx,xt(1_ip),xt(lxt),extrap)

    ihi = ilo + 1_ip
    if ( ihi>=lxt ) then
        if ( x>=xt(lxt) ) then
            mflag = -2_ip
            ileft = lxt
            return
        end if
        if ( lxt<=1 ) then
            mflag = -1_ip
            ileft = 1_ip
            return
        end if
        ilo = lxt - 1_ip
        ihi = lxt
    end if

    if ( x>=xt(ihi) ) then

        ! now x >= xt(ilo). find upper bound
        istep = 1_ip
        do
            ilo = ihi
            ihi = ilo + istep
            if ( ihi>=lxt ) then
                if ( x>=xt(lxt) ) then
                    mflag = -2_ip
                    ileft = lxt
                    return
                end if
                ihi = lxt
            else if ( x>=xt(ihi) ) then
                istep = istep*2_ip
                cycle
            end if
            exit
        end do

    else

        if ( x>=xt(ilo) ) then
            mflag = 0_ip
            ileft = ilo
            return
        end if
        ! now x <= xt(ihi). find lower bound
        istep = 1_ip
        do
            ihi = ilo
            ilo = ihi - istep
            if ( ilo<=1_ip ) then
                ilo = 1_ip
                if ( x<xt(1_ip) ) then
                    mflag = -1_ip
                    ileft = 1_ip
                    return
                end if
            else if ( x<xt(ilo) ) then
                istep = istep*2_ip
                cycle
            end if
            exit
        end do

    end if

    ! now xt(ilo) <= x < xt(ihi). narrow the interval
    do
        middle = (ilo+ihi)/2_ip
        if ( middle==ilo ) then
            mflag = 0_ip
            ileft = ilo
            return
        end if
        ! note. it is assumed that middle = ilo in case ihi = ilo+1
        if ( x<xt(middle) ) then
            ihi = middle
        else
            ilo = middle
        end if
    end do

    end subroutine dintrv