splde Function

private function splde(me, ndim, x, nderiv, coef, xmin, xmax, nodes, ierror)

N-dimensional cubic spline derivative evaluation.

A grid of evenly spaced nodes in ndim space is defined by the arguments xmin, xmax and nodes. A linear basis for the class of natural splines on these nodes is formed, and to each basis function corresponds a coefficient in the array coef (computed in splcc or splcw). Using nderiv to indicate the appropriate partial derivatives, each basis function is evaluated at the point x in ndim space. These values are then multiplied by the corresponding coefficient and summed to form the function result.

See also

Note

The original version of this routine would stop for an error. Now it just returns.

Note

coef, xmin, xmax and nodes must be exactly retained from the call to splcc (or splcw).

Type Bound

splpak_type

Arguments

Type IntentOptional Attributes Name
class(splpak_type), intent(inout) :: me
integer, intent(in) :: ndim

the dimensionality of the problem. the spline is a function of ndim variables or coordinates and thus a point in the independent variable space is an ndim vector. ndim must be in the range 1 <= ndim <= 4.

real(kind=wp), intent(in) :: x(ndim)

an ndim vector describing the point in the independent variable space at which the spline is to be evaluated.

integer, intent(in) :: nderiv(ndim)

an ndim vector of integers specifying the partial derivative to be evaluated. the order of the derivative along the idim axis is nderiv(idim). these integers must be in the range 0 <= nderiv(idim) <= 2.

real(kind=wp), intent(out) :: coef(*)

the array of coefficients which determine the spline. each coefficient corresponds to a particular basis function which in turn corresponds to a node in the node grid. this correspondence between the node grid and the array coef is as if coef were an ndim-dimensional fortran array with dimensions nodes(1),...,nodes(ndim), i.e., to store the array linearly, the leftmost indices are incremented most frequently. coef may be computed by using routines splcc or splcw.

the dimension is assumed to be coef(nodes(1)*...*nodes(ndim)).

real(kind=wp), intent(in) :: xmin(ndim)

a vector describing the lower extreme corner of the node grid. a set of evenly spaced nodes is formed along each coordinate axis and xmin(idim) is the location of the first node along the idim axis.

real(kind=wp), intent(in) :: xmax(ndim)

a vector describing the upper extreme corner of the node grid. a set of evenly spaced nodes is formed along each coordinate axis and xmax(idim) is the location of the last node along the idim axis.

integer, intent(in) :: nodes(ndim)

a vector of integers describing the number of nodes along each axis. nodes(idim) is the the number of nodes (counting endpoints) along the idim axis and determines the flexibility of the spline in that coordinate direction. nodes(idim) must be >= 4 but may be as large as the arrays coef and work allow.

note: the node grid is completely defined by the arguments xmin, xmax and nodes. the spacing of this grid in the idim coordinate direction is dx(idim) = (xmax(idim)-xmin(idim)) / (nodes(idim)-1). a node in this grid may be indexed by an ndim vector of integers (in(1),...,in(ndim)) where 1 <= in(idim) <= nodes(idim). the location of such a node may be represented by an ndim vector (x(1),...,x(ndim)) where x(idim) = xmin(idim)+(in(idim)-1) * dx(idim).

integer, intent(out) :: ierror

an error flag with the following meanings:

  • 0 no error.
  • 101 ndim is < 1 or is > 4.
  • 102 nodes(idim) is < 4 for some idim.
  • 103 xmin(idim) = xmax(idim) for some idim.
  • 104 nderiv(idim) is < 0 or is > 2 for some idim.

Return Value real(kind=wp)

the function value returned is the partial derivative (indicated by nderiv) of the spline evaluated at x.


Calls

proc~~splde~~CallsGraph proc~splde splpak_module::splpak_type%splde proc~bascmp splpak_module::splpak_type%bascmp proc~splde->proc~bascmp proc~cfaerr splpak_module::cfaerr proc~splde->proc~cfaerr

Called by

proc~~splde~~CalledByGraph proc~splde splpak_module::splpak_type%splde proc~splfe splpak_module::splpak_type%splfe proc~splfe->proc~splde

Source Code

