Processing math: 100%

db5val Subroutine

public pure subroutine db5val(xval, yval, zval, qval, rval, idx, idy, idz, idq, idr, tx, ty, tz, tq, tr, nx, ny, nz, nq, nr, kx, ky, kz, kq, kr, bcoef, f, iflag, inbvx, inbvy, inbvz, inbvq, inbvr, iloy, iloz, iloq, ilor, w4, w3, w2, w1, w0, extrap)

Evaluates the tensor product piecewise polynomial interpolant constructed by the routine db5ink or one of its derivatives at the point (xval,yval,zval,qval,rval).

To evaluate the interpolant itself, set idx=idy=idz=idq=idr=0, to evaluate the first partial with respect to x, set idx=1,idy=idz=idq=idr=0, and so on.

See db3val header for more information.

History

  • Jacob Williams, 2/24/2015 : Created this routine.

Arguments

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

x coordinate of evaluation point.

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

y coordinate of evaluation point.

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

z coordinate of evaluation point.

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

q coordinate of evaluation point.

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

r coordinate of evaluation point.

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

x derivative of piecewise polynomial to evaluate.

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

y derivative of piecewise polynomial to evaluate.

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

z derivative of piecewise polynomial to evaluate.

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

q derivative of piecewise polynomial to evaluate.

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

r derivative of piecewise polynomial to evaluate.

real(kind=wp), intent(in), dimension(nx+kx) :: tx

sequence of knots defining the piecewise polynomial in the x direction. (same as in last call to db5ink)

real(kind=wp), intent(in), dimension(ny+ky) :: ty

sequence of knots defining the piecewise polynomial in the y direction. (same as in last call to db5ink)

real(kind=wp), intent(in), dimension(nz+kz) :: tz

sequence of knots defining the piecewise polynomial in the z direction. (same as in last call to db5ink)

real(kind=wp), intent(in), dimension(nq+kq) :: tq

sequence of knots defining the piecewise polynomial in the q direction. (same as in last call to db5ink)

real(kind=wp), intent(in), dimension(nr+kr) :: tr

sequence of knots defining the piecewise polynomial in the r direction. (same as in last call to db5ink)

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

the number of interpolation points in x. (same as in last call to db5ink)

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

the number of interpolation points in y. (same as in last call to db5ink)

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

the number of interpolation points in z. (same as in last call to db5ink)

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

the number of interpolation points in q. (same as in last call to db5ink)

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

the number of interpolation points in r. (same as in last call to db5ink)

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

order of polynomial pieces in x. (same as in last call to db5ink)

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

order of polynomial pieces in y. (same as in last call to db5ink)

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

order of polynomial pieces in z. (same as in last call to db5ink)

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

order of polynomial pieces in q. (same as in last call to db5ink)

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

order of polynomial pieces in r. (same as in last call to db5ink)

real(kind=wp), intent(in), dimension(nx,ny,nz,nq,nr) :: bcoef

the b-spline coefficients computed by db5ink.

real(kind=wp), intent(out) :: f

interpolated value

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

status flag:

  • =0 : no errors
  • 0 : error
integer(kind=ip), intent(inout) :: inbvx

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

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

initialization parameter which must be set to 1 the first time this routine is called, and must not be changed by the user.

real(kind=wp), intent(inout), dimension(ky,kz,kq,kr) :: w4

work array

real(kind=wp), intent(inout), dimension(kz,kq,kr) :: w3

work array

real(kind=wp), intent(inout), dimension(kq,kr) :: w2

work array

real(kind=wp), intent(inout), dimension(kr) :: w1

work array

real(kind=wp), intent(inout), dimension(3_ip*max(kx,ky,kz,kq,kr)) :: w0

work array

logical, intent(in), optional :: extrap

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


Calls

proc~~db5val~~CallsGraph proc~db5val db5val proc~check_value check_value proc~db5val->proc~check_value proc~dbvalu dbvalu proc~db5val->proc~dbvalu proc~dintrv dintrv proc~db5val->proc~dintrv proc~dbvalu->proc~dintrv proc~get_temp_x_for_extrap get_temp_x_for_extrap proc~dintrv->proc~get_temp_x_for_extrap

