ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
splder.f
Go to the documentation of this file.
1  subroutine splder(t,n,c,k,nu,x,y,m,wrk,ier)
2  implicit none
3 c subroutine splder evaluates in a number of points x(i),i=1,2,...,m
4 c the derivative of order nu of a spline s(x) of degree k,given in
5 c its b-spline representation.
6 c
7 c calling sequence:
8 c call splder(t,n,c,k,nu,x,y,m,wrk,ier)
9 c
10 c input parameters:
11 c t : array,length n, which contains the position of the knots.
12 c n : integer, giving the total number of knots of s(x).
13 c c : array,length n, which contains the b-spline coefficients.
14 c k : integer, giving the degree of s(x).
15 c nu : integer, specifying the order of the derivative. 0<=nu<=k
16 c x : array,length m, which contains the points where the deriv-
17 c ative of s(x) must be evaluated.
18 c m : integer, giving the number of points where the derivative
19 c of s(x) must be evaluated
20 c wrk : real array of dimension n. used as working space.
21 c
22 c output parameters:
23 c y : array,length m, giving the value of the derivative of s(x)
24 c at the different points.
25 c ier : error flag
26 c ier = 0 : normal return
27 c ier =10 : invalid input data (see restrictions)
28 c
29 c restrictions:
30 c 0 <= nu <= k
31 c m >= 1
32 c t(k+1) <= x(i) <= x(i+1) <= t(n-k) , i=1,2,...,m-1.
33 c
34 c other subroutines required: fpbspl
35 c
36 c references :
37 c de boor c : on calculating with b-splines, j. approximation theory
38 c 6 (1972) 50-62.
39 c cox m.g. : the numerical evaluation of b-splines, j. inst. maths
40 c applics 10 (1972) 134-149.
41 c dierckx p. : curve and surface fitting with splines, monographs on
42 c numerical analysis, oxford university press, 1993.
43 c
44 c author :
45 c p.dierckx
46 c dept. computer science, k.u.leuven
47 c celestijnenlaan 200a, b-3001 heverlee, belgium.
48 c e-mail : Paul.Dierckx@cs.kuleuven.ac.be
49 c
50 c latest update : march 1987
51 c
52 c ..scalar arguments..
53  integer n,k,nu,m,ier
54 c ..array arguments..
55  real*8 t(n),c(n),x(m),y(m),wrk(n)
56 c ..local scalars..
57  integer i,j,kk,k1,k2,l,ll,l1,l2,nk1,nk2,nn
58  real*8 ak,arg,fac,sp,tb,te
59 c ..local arrays ..
60  real*8 h(6)
61 c before starting computations a data check is made. if the input data
62 c are invalid control is immediately repassed to the calling program.
63  ier = 10
64  if(nu.lt.0 .or. nu.gt.k) go to 200
65  if(m-1) 200,30,10
66  10 do 20 i=2,m
67  if(x(i).lt.x(i-1)) go to 200
68  20 continue
69  30 ier = 0
70 c fetch tb and te, the boundaries of the approximation interval.
71  k1 = k+1
72  nk1 = n-k1
73  tb = t(k1)
74  te = t(nk1+1)
75 c the derivative of order nu of a spline of degree k is a spline of
76 c degree k-nu,the b-spline coefficients wrk(i) of which can be found
77 c using the recurrence scheme of de boor.
78  l = 1
79  kk = k
80  nn = n
81  do 40 i=1,nk1
82  wrk(i) = c(i)
83  40 continue
84  if(nu.eq.0) go to 100
85  nk2 = nk1
86  do 60 j=1,nu
87  ak = kk
88  nk2 = nk2-1
89  l1 = l
90  do 50 i=1,nk2
91  l1 = l1+1
92  l2 = l1+kk
93  fac = t(l2)-t(l1)
94  if(fac.le.0.) go to 50
95  wrk(i) = ak*(wrk(i+1)-wrk(i))/fac
96  50 continue
97  l = l+1
98  kk = kk-1
99  60 continue
100  if(kk.ne.0) go to 100
101 c if nu=k the derivative is a piecewise constant function
102  j = 1
103  do 90 i=1,m
104  arg = x(i)
105  70 if(arg.lt.t(l+1) .or. l.eq.nk1) go to 80
106  l = l+1
107  j = j+1
108  go to 70
109  80 y(i) = wrk(j)
110  90 continue
111  go to 200
112  100 l = k1
113  l1 = l+1
114  k2 = k1-nu
115 c main loop for the different points.
116  do 180 i=1,m
117 c fetch a new x-value arg.
118  arg = x(i)
119  if(arg.lt.tb) arg = tb
120  if(arg.gt.te) arg = te
121 c search for knot interval t(l) <= arg < t(l+1)
122  140 if(arg.lt.t(l1) .or. l.eq.nk1) go to 150
123  l = l1
124  l1 = l+1
125  go to 140
126 c evaluate the non-zero b-splines of degree k-nu at arg.
127  150 call fpbspl(t,n,kk,arg,l,h)
128 c find the value of the derivative at x=arg.
129  sp = 0.
130  ll = l-k1
131  do 160 j=1,k2
132  ll = ll+1
133  sp = sp+wrk(ll)*h(j)
134  160 continue
135  y(i) = sp
136  180 continue
137  200 return
138  end
subroutine splder(t, n, c, k, nu, x, y, m, wrk, ier)
Definition: splder.f:1
subroutine fpbspl(t, n, k, x, l, h)
Definition: fpbspl.f:1