function splde(me,ndim,x,nderiv,coef,xmin,xmax,nodes,ierror)

    class(splpak_type),intent(inout) :: me
    real(wp) :: splde !! the function value returned is the partial
                      !! derivative (indicated by `nderiv`) of the
                      !! spline evaluated at `x`.
    integer,intent(in) :: ndim !! the dimensionality of the problem.  the
                               !! spline is a function of `ndim` variables or
                               !! coordinates and thus a point in the
                               !! independent variable space is an `ndim` vector.
                               !! `ndim` must be in the range `1 <= ndim <= 4`.
    real(wp),intent(in) :: x(ndim) !! an `ndim` vector describing the point in the
                                   !! independent variable space at which the
                                   !! spline is to be evaluated.
    real(wp),intent(out) :: coef(*) !! the array of coefficients which determine the
                                    !! spline.  each coefficient corresponds to a
                                    !! particular basis function which in turn
                                    !! corresponds to a node in the node grid.  this
                                    !! correspondence between the node grid and the
                                    !! array `coef` is as if `coef` were an
                                    !! ndim-dimensional fortran array with
                                    !! dimensions `nodes(1),...,nodes(ndim)`, i.e., to
                                    !! store the array linearly, the leftmost
                                    !! indices are incremented most frequently.
                                    !! coef may be computed by using routines [[splcc]]
                                    !! or [[splcw]].
                                    !!
                                    !! the dimension is assumed to be
                                    !! `coef(nodes(1)*...*nodes(ndim))`.
    real(wp),intent(in) :: xmin(ndim) !! a vector describing the lower extreme corner
                                      !! of the node grid.  a set of evenly spaced
                                      !! nodes is formed along each coordinate axis
                                      !! and `xmin(idim)` is the location of the first
                                      !! node along the `idim` axis.
    real(wp),intent(in) :: xmax(ndim) !! a vector describing the upper extreme corner
                                      !! of the node grid.  a set of evenly spaced
                                      !! nodes is formed along each coordinate axis
                                      !! and `xmax(idim)` is the location of the last
                                      !! node along the `idim` axis.
    integer,intent(in) :: nderiv(ndim) !! an `ndim` vector of integers specifying the
                                       !! partial derivative to be evaluated.  the
                                       !! order of the derivative along the `idim` axis
                                       !! is `nderiv(idim)`.  these integers must be in
                                       !! the range `0 <= nderiv(idim) <= 2`.
    integer,intent(in) :: nodes(ndim) !! a vector of integers describing the number of
                                      !! nodes along each axis.  `nodes(idim)` is the
                                      !! the number of nodes (counting endpoints)
                                      !! along the `idim` axis and determines the
                                      !! flexibility of the spline in that coordinate
                                      !! direction.  `nodes(idim)` must be >= 4 but
                                      !! may be as large as the arrays `coef` and `work`
                                      !! allow.
                                      !!
                                      !! *note:*  the node grid is completely defined by
                                      !! the arguments `xmin`, `xmax` and `nodes`.
                                      !! the spacing of this grid in the `idim`
                                      !! coordinate direction is
                                      !! `dx(idim) = (xmax(idim)-xmin(idim)) / (nodes(idim)-1)`.
                                      !! a node in this grid may be indexed by
                                      !! an `ndim` vector of integers
                                      !! `(in(1),...,in(ndim))` where
                                      !! `1 <= in(idim) <= nodes(idim)`.  the
                                      !! location of such a node may be
                                      !! represented by an `ndim` vector
                                      !! `(x(1),...,x(ndim)) ` where
                                      !! `x(idim) = xmin(idim)+(in(idim)-1) * dx(idim)`.
    integer,intent(out) :: ierror !! an error flag with the following meanings:
                                  !!
                                  !! *   0  no error.
                                  !! * 101  ndim is < 1 or is > 4.
                                  !! * 102  nodes(idim) is < 4 for some idim.
                                  !! * 103  xmin(idim) = xmax(idim) for some idim.
                                  !! * 104  nderiv(idim) is < 0 or is > 2 for some idim.

    real(wp) :: xrng,sum,basm
    integer :: iibmx,idim,nod,it,iib,icof

    ierror = 0
    me%mdim = ndim
    if (me%mdim<1) then
        ierror = 101
        call cfaerr(ierror, &
            ' splfe or splde - NDIM is less than 1')
        return
    end if
    iibmx = 1
    do idim = 1,me%mdim
        nod = nodes(idim)
        if (nod<4) then
            ierror = 102
            call cfaerr(ierror, &
                ' splfe or splde - NODES(IDIM) is less than  4for some IDIM')
            return
        end if
        xrng = xmax(idim) - xmin(idim)
        if (xrng==0.0_wp) then
            ierror = 103
            call cfaerr(ierror, &
                ' splfe or splde - XMIN(IDIM) = XMAX(IDIM) for some IDIM')
            return
        end if
        if (nderiv(idim)<0 .or. nderiv(idim)>2) then
            ierror = 104
            call cfaerr(ierror, &
                ' splde - NDERIV(IDIM) IS less than 0 or greater than 2 for some IDIM')
        end if

        !  DX(IDIM) is the node spacing along the IDIM coordinate.
        me%dx(idim) = xrng/real(nod-1,wp)
        me%dxin(idim) = 1.0_wp/me%dx(idim)

        !  Compute indices of basis functions which are nonzero at X.
        it = me%dxin(idim)*(x(idim)-xmin(idim))

        !  IBMN must be in the range 0 to NODES-2.
        me%ibmn(idim) = min(max(it-1,0),nod-2)

        !  IBMX must be in the range 1 to NODES-1.
        me%ibmx(idim) = max(min(it+2,nod-1),1)
        iibmx = iibmx* (me%ibmx(idim)-me%ibmn(idim)+1)
        me%ib(idim) = me%ibmn(idim)
    end do

    sum = 0.0_wp
    iib = 0

    basis_index : do
        !  Begining of basis index loop - traverse all indices corresponding
        !  to basis functions which are nonzero at X.
        iib = iib + 1

        !  The indices are in IB and are passed through common to BASCMP.
        call me%bascmp(x,nderiv,xmin,nodes,icof,basm)

        !  BASCMP computes ICOF and BASM where BASM is the value at X of the
        !  N-dimensional basis function corresponding to COEF(ICOF).
        sum = sum + coef(icof)*basm
        if (iib<iibmx) then
            !  Increment the basis indices.
            do idim = 1,me%mdim
                me%ib(idim) = me%ib(idim) + 1
                if (me%ib(idim)<=me%ibmx(idim)) cycle basis_index
                me%ib(idim) = me%ibmn(idim)
            end do
        end if

        exit basis_index !  End of basis index loop.
    end do basis_index

    splde = sum

end function splde