gtkIOStream  1.7.0
GTK+ << C++ IOStream operators for GTK+. Now with ORBing, numerical computation, audio client and more ...
bsplpp.f
Go to the documentation of this file.
1  subroutine bsplpp ( t, bcoef, n, k, scrtch, break, coef, l )
2 
3 c*********************************************************************72
4 c
5 cc BSPLPP converts from B-spline to piecewise polynomial form.
6 c
7 c from * a practical guide to splines * by c. de Boor (7 may 92)
8 calls bsplvb
9 c
10 converts the b-representation t, bcoef, n, k of some spline into its
11 c pp-representation break, coef, l, k .
12 c
13 c****** i n p u t ******
14 c t.....knot sequence, of length n+k
15 c bcoef.....b-spline coefficient sequence, of length n
16 c n.....length of bcoef and dimension of spline space spline(k,t)
17 c k.....order of the spline
18 c
19 c w a r n i n g . . . the restriction k .le. kmax (= 20) is impo-
20 c sed by the arbitrary dimension statement for biatx below, but
21 c is n o w h e r e c h e c k e d for.
22 c
23 c****** w o r k a r e a ******
24 c scrtch......of size (k,k) , needed to contain bcoeffs of a piece of
25 c the spline and its k-1 derivatives
26 c
27 c****** o u t p u t ******
28 c break.....breakpoint sequence, of length l+1, contains (in increas-
29 c ing order) the distinct points in the sequence t(k),...,t(n+1)
30 c coef.....array of size (k,l), with coef(i,j) = (i-1)st derivative of
31 c spline at break(j) from the right
32 c l.....number of polynomial pieces which make up the spline in the in-
33 c terval (t(k), t(n+1))
34 c
35 c****** m e t h o d ******
36 c for each breakpoint interval, the k relevant b-coeffs of the
37 c spline are found and then differenced repeatedly to get the b-coeffs
38 c of all the derivatives of the spline on that interval. the spline and
39 c its first k-1 derivatives are then evaluated at the left end point
40 c of that interval, using bsplvb repeatedly to obtain the values of
41 c all b-splines of the appropriate order at that point.
42 c
43  implicit none
44 
45  integer k,l,n, i,j,jp1,kmax,kmj,left,lsofar
46  parameter(kmax = 20)
47  double precision bcoef(n),break(l+1),coef(k,l),t(n+k), scrtch(k,k)
48  & ,biatx(kmax),diff,factor,sum
49 c
50  lsofar = 0
51  break(1) = t(k)
52  do 50 left=k,n
53 c find the next nontrivial knot interval.
54  if (t(left+1) .eq. t(left)) go to 50
55  lsofar = lsofar + 1
56  break(lsofar+1) = t(left+1)
57  if (k .gt. 1) go to 9
58  coef(1,lsofar) = bcoef(left)
59  go to 50
60 c store the k b-spline coeff.s relevant to current knot interval
61 c in scrtch(.,1) .
62  9 do 10 i=1,k
63  10 scrtch(i,1) = bcoef(left-k+i)
64 c
65 c for j=1,...,k-1, compute the k-j b-spline coeff.s relevant to
66 c current knot interval for the j-th derivative by differencing
67 c those for the (j-1)st derivative, and store in scrtch(.,j+1) .
68  do 20 jp1=2,k
69  j = jp1 - 1
70  kmj = k - j
71  do 20 i=1,kmj
72  diff = t(left+i) - t(left+i - kmj)
73  if (diff .gt. 0.0d+00) scrtch(i,jp1) =
74  & (scrtch(i+1,j)-scrtch(i,j))/diff
75  20 continue
76 c
77 c for j = 0, ..., k-1, find the values at t(left) of the j+1
78 c b-splines of order j+1 whose support contains the current
79 c knot interval from those of order j (in biatx ), then comb-
80 c ine with the b-spline coeff.s (in scrtch(.,k-j) ) found earlier
81 c to compute the (k-j-1)st derivative at t(left) of the given
82 c spline.
83 c note. if the repeated calls to bsplvb are thought to gene-
84 c rate too much overhead, then replace the first call by
85 c biatx(1) = 1.
86 c and the subsequent call by the statement
87 c j = jp1 - 1
88 c followed by a direct copy of the lines
89 c deltar(j) = t(left+j) - x
90 c ......
91 c biatx(j+1) = saved
92 c from bsplvb . deltal(kmax) and deltar(kmax) would have to
93 c appear in a dimension statement, of course.
94 c
95  call bsplvb ( t, 1, 1, t(left), left, biatx )
96  coef(k,lsofar) = scrtch(1,k)
97  do 30 jp1=2,k
98  call bsplvb ( t, jp1, 2, t(left), left, biatx )
99  kmj = k+1 - jp1
100  sum = 0.0d+00
101  do 28 i=1,jp1
102  28 sum = biatx(i)*scrtch(i,kmj) + sum
103  30 coef(kmj,lsofar) = sum
104  50 continue
105  l = lsofar
106  if (k .eq. 1) return
107  factor = 1.0d+00
108  do 60 i=2,k
109  factor = factor*dble(k+1-i)
110  do 60 j=1,lsofar
111  60 coef(i,j) = coef(i,j)*factor
112  return
113  end
subroutine bsplpp(t, bcoef, n, k, scrtch, break, coef, l)
Definition: bsplpp.f:2
subroutine bsplvb(t, jhigh, index, x, left, biatx)
Definition: bsplvb.f:2
gtkIOStream: /tmp/gtkiostream/src/deBoor/bsplpp.f Source File
GTK+ IOStream  Beta