************80
! PUTIT puts together one block of the collocation equation system.
Method:
The K collocation equations for the interval ( T(LEFT), T(LEFT+1) )
are constructed with the aid of the subroutine DIFEQU( 2, ., . )
and interspersed (in order) with the side conditions, if any, in
this interval, using DIFEQU ( 3, ., . ) for the information.
The block Q has KPM columns, corresponding to the KPM B-splines of order
KPM which have the interval ( T(LEFT), T(LEFT+1) ) in their support.
The block's diagonal is part of the diagonal of the total system.
The first equation in this block not overlapped by the preceding block
is therefore equation LOWROW, with LOWROW = number of side conditions
in preceding intervals (or blocks).
Modified:
14 February 2007
Author:
Carl DeBoor
Reference:
Carl DeBoor,
A Practical Guide to Splines,
Springer, 2001,
ISBN: 0387953663,
LC: QA1.A647.v27.
Parameters:
Input, real ( kind = 8 ) T(LEFT+KPM), the knot sequence.
Input, integer ( kind = 4 ) KPM, the order of the spline.
Input, integer ( kind = 4 ) LEFT, indicates the interval of interest,
that is, the interval ( T(LEFT), T(LEFT+1) ).
Workspace, real ( kind = 8 ) SCRTCH(KPM,KPM).
Workspace, real ( kind = 8 ) DBIATX(KPM,M+1), derivatives of B-splines,
with DBIATX(J,I+1) containing the I-th derivative of the J-th B-spline
of interest.
Output, real ( kind = 8 ) Q(NROW,KPM), the block.
Input, integer ( kind = 4 ) NROW, number of rows in block to be
! put together.
Output, real ( kind = 8 ) B(NROW), the corresponding piece of
the right hand side.
Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
real(kind=8) | :: | t(left+kpm) | ||||
integer(kind=4) | :: | kpm | ||||
integer(kind=4) | :: | left | ||||
real(kind=8) | :: | scrtch(kpm,kpm) | ||||
real(kind=8) | :: | dbiatx(kpm,*) | ||||
real(kind=8) | :: | q(nrow,kpm) | ||||
integer(kind=4) | :: | nrow | ||||
real(kind=8) | :: | b(nrow) |
subroutine putit ( t, kpm, left, scrtch, dbiatx, q, nrow, b ) !*****************************************************************************80 ! !! PUTIT puts together one block of the collocation equation system. ! ! Method: ! ! The K collocation equations for the interval ( T(LEFT), T(LEFT+1) ) ! are constructed with the aid of the subroutine DIFEQU( 2, ., . ) ! and interspersed (in order) with the side conditions, if any, in ! this interval, using DIFEQU ( 3, ., . ) for the information. ! ! The block Q has KPM columns, corresponding to the KPM B-splines of order ! KPM which have the interval ( T(LEFT), T(LEFT+1) ) in their support. ! ! The block's diagonal is part of the diagonal of the total system. ! ! The first equation in this block not overlapped by the preceding block ! is therefore equation LOWROW, with LOWROW = number of side conditions ! in preceding intervals (or blocks). ! ! Modified: ! ! 14 February 2007 ! ! Author: ! ! Carl DeBoor ! ! Reference: ! ! Carl DeBoor, ! A Practical Guide to Splines, ! Springer, 2001, ! ISBN: 0387953663, ! LC: QA1.A647.v27. ! ! Parameters: ! ! Input, real ( kind = 8 ) T(LEFT+KPM), the knot sequence. ! ! Input, integer ( kind = 4 ) KPM, the order of the spline. ! ! Input, integer ( kind = 4 ) LEFT, indicates the interval of interest, ! that is, the interval ( T(LEFT), T(LEFT+1) ). ! ! Workspace, real ( kind = 8 ) SCRTCH(KPM,KPM). ! ! Workspace, real ( kind = 8 ) DBIATX(KPM,M+1), derivatives of B-splines, ! with DBIATX(J,I+1) containing the I-th derivative of the J-th B-spline ! of interest. ! ! Output, real ( kind = 8 ) Q(NROW,KPM), the block. ! ! Input, integer ( kind = 4 ) NROW, number of rows in block to be !! put together. ! ! Output, real ( kind = 8 ) B(NROW), the corresponding piece of ! the right hand side. ! implicit none integer ( kind = 4 ) kpm integer ( kind = 4 ) left integer ( kind = 4 ) nrow real ( kind = 8 ) b(nrow) real ( kind = 8 ) dbiatx(kpm,*) real ( kind = 8 ) dx integer ( kind = 4 ) i integer ( kind = 4 ) irow integer ( kind = 4 ) iside integer ( kind = 4 ) itermx integer ( kind = 4 ) j integer ( kind = 4 ) k integer ( kind = 4 ) ll integer ( kind = 4 ) lowrow integer ( kind = 4 ) m integer ( kind = 4 ) mode integer ( kind = 4 ) mp1 real ( kind = 8 ) q(nrow,kpm) real ( kind = 8 ) rho real ( kind = 8 ) scrtch(kpm,kpm) real ( kind = 8 ) sum1 real ( kind = 8 ) t(left+kpm) real ( kind = 8 ) v(20) real ( kind = 8 ) xm real ( kind = 8 ) xside real ( kind = 8 ) xx save / other / save / side / common / other / itermx, k, rho(19) common / side / m, iside, xside(10) mp1 = m + 1 q(1:nrow,1:kpm) = 0.0D+00 xm = ( t(left+1) + t(left) ) / 2.0D+00 dx = ( t(left+1) - t(left) ) / 2.0D+00 ll = 1 lowrow = iside do irow = lowrow, nrow if ( k < ll ) then go to 20 end if mode = 2 ! ! Next collocation point: ! xx = xm + dx * rho(ll) ll = ll + 1 ! ! The corresponding collocation equation is next unless the next side ! condition occurs at a point at, or to the left of, the next ! collocation point. ! if ( m < iside ) then go to 30 end if if ( xx < xside(iside) ) then go to 30 end if ll = ll - 1 20 continue mode = 3 xx = xside(iside) 30 continue call difequ ( mode, xx, v ) ! ! The next equation, a collocation equation (MODE=2) or a side ! condition (MODE=3), reads ! ! (*) (V(M+1)*D**M+V(M)*D**(M-1) +...+ V(1)*D**0)F(XX) = V(M+2) ! ! in terms of the information supplied by DIFEQU. ! ! The corresponding equation for the B-spline coefficients of F therefore ! has the left side of (*), evaluated at each of the KPM B-splines having ! XX in their support, as its KPM possibly nonzero coefficients. ! call bsplvd ( t, kpm, xx, left, scrtch, dbiatx, mp1 ) do j = 1, kpm q(irow,j) = dot_product ( dbiatx(j,1:mp1), v(1:mp1) ) end do b(irow) = v(m+2) end do return end subroutine putit