dbtpcf Subroutine

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

Arguments

Type IntentOptional AttributesName
real(kind=wp) :: x(n)
integer, intent(in) :: n
real(kind=wp) :: fcn(ldf,nf)
integer, intent(in) :: ldf
integer, intent(in) :: nf
real(kind=wp) :: t(*)
integer, intent(in) :: k
real(kind=wp) :: bcoef(nf,n)
real(kind=wp) :: work(*)
integer, intent(out) :: iflag

Called by

proc~~dbtpcf~~CalledByGraph proc~dbtpcf dbtpcf proc~db1ink db1ink proc~db1ink->proc~dbtpcf proc~db2ink db2ink proc~db2ink->proc~dbtpcf

Contents

Source Code


Source Code

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

    integer,intent(in)  :: n
    integer,intent(in)  :: nf
    integer,intent(in)  :: ldf
    integer,intent(in)  :: k
    real(wp)            :: x(n)
    real(wp)            :: fcn(ldf,nf)
    real(wp)            :: t(*)
    real(wp)            :: bcoef(nf,n)
    real(wp)            :: work(*)
    integer,intent(out) :: iflag  !!   0: no errors
                                  !! 301: n should be >0

    integer :: i, j, m1, m2, iq, iw

    ! check for null input

    if (nf > 0)  then

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

        ! compute b-spline coefficients

        ! first data set

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

            !  all remaining data sets by back-substitution

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

    else
        write(error_unit,'(A)') 'dbtpcf - n should be >0'
        iflag = 301
    endif

    endsubroutine dbtpcf