Called by

proc~~db5val~~CalledByGraph proc~db5val db5val proc~evaluate_5d bspline_5d%evaluate_5d proc~evaluate_5d->proc~db5val

Source Code

    pure subroutine db5val(xval,yval,zval,qval,rval,&
                           idx,idy,idz,idq,idr,&
                           tx,ty,tz,tq,tr,&
                           nx,ny,nz,nq,nr,&
                           kx,ky,kz,kq,kr,&
                           bcoef,f,iflag,&
                           inbvx,inbvy,inbvz,inbvq,inbvr,&
                           iloy,iloz,iloq,ilor,&
                           w4,w3,w2,w1,w0,extrap)

    implicit none

    integer(ip),intent(in)                        :: idx      !! \(x\) derivative of piecewise polynomial to evaluate.
    integer(ip),intent(in)                        :: idy      !! \(y\) derivative of piecewise polynomial to evaluate.
    integer(ip),intent(in)                        :: idz      !! \(z\) derivative of piecewise polynomial to evaluate.
    integer(ip),intent(in)                        :: idq      !! \(q\) derivative of piecewise polynomial to evaluate.
    integer(ip),intent(in)                        :: idr      !! \(r\) derivative of piecewise polynomial to evaluate.
    integer(ip),intent(in)                        :: nx       !! the number of interpolation points in \(x\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: ny       !! the number of interpolation points in \(y\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: nz       !! the number of interpolation points in \(z\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: nq       !! the number of interpolation points in \(q\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: nr       !! the number of interpolation points in \(r\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: kx       !! order of polynomial pieces in \(x\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: ky       !! order of polynomial pieces in \(y\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: kz       !! order of polynomial pieces in \(z\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: kq       !! order of polynomial pieces in \(q\).
                                                              !! (same as in last call to [[db5ink]])
    integer(ip),intent(in)                        :: kr       !! order of polynomial pieces in \(r\).
                                                              !! (same as in last call to [[db5ink]])
    real(wp),intent(in)                           :: xval     !! \(x\) coordinate of evaluation point.
    real(wp),intent(in)                           :: yval     !! \(y\) coordinate of evaluation point.
    real(wp),intent(in)                           :: zval     !! \(z\) coordinate of evaluation point.
    real(wp),intent(in)                           :: qval     !! \(q\) coordinate of evaluation point.
    real(wp),intent(in)                           :: rval     !! \(r\) coordinate of evaluation point.
    real(wp),dimension(nx+kx),intent(in)          :: tx       !! sequence of knots defining the piecewise polynomial
                                                              !! in the \(x\) direction.
                                                              !! (same as in last call to [[db5ink]])
    real(wp),dimension(ny+ky),intent(in)          :: ty       !! sequence of knots defining the piecewise polynomial
                                                              !! in the \(y\) direction.
                                                              !! (same as in last call to [[db5ink]])
    real(wp),dimension(nz+kz),intent(in)          :: tz       !! sequence of knots defining the piecewise polynomial
                                                              !! in the \(z\) direction.
                                                              !! (same as in last call to [[db5ink]])
    real(wp),dimension(nq+kq),intent(in)          :: tq       !! sequence of knots defining the piecewise polynomial
                                                              !! in the \(q\) direction.
                                                              !! (same as in last call to [[db5ink]])
    real(wp),dimension(nr+kr),intent(in)          :: tr       !! sequence of knots defining the piecewise polynomial
                                                              !! in the \(r\) direction.
                                                              !! (same as in last call to [[db5ink]])
    real(wp),dimension(nx,ny,nz,nq,nr),intent(in) :: bcoef    !! the b-spline coefficients computed by [[db5ink]].
    real(wp),intent(out)                          :: f        !! interpolated value
    integer(ip),intent(out)                       :: iflag    !! status flag:
                                                              !!
                                                              !! * \( = 0 \)   : no errors
                                                              !! * \( \ne 0 \) : error
    integer(ip),intent(inout)                     :: inbvx    !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: inbvy    !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: inbvz    !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: inbvq    !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: inbvr    !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: iloy     !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: iloz     !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: iloq     !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    integer(ip),intent(inout)                     :: ilor     !! initialization parameter which must be set
                                                              !! to 1 the first time this routine is called,
                                                              !! and must not be changed by the user.
    real(wp),dimension(ky,kz,kq,kr),intent(inout)           :: w4  !! work array
    real(wp),dimension(kz,kq,kr),intent(inout)              :: w3  !! work array
    real(wp),dimension(kq,kr),intent(inout)                 :: w2  !! work array
    real(wp),dimension(kr),intent(inout)                    :: w1  !! work array
    real(wp),dimension(3_ip*max(kx,ky,kz,kq,kr)),intent(inout) :: w0  !! work array
    logical,intent(in),optional                   :: extrap   !! if extrapolation is allowed
                                                              !! (if not present, default is False)

    integer(ip) :: lefty, leftz, leftq, leftr, &
               kcoly, kcolz, kcolq, kcolr, j, k, q, r

    f = 0.0_wp

    iflag = check_value(xval,tx,1_ip,extrap); if (iflag/=0_ip) return
    iflag = check_value(yval,ty,2_ip,extrap); if (iflag/=0_ip) return
    iflag = check_value(zval,tz,3_ip,extrap); if (iflag/=0_ip) return
    iflag = check_value(qval,tq,4_ip,extrap); if (iflag/=0_ip) return
    iflag = check_value(rval,tr,5_ip,extrap); if (iflag/=0_ip) return

    call dintrv(ty,ny+ky,yval,iloy,lefty,iflag,extrap); if (iflag/=0_ip) return
    call dintrv(tz,nz+kz,zval,iloz,leftz,iflag,extrap); if (iflag/=0_ip) return
    call dintrv(tq,nq+kq,qval,iloq,leftq,iflag,extrap); if (iflag/=0_ip) return
    call dintrv(tr,nr+kr,rval,ilor,leftr,iflag,extrap); if (iflag/=0_ip) return

    iflag = 0_ip

    ! x -> y, z, q, r
    kcolr = leftr - kr
    do r=1_ip,kr
        kcolr = kcolr + 1_ip
        kcolq = leftq - kq
        do q=1_ip,kq
            kcolq = kcolq + 1_ip
            kcolz = leftz - kz
            do k=1_ip,kz
                kcolz = kcolz + 1_ip
                kcoly = lefty - ky
                do j=1_ip,ky
                    kcoly = kcoly + 1_ip
                    call dbvalu(tx,bcoef(:,kcoly,kcolz,kcolq,kcolr),&
                                nx,kx,idx,xval,inbvx,w0,iflag,w4(j,k,q,r),&
                                extrap)
                    if (iflag/=0_ip) return
                end do
            end do
        end do
    end do

    ! y -> z, q, r
    kcoly = lefty - ky + 1_ip
    do r=1_ip,kr
        do q=1_ip,kq
            do k=1_ip,kz
                call dbvalu(ty(kcoly:),w4(:,k,q,r),ky,ky,idy,yval,inbvy,&
                            w0,iflag,w3(k,q,r),extrap)
                if (iflag/=0_ip) return
            end do
        end do
    end do

    ! z -> q, r
    kcolz = leftz - kz + 1_ip
    do r=1_ip,kr
        do q=1_ip,kq
            call dbvalu(tz(kcolz:),w3(:,q,r),kz,kz,idz,zval,inbvz,&
                        w0,iflag,w2(q,r),extrap)
            if (iflag/=0_ip) return
        end do
    end do

    ! q -> r
    kcolq = leftq - kq + 1_ip
    do r=1_ip,kr
        call dbvalu(tq(kcolq:),w2(:,r),kq,kq,idq,qval,inbvq,&
                    w0,iflag,w1(r),extrap)
        if (iflag/=0_ip) return
    end do

    ! r
    kcolr = leftr - kr + 1_ip
    call dbvalu(tr(kcolr:),w1,kr,kr,idr,rval,inbvr,w0,iflag,f,extrap)

    end subroutine db5val