dbtpcf Subroutine

private pure subroutine dbtpcf(x, n, fcn, ldf, nf, t, k, bcoef, work, iflag)

dbtpcf computes b-spline interpolation coefficients for nf sets of data stored in the columns of the array fcn. the b-spline coefficients are stored in the rows of bcoef however. each interpolation is based on the n abcissa stored in the array x, and the n+k knots stored in the array t. the order of each interpolation is k.

History

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

Arguments

Type IntentOptional Attributes Name
real(kind=wp), intent(in), dimension(:) :: x
integer(kind=ip), intent(in) :: n

dimension of x

real(kind=wp), intent(in), dimension(ldf,nf) :: fcn
integer(kind=ip), intent(in) :: ldf
integer(kind=ip), intent(in) :: nf
real(kind=wp), intent(in), dimension(:) :: t
integer(kind=ip), intent(in) :: k
real(kind=wp), intent(out), dimension(nf,n) :: bcoef
real(kind=wp), intent(out), dimension(*) :: work

work array of size >= 2*k*(n+1)

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

status flag:

  • 0: no errors
  • 301: n should be >0

Calls

proc~~dbtpcf~~CallsGraph proc~dbtpcf bspline_sub_module::dbtpcf proc~dbintk bspline_sub_module::dbintk proc~dbtpcf->proc~dbintk proc~dbnslv bspline_sub_module::dbnslv proc~dbtpcf->proc~dbnslv proc~dbintk->proc~dbnslv proc~dbnfac bspline_sub_module::dbnfac proc~dbintk->proc~dbnfac proc~dbspvn bspline_sub_module::dbspvn proc~dbintk->proc~dbspvn

Called by

proc~~dbtpcf~~CalledByGraph proc~dbtpcf bspline_sub_module::dbtpcf proc~db1ink_default bspline_sub_module::db1ink_default proc~db1ink_default->proc~dbtpcf proc~db2ink bspline_sub_module::db2ink proc~db2ink->proc~dbtpcf proc~db3ink bspline_sub_module::db3ink proc~db3ink->proc~dbtpcf proc~db4ink bspline_sub_module::db4ink proc~db4ink->proc~dbtpcf proc~db5ink bspline_sub_module::db5ink proc~db5ink->proc~dbtpcf proc~db6ink bspline_sub_module::db6ink proc~db6ink->proc~dbtpcf interface~db1ink bspline_sub_module::db1ink interface~db1ink->proc~db1ink_default proc~initialize_2d_auto_knots bspline_oo_module::bspline_2d%initialize_2d_auto_knots proc~initialize_2d_auto_knots->proc~db2ink proc~initialize_2d_specify_knots bspline_oo_module::bspline_2d%initialize_2d_specify_knots proc~initialize_2d_specify_knots->proc~db2ink proc~initialize_3d_auto_knots bspline_oo_module::bspline_3d%initialize_3d_auto_knots proc~initialize_3d_auto_knots->proc~db3ink proc~initialize_3d_specify_knots bspline_oo_module::bspline_3d%initialize_3d_specify_knots proc~initialize_3d_specify_knots->proc~db3ink proc~initialize_4d_auto_knots bspline_oo_module::bspline_4d%initialize_4d_auto_knots proc~initialize_4d_auto_knots->proc~db4ink proc~initialize_4d_specify_knots bspline_oo_module::bspline_4d%initialize_4d_specify_knots proc~initialize_4d_specify_knots->proc~db4ink proc~initialize_5d_auto_knots bspline_oo_module::bspline_5d%initialize_5d_auto_knots proc~initialize_5d_auto_knots->proc~db5ink proc~initialize_5d_specify_knots bspline_oo_module::bspline_5d%initialize_5d_specify_knots proc~initialize_5d_specify_knots->proc~db5ink proc~initialize_6d_auto_knots bspline_oo_module::bspline_6d%initialize_6d_auto_knots proc~initialize_6d_auto_knots->proc~db6ink proc~initialize_6d_specify_knots bspline_oo_module::bspline_6d%initialize_6d_specify_knots proc~initialize_6d_specify_knots->proc~db6ink proc~bspline_2d_constructor_auto_knots bspline_oo_module::bspline_2d_constructor_auto_knots proc~bspline_2d_constructor_auto_knots->proc~initialize_2d_auto_knots proc~bspline_2d_constructor_specify_knots bspline_oo_module::bspline_2d_constructor_specify_knots proc~bspline_2d_constructor_specify_knots->proc~initialize_2d_specify_knots proc~bspline_3d_constructor_auto_knots bspline_oo_module::bspline_3d_constructor_auto_knots proc~bspline_3d_constructor_auto_knots->proc~initialize_3d_auto_knots proc~bspline_3d_constructor_specify_knots bspline_oo_module::bspline_3d_constructor_specify_knots proc~bspline_3d_constructor_specify_knots->proc~initialize_3d_specify_knots proc~bspline_4d_constructor_auto_knots bspline_oo_module::bspline_4d_constructor_auto_knots proc~bspline_4d_constructor_auto_knots->proc~initialize_4d_auto_knots proc~bspline_4d_constructor_specify_knots bspline_oo_module::bspline_4d_constructor_specify_knots proc~bspline_4d_constructor_specify_knots->proc~initialize_4d_specify_knots proc~bspline_5d_constructor_auto_knots bspline_oo_module::bspline_5d_constructor_auto_knots proc~bspline_5d_constructor_auto_knots->proc~initialize_5d_auto_knots proc~bspline_5d_constructor_specify_knots bspline_oo_module::bspline_5d_constructor_specify_knots proc~bspline_5d_constructor_specify_knots->proc~initialize_5d_specify_knots proc~bspline_6d_constructor_auto_knots bspline_oo_module::bspline_6d_constructor_auto_knots proc~bspline_6d_constructor_auto_knots->proc~initialize_6d_auto_knots proc~bspline_6d_constructor_specify_knots bspline_oo_module::bspline_6d_constructor_specify_knots proc~bspline_6d_constructor_specify_knots->proc~initialize_6d_specify_knots proc~initialize_1d_auto_knots bspline_oo_module::bspline_1d%initialize_1d_auto_knots proc~initialize_1d_auto_knots->interface~db1ink proc~initialize_1d_specify_knots bspline_oo_module::bspline_1d%initialize_1d_specify_knots proc~initialize_1d_specify_knots->interface~db1ink interface~bspline_2d bspline_oo_module::bspline_2d interface~bspline_2d->proc~bspline_2d_constructor_auto_knots interface~bspline_2d->proc~bspline_2d_constructor_specify_knots interface~bspline_3d bspline_oo_module::bspline_3d interface~bspline_3d->proc~bspline_3d_constructor_auto_knots interface~bspline_3d->proc~bspline_3d_constructor_specify_knots interface~bspline_4d bspline_oo_module::bspline_4d interface~bspline_4d->proc~bspline_4d_constructor_auto_knots interface~bspline_4d->proc~bspline_4d_constructor_specify_knots interface~bspline_5d bspline_oo_module::bspline_5d interface~bspline_5d->proc~bspline_5d_constructor_auto_knots interface~bspline_5d->proc~bspline_5d_constructor_specify_knots interface~bspline_6d bspline_oo_module::bspline_6d interface~bspline_6d->proc~bspline_6d_constructor_auto_knots interface~bspline_6d->proc~bspline_6d_constructor_specify_knots proc~bspline_1d_constructor_auto_knots bspline_oo_module::bspline_1d_constructor_auto_knots proc~bspline_1d_constructor_auto_knots->proc~initialize_1d_auto_knots proc~bspline_1d_constructor_specify_knots bspline_oo_module::bspline_1d_constructor_specify_knots proc~bspline_1d_constructor_specify_knots->proc~initialize_1d_specify_knots interface~bspline_1d bspline_oo_module::bspline_1d interface~bspline_1d->proc~bspline_1d_constructor_auto_knots interface~bspline_1d->proc~bspline_1d_constructor_specify_knots

