Free cookie consent management tool by TermsFeed Policy Generator

source: branches/2789_MathNetNumerics-Exploration/HeuristicLab.Algorithms.DataAnalysis.Experimental/sbart/bvalue.f @ 16409

Last change on this file since 16409 was 15467, checked in by gkronber, 7 years ago

#2789 dllexport bvalue function

File size: 5.2 KB
Line 
1      module modbvalue
2      contains
3      real function bvalue ( t, bcoef, n, k, x, jderiv )
4     & bind(C,name='bvalue')
5      USE ISO_C_BINDING
6      !DEC$ ATTRIBUTES DLLEXPORT::bvalue
7c  from  * a practical guide to splines *  by c. de boor   
8calls  interv
9c
10calculates value at  x  of  jderiv-th derivative of spline from b-repr.
11c  the spline is taken to be continuous from the right, EXCEPT at the
12c  rightmost knot, where it is taken to be continuous from the left.
13c
14c******  i n p u t ******
15c  t, bcoef, n, k......forms the b-representation of the spline  f  to
16c        be evaluated. specifically,
17c  t.....knot sequence, of length  n+k, assumed nondecreasing.
18c  bcoef.....b-coefficient sequence, of length  n .
19c  n.....length of  bcoef  and dimension of spline(k,t),
20c        a s s u m e d  positive .
21c  k.....order of the spline .
22c
23c  w a r n i n g . . .   the restriction  k .le. kmax (=20)  is imposed
24c        arbitrarily by the dimension statement for  aj, dl, dr  below,
25c        but is  n o w h e r e  c h e c k e d  for.
26c
27c  x.....the point at which to evaluate .
28c  jderiv.....integer giving the order of the derivative to be evaluated
29c        a s s u m e d  to be zero or positive.
30c
31c******  o u t p u t  ******
32c  bvalue.....the value of the (jderiv)-th derivative of  f  at  x .
33c
34c******  m e t h o d  ******
35c     The nontrivial knot interval  (t(i),t(i+1))  containing  x  is lo-
36c  cated with the aid of  interv . The  k  b-coeffs of  f  relevant for
37c  this interval are then obtained from  bcoef (or taken to be zero if
38c  not explicitly available) and are then differenced  jderiv  times to
39c  obtain the b-coeffs of  (d**jderiv)f  relevant for that interval.
40c  Precisely, with  j = jderiv, we have from x.(12) of the text that
41c
42c     (d**j)f  =  sum ( bcoef(.,j)*b(.,k-j,t) )
43c
44c  where
45c                   / bcoef(.),                     ,  j .eq. 0
46c                   /
47c    bcoef(.,j)  =  / bcoef(.,j-1) - bcoef(.-1,j-1)
48c                   / ----------------------------- ,  j .gt. 0
49c                   /    (t(.+k-j) - t(.))/(k-j)
50c
51c     Then, we use repeatedly the fact that
52c
53c    sum ( a(.)*b(.,m,t)(x) )  =  sum ( a(.,x)*b(.,m-1,t)(x) )
54c  with
55c                 (x - t(.))*a(.) + (t(.+m-1) - x)*a(.-1)
56c    a(.,x)  =    ---------------------------------------
57c                 (x - t(.))      + (t(.+m-1) - x)
58c
59c  to write  (d**j)f(x)  eventually as a linear combination of b-splines
60c  of order  1 , and the coefficient for  b(i,1,t)(x)  must then be the
61c  desired number  (d**j)f(x). (see x.(17)-(19) of text).
62c
63      integer jderiv,k,n,   i,ilo,imk,j,jc,jcmin,jcmax,jj,kmax,kmj,km1
64     *                     ,mflag,nmi,jdrvp1
65      parameter (kmax = 20)
66      real bcoef(n),t(n+k),x,   aj(kmax),dl(kmax),dr(kmax),fkmj
67      bvalue = 0.
68      if (jderiv .ge. k)                go to 99
69c
70c  *** Find  i   s.t.   1 .le. i .lt. n+k   and   t(i) .lt. t(i+1)   and
71c      t(i) .le. x .lt. t(i+1) . If no such i can be found,  x  lies
72c      outside the support of  the spline  f , hence  bvalue = 0.
73c      (The asymmetry in this choice of  i  makes  f  rightcontinuous, except
74c      at  t(n+k) where it is leftcontinuous.)
75      call interv ( t, n+k, x, i, mflag )
76      if (mflag .ne. 0)                 go to 99
77c  *** if k = 1 (and jderiv = 0), bvalue = bcoef(i).
78      km1 = k - 1
79      if (km1 .gt. 0)                   go to 1
80      bvalue = bcoef(i)
81                                        go to 99
82c
83c  *** store the k b-spline coefficients relevant for the knot interval
84c     (t(i),t(i+1)) in aj(1),...,aj(k) and compute dl(j) = x - t(i+1-j),
85c     dr(j) = t(i+j) - x, j=1,...,k-1 . set any of the aj not obtainable
86c     from input to zero. set any t.s not obtainable equal to t(1) or
87c     to t(n+k) appropriately.
88    1 jcmin = 1
89      imk = i - k
90      if (imk .ge. 0)                   go to 8
91      jcmin = 1 - imk
92      do 5 j=1,i
93    5    dl(j) = x - t(i+1-j)
94      do 6 j=i,km1
95         aj(k-j) = 0.
96    6    dl(j) = dl(i)
97                                        go to 10
98    8 do 9 j=1,km1
99    9    dl(j) = x - t(i+1-j)
100c
101   10 jcmax = k
102      nmi = n - i
103      if (nmi .ge. 0)                   go to 18
104      jcmax = k + nmi
105      do 15 j=1,jcmax
106   15    dr(j) = t(i+j) - x
107      do 16 j=jcmax,km1
108         aj(j+1) = 0.
109   16    dr(j) = dr(jcmax)
110                                        go to 20
111   18 do 19 j=1,km1
112   19    dr(j) = t(i+j) - x
113c
114   20 do 21 jc=jcmin,jcmax
115   21    aj(jc) = bcoef(imk + jc)
116c
117c               *** difference the coefficients  jderiv  times.
118      if (jderiv .eq. 0)                go to 30
119      do 23 j=1,jderiv
120         kmj = k-j
121         fkmj = float(kmj)
122         ilo = kmj
123         do 23 jj=1,kmj
124            aj(jj) = ((aj(jj+1) - aj(jj))/(dl(ilo) + dr(jj)))*fkmj
125   23       ilo = ilo - 1
126c
127c  *** compute value at  x  in (t(i),t(i+1)) of jderiv-th derivative,
128c     given its relevant b-spline coeffs in aj(1),...,aj(k-jderiv).
129   30 if (jderiv .eq. km1)              go to 39
130      jdrvp1 = jderiv + 1     
131      do 33 j=jdrvp1,km1
132         kmj = k-j
133         ilo = kmj
134         do 33 jj=1,kmj
135            aj(jj) = (aj(jj+1)*dl(ilo) + aj(jj)*dr(jj))/(dl(ilo)+dr(jj))
136   33       ilo = ilo - 1
137   39 bvalue = aj(1)
138c
139   99                                   return
140      end
141      end
Note: See TracBrowser for help on using the repository browser.