Source Code

    pure subroutine dbtpcf(x,n,fcn,ldf,nf,t,k,bcoef,work,iflag)

    integer(ip),intent(in)                :: n  !! dimension of `x`
    integer(ip),intent(in)                :: nf
    integer(ip),intent(in)                :: ldf
    integer(ip),intent(in)                :: k
    real(wp),dimension(:),intent(in)      :: x
    real(wp),dimension(ldf,nf),intent(in) :: fcn
    real(wp),dimension(:),intent(in)      :: t
    real(wp),dimension(nf,n),intent(out)  :: bcoef
    real(wp),dimension(*),intent(out)     :: work   !! work array of size >= `2*k*(n+1)`
    integer(ip),intent(out)               :: iflag  !! status flag:
                                                    !!
                                                    !! * 0: no errors
                                                    !! * 301: n should be >0

    integer(ip) :: i, j, m1, m2, iq, iw

    ! check for null input

    if (nf > 0_ip)  then

        ! partition work array
        m1 = k - 1_ip
        m2 = m1 + k
        iq = 1_ip + n
        iw = iq + m2*n+1_ip

        ! compute b-spline coefficients

        ! first data set

        call dbintk(x,fcn,t,n,k,work,work(iq),work(iw),iflag)
        if (iflag == 0_ip) then
            do i=1_ip,n
                bcoef(1_ip,i) = work(i)
            end do

            !  all remaining data sets by back-substitution

            if (nf == 1_ip)  return
            do j=2_ip,nf
                do i=1_ip,n
                    work(i) = fcn(i,j)
                end do
                call dbnslv(work(iq),m2,n,m1,m1,work)
                do i=1_ip,n
                    bcoef(j,i) = work(i)
                end do
            end do
        end if

    else
        !write(error_unit,'(A)') 'dbtpcf - n should be >0'
        iflag = 301_ip
    end if

    end subroutine dbtpcf