1 | *** For comments, see sbart.r, written by Finbarr O'Sullivan |
---|
2 | *** ratfor output of 14 May 1987 version |
---|
3 | *** with bug fix in setreg, 27 Sep 1990 |
---|
4 | module modsbart |
---|
5 | contains |
---|
6 | subroutine sbart(xs,ys,ws,n,knot,nk,coef,sz,lev,crit,icrit,spar, |
---|
7 | &ispar,lspar,uspar,tol,isetup,xwy,hs0,hs1,hs2,hs3,sg0,sg1,sg2,sg3, |
---|
8 | &abd,p1ip,p2ip,ld4,ldnk,ier) bind(c, name='sbart') |
---|
9 | USE ISO_C_BINDING |
---|
10 | !DEC$ ATTRIBUTES DLLEXPORT::sbart |
---|
11 | real xs(n),ys(n),ws(n),knot(nk+4),coef(nk),sz(n),lev(n),crit,spar, |
---|
12 | &lspar,uspar,tol,xwy(nk),hs0(nk),hs1(nk),hs2(nk),hs3(nk), sg0(nk), |
---|
13 | &sg1(nk),sg2(nk),sg3(nk),abd(ld4,nk),p1ip(ld4,nk),p2ip(ldnk,nk) |
---|
14 | integer n,nk,isetup,icrit,ispar,ld4,ldnk,ier |
---|
15 | realt1,t2,ratio, a,b,c,d,e,eps,xm,p,q,r,tol1,tol2,u,v,w, fu,fv,fw, |
---|
16 | &fx,x,ax,bx |
---|
17 | integer i |
---|
18 | if(.not.(isetup.eq.0))goto 23000 |
---|
19 | call sgram(sg0,sg1,sg2,sg3,knot,nk) |
---|
20 | call stxwx(xs,ys,ws,n,knot,nk,xwy,hs0,hs1,hs2,hs3) |
---|
21 | t1=0. |
---|
22 | t2=0. |
---|
23 | do 23002 i=3,nk-3 |
---|
24 | t1 = t1 + hs0(i) |
---|
25 | 23002 continue |
---|
26 | do 23004 i=3,nk-3 |
---|
27 | t2 = t2 + sg0(i) |
---|
28 | 23004 continue |
---|
29 | ratio = t1/t2 |
---|
30 | isetup = 1 |
---|
31 | 23000 continue |
---|
32 | if(.not.(ispar.eq.1))goto 23006 |
---|
33 | call sslvrg(xs,ys,ws,n,knot,nk,coef,sz,lev,crit,icrit,spar,ratio, |
---|
34 | &xwy,hs0,hs1,hs2,hs3,sg0,sg1,sg2,sg3,abd,p1ip,p2ip,ld4,ldnk,ier) |
---|
35 | return |
---|
36 | 23006 continue |
---|
37 | ax=lspar |
---|
38 | bx=uspar |
---|
39 | c = 0.5*(3. - sqrt(5.0)) |
---|
40 | eps = 1.00 |
---|
41 | 10 eps = eps/2.00 |
---|
42 | tol1 = 1.0 + eps |
---|
43 | if(.not.(tol1 .gt. 1.00))goto 23008 |
---|
44 | go to 10 |
---|
45 | 23008 continue |
---|
46 | eps = sqrt(eps) |
---|
47 | a = ax |
---|
48 | b = bx |
---|
49 | v = a + c*(b - a) |
---|
50 | w = v |
---|
51 | x = v |
---|
52 | e = 0.0 |
---|
53 | spar = x |
---|
54 | call sslvrg(xs,ys,ws,n,knot,nk,coef,sz,lev,crit,icrit,spar,ratio, |
---|
55 | &xwy,hs0,hs1,hs2,hs3,sg0,sg1,sg2,sg3,abd,p1ip,p2ip,ld4,ldnk,ier) |
---|
56 | fx = crit |
---|
57 | fv = fx |
---|
58 | fw = fx |
---|
59 | 20 xm = 0.5*(a + b) |
---|
60 | tol1 = eps*abs(x) + tol/3.0 |
---|
61 | tol2 = 2.0*tol1 |
---|
62 | if(.not.(abs(x - xm) .le. (tol2 - 0.5*(b - a))))goto 23010 |
---|
63 | go to 90 |
---|
64 | 23010 continue |
---|
65 | if(.not.(abs(e) .le. tol1))goto 23012 |
---|
66 | go to 40 |
---|
67 | 23012 continue |
---|
68 | r = (x - w)*(fx - fv) |
---|
69 | q = (x - v)*(fx - fw) |
---|
70 | p = (x - v)*q - (x - w)*r |
---|
71 | q = 2.00*(q - r) |
---|
72 | if(.not.(q .gt. 0.0))goto 23014 |
---|
73 | p = -p |
---|
74 | 23014 continue |
---|
75 | q = abs(q) |
---|
76 | r = e |
---|
77 | e = d |
---|
78 | 30 if(.not.(abs(p) .ge. abs(0.5*q*r)))goto 23016 |
---|
79 | go to 40 |
---|
80 | 23016 continue |
---|
81 | if(.not.(p .le. q*(a - x)))goto 23018 |
---|
82 | go to 40 |
---|
83 | 23018 continue |
---|
84 | if(.not.(p .ge. q*(b - x)))goto 23020 |
---|
85 | go to 40 |
---|
86 | 23020 continue |
---|
87 | d = p/q |
---|
88 | u = x + d |
---|
89 | if(.not.((u - a) .lt. tol2))goto 23022 |
---|
90 | d = sign(tol1, xm - x) |
---|
91 | 23022 continue |
---|
92 | if(.not.((b - u) .lt. tol2))goto 23024 |
---|
93 | d = sign(tol1, xm - x) |
---|
94 | 23024 continue |
---|
95 | go to 50 |
---|
96 | 40 if(.not.(x .ge. xm))goto 23026 |
---|
97 | e = a - x |
---|
98 | 23026 continue |
---|
99 | if(.not.(x .lt. xm))goto 23028 |
---|
100 | e = b - x |
---|
101 | 23028 continue |
---|
102 | d = c*e |
---|
103 | 50 if(.not.(abs(d) .ge. tol1))goto 23030 |
---|
104 | u = x + d |
---|
105 | 23030 continue |
---|
106 | if(.not.(abs(d) .lt. tol1))goto 23032 |
---|
107 | u = x + sign(tol1, d) |
---|
108 | 23032 continue |
---|
109 | spar = u |
---|
110 | call sslvrg(xs,ys,ws,n,knot,nk,coef,sz,lev,crit,icrit,spar,ratio, |
---|
111 | &xwy,hs0,hs1,hs2,hs3,sg0,sg1,sg2,sg3,abd,p1ip,p2ip,ld4,ldnk,ier) |
---|
112 | fu = crit |
---|
113 | if(.not.(fu .gt. fx))goto 23034 |
---|
114 | go to 60 |
---|
115 | 23034 continue |
---|
116 | if(.not.(u .ge. x))goto 23036 |
---|
117 | a = x |
---|
118 | 23036 continue |
---|
119 | if(.not.(u .lt. x))goto 23038 |
---|
120 | b = x |
---|
121 | 23038 continue |
---|
122 | v = w |
---|
123 | fv = fw |
---|
124 | w = x |
---|
125 | fw = fx |
---|
126 | x = u |
---|
127 | fx = fu |
---|
128 | go to 20 |
---|
129 | 60 if(.not.(u .lt. x))goto 23040 |
---|
130 | a = u |
---|
131 | 23040 continue |
---|
132 | if(.not.(u .ge. x))goto 23042 |
---|
133 | b = u |
---|
134 | 23042 continue |
---|
135 | if(.not.(fu .le. fw))goto 23044 |
---|
136 | go to 70 |
---|
137 | 23044 continue |
---|
138 | if(.not.(w .eq. x))goto 23046 |
---|
139 | go to 70 |
---|
140 | 23046 continue |
---|
141 | if(.not.(fu .le. fv))goto 23048 |
---|
142 | go to 80 |
---|
143 | 23048 continue |
---|
144 | if(.not.(v .eq. x))goto 23050 |
---|
145 | go to 80 |
---|
146 | 23050 continue |
---|
147 | if(.not.(v .eq. w))goto 23052 |
---|
148 | go to 80 |
---|
149 | 23052 continue |
---|
150 | go to 20 |
---|
151 | 70 v = w |
---|
152 | fv = fw |
---|
153 | w = u |
---|
154 | fw = fu |
---|
155 | go to 20 |
---|
156 | 80 v = u |
---|
157 | fv = fu |
---|
158 | go to 20 |
---|
159 | 90 continue |
---|
160 | spar = x |
---|
161 | crit = fx |
---|
162 | return |
---|
163 | 23007 continue |
---|
164 | return |
---|
165 | end |
---|
166 | subroutine sgram(sg0,sg1,sg2,sg3,tb,nb) |
---|
167 | real sg0(nb),sg1(nb),sg2(nb),sg3(nb),tb(nb+4),vnikx(4,3),work(16), |
---|
168 | &yw1(4),yw2(4),wpt |
---|
169 | integer nb,ileft,ilo,mflag,i,ii,jj |
---|
170 | do 23054 i=1,nb |
---|
171 | sg0(i)=0. |
---|
172 | sg1(i)=0. |
---|
173 | sg2(i)=0. |
---|
174 | sg3(i)=0. |
---|
175 | 23054 continue |
---|
176 | ilo = 1 |
---|
177 | do 23056 i=1,nb |
---|
178 | call interv(tb(1),(nb+1),tb(i),ileft,mflag) |
---|
179 | call bsplvd (tb,4,tb(i),ileft,work,vnikx,3) |
---|
180 | do 23058 ii=1,4 |
---|
181 | yw1(ii) = vnikx(ii,3) |
---|
182 | 23058 continue |
---|
183 | call bsplvd (tb,4,tb(i+1),ileft,work,vnikx,3) |
---|
184 | do 23060 ii=1,4 |
---|
185 | yw2(ii) = vnikx(ii,3) - yw1(ii) |
---|
186 | 23060 continue |
---|
187 | wpt = tb(i+1) - tb(i) |
---|
188 | if(.not.(ileft.ge.4))goto 23062 |
---|
189 | do 23064 ii=1,4 |
---|
190 | jj=ii |
---|
191 | sg0(ileft-4+ii) = sg0(ileft-4+ii) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
192 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
193 | jj=ii+1 |
---|
194 | if(.not.(jj.le.4))goto 23066 |
---|
195 | sg1(ileft+ii-4) = sg1(ileft+ii-4) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
196 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
197 | 23066 continue |
---|
198 | jj=ii+2 |
---|
199 | if(.not.(jj.le.4))goto 23068 |
---|
200 | sg2(ileft+ii-4) = sg2(ileft+ii-4) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
201 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
202 | 23068 continue |
---|
203 | jj=ii+3 |
---|
204 | if(.not.(jj.le.4))goto 23070 |
---|
205 | sg3(ileft+ii-4) = sg3(ileft+ii-4) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
206 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
207 | 23070 continue |
---|
208 | 23064 continue |
---|
209 | goto 23063 |
---|
210 | 23062 continue |
---|
211 | if(.not.(ileft.eq.3))goto 23072 |
---|
212 | do 23074 ii=1,3 |
---|
213 | jj=ii |
---|
214 | sg0(ileft-3+ii) = sg0(ileft-3+ii) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
215 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
216 | jj=ii+1 |
---|
217 | if(.not.(jj.le.3))goto 23076 |
---|
218 | sg1(ileft+ii-3) = sg1(ileft+ii-3) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
219 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
220 | 23076 continue |
---|
221 | jj=ii+2 |
---|
222 | if(.not.(jj.le.3))goto 23078 |
---|
223 | sg2(ileft+ii-3) = sg2(ileft+ii-3) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
224 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
225 | 23078 continue |
---|
226 | 23074 continue |
---|
227 | goto 23073 |
---|
228 | 23072 continue |
---|
229 | if(.not.(ileft.eq.2))goto 23080 |
---|
230 | do 23082 ii=1,2 |
---|
231 | jj=ii |
---|
232 | sg0(ileft-2+ii) = sg0(ileft-2+ii) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
233 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
234 | jj=ii+1 |
---|
235 | if(.not.(jj.le.2))goto 23084 |
---|
236 | sg1(ileft+ii-2) = sg1(ileft+ii-2) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
237 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
238 | 23084 continue |
---|
239 | 23082 continue |
---|
240 | goto 23081 |
---|
241 | 23080 continue |
---|
242 | if(.not.(ileft.eq.1))goto 23086 |
---|
243 | do 23088 ii=1,1 |
---|
244 | jj=ii |
---|
245 | sg0(ileft-1+ii) = sg0(ileft-1+ii) +wpt* (yw1(ii)*yw1(jj) + (yw2( |
---|
246 | &ii)*yw1(jj) + yw2(jj)*yw1(ii))*.50 +yw2(ii)*yw2(jj)*.3330 ) |
---|
247 | 23088 continue |
---|
248 | 23086 continue |
---|
249 | 23081 continue |
---|
250 | 23073 continue |
---|
251 | 23063 continue |
---|
252 | 23056 continue |
---|
253 | return |
---|
254 | end |
---|
255 | subroutine sinerp(abd,ld4,nk,p1ip,p2ip,ldnk,flag) |
---|
256 | realabd(ld4,nk),p1ip(ld4,nk),p2ip(ldnk,nk),wjm3(3),wjm2(2),wjm1(1) |
---|
257 | &,c0,c1,c2,c3 |
---|
258 | integerflag,ld4,nk,ldnk,i,j,k |
---|
259 | wjm3(1)=0. |
---|
260 | wjm3(2)=0. |
---|
261 | wjm3(1)=0. |
---|
262 | wjm2(1)=0. |
---|
263 | wjm2(2)=0. |
---|
264 | wjm1(1)=0. |
---|
265 | do 23090 i=1,nk |
---|
266 | j=nk-i+1 |
---|
267 | c0 = 1./abd(4,j) |
---|
268 | if(.not.(j.le.nk-3))goto 23092 |
---|
269 | c1 = abd(1,j+3)*c0 |
---|
270 | c2 = abd(2,j+2)*c0 |
---|
271 | c3 = abd(3,j+1)*c0 |
---|
272 | goto 23093 |
---|
273 | 23092 continue |
---|
274 | if(.not.(j.eq.nk-2))goto 23094 |
---|
275 | c1 = 0. |
---|
276 | c2 = abd(2,j+2)*c0 |
---|
277 | c3 = abd(3,j+1)*c0 |
---|
278 | goto 23095 |
---|
279 | 23094 continue |
---|
280 | if(.not.(j.eq.nk-1))goto 23096 |
---|
281 | c1 = 0. |
---|
282 | c2 = 0. |
---|
283 | c3 = abd(3,j+1)*c0 |
---|
284 | goto 23097 |
---|
285 | 23096 continue |
---|
286 | if(.not.(j.eq.nk))goto 23098 |
---|
287 | c1 = 0. |
---|
288 | c2 = 0. |
---|
289 | c3 = 0. |
---|
290 | 23098 continue |
---|
291 | 23097 continue |
---|
292 | 23095 continue |
---|
293 | 23093 continue |
---|
294 | p1ip(1,j) = 0. - (c1*wjm3(1)+c2*wjm3(2)+c3*wjm3(3)) |
---|
295 | p1ip(2,j) = 0. - (c1*wjm3(2)+c2*wjm2(1)+c3*wjm2(2)) |
---|
296 | p1ip(3,j) = 0. - (c1*wjm3(3)+c2*wjm2(2)+c3*wjm1(1)) |
---|
297 | p1ip(4,j) = c0**2 +c1**2*wjm3(1)+2.*c1*c2*wjm3(2)+2.*c1*c3*wjm3(3) |
---|
298 | & +c2**2*wjm2(1)+2.*c2*c3*wjm2(2) +c3**2*wjm1(1) |
---|
299 | wjm3(1)=wjm2(1) |
---|
300 | wjm3(2)=wjm2(2) |
---|
301 | wjm3(3)=p1ip(2,j) |
---|
302 | wjm2(1)=wjm1(1) |
---|
303 | wjm2(2)=p1ip(3,j) |
---|
304 | wjm1(1)=p1ip(4,j) |
---|
305 | 23090 continue |
---|
306 | if(.not.(flag.eq.0))goto 23100 |
---|
307 | return |
---|
308 | 23100 continue |
---|
309 | do 23102 i=1,nk |
---|
310 | j=nk-i+1 |
---|
311 | k=1 |
---|
312 | 23104 if(.not.(k.le.4.and.j+k-1.le.nk))goto 23106 |
---|
313 | p2ip(j,j+k-1) = p1ip(5-k,j) |
---|
314 | k=k+1 |
---|
315 | goto 23104 |
---|
316 | 23106 continue |
---|
317 | 23102 continue |
---|
318 | do 23107 i=1,nk |
---|
319 | j=nk-i+1 |
---|
320 | k=j-4 |
---|
321 | 23109 if(.not.(k.ge.1))goto 23111 |
---|
322 | c0 = 1./abd(4,k) |
---|
323 | c1 = abd(1,k+3)*c0 |
---|
324 | c2 = abd(2,k+2)*c0 |
---|
325 | c3 = abd(3,k+1)*c0 |
---|
326 | p2ip(k,j) = 0. - ( c1*p2ip(k+3,j) +c2*p2ip(k+2,j) +c3*p2ip(k+1,j) |
---|
327 | &) |
---|
328 | k=k-1 |
---|
329 | goto 23109 |
---|
330 | 23111 continue |
---|
331 | 23107 continue |
---|
332 | return |
---|
333 | 23101 continue |
---|
334 | end |
---|
335 | subroutine sslvrg(x,y,w,n,knot,nk,coef,sz,lev,crit,icrit,spar, |
---|
336 | &ratio,xwy,hs0,hs1,hs2,hs3,sg0,sg1,sg2,sg3,abd,p1ip,p2ip,ld4,ldnk, |
---|
337 | &info) |
---|
338 | USE MODBVALUE |
---|
339 | real x(n),y(n),w(n),knot(nk+4),coef(nk),sz(n),lev(n),crit,ratio, |
---|
340 | &spar,xwy(nk),hs0(nk),hs1(nk),hs2(nk),hs3(nk), sg0(nk),sg1(nk),sg2( |
---|
341 | &nk),sg3(nk),abd(ld4,nk),p1ip(ld4,nk),p2ip(ldnk,nk),lambda,b0,b1, |
---|
342 | &b2,b3,eps,vnikx(4,1),work(16),xv,rss,df |
---|
343 | integer n,nk,icrit,ld4,ldnk,i,icoef,ileft,ilo,info,j,mflag |
---|
344 | ilo = 1 |
---|
345 | eps = .1e-10 |
---|
346 | lambda = ratio*16.**(-2. + spar*(6.)) |
---|
347 | do 23112 i=1,nk |
---|
348 | coef(i) = xwy(i) |
---|
349 | 23112 continue |
---|
350 | do 23114 i=1,nk |
---|
351 | abd(4,i) = hs0(i)+lambda*sg0(i) |
---|
352 | 23114 continue |
---|
353 | do 23116 i=1,(nk-1) |
---|
354 | abd(3,i+1) = hs1(i)+lambda*sg1(i) |
---|
355 | 23116 continue |
---|
356 | do 23118 i=1,(nk-2) |
---|
357 | abd(2,i+2) = hs2(i)+lambda*sg2(i) |
---|
358 | 23118 continue |
---|
359 | do 23120 i=1,(nk-3) |
---|
360 | abd(1,i+3) = hs3(i)+lambda*sg3(i) |
---|
361 | 23120 continue |
---|
362 | call spbfa(abd,ld4,nk,3,info) |
---|
363 | if(.not.(info.ne.0))goto 23122 |
---|
364 | return |
---|
365 | 23122 continue |
---|
366 | call spbsl(abd,ld4,nk,3,coef) |
---|
367 | icoef = 1 |
---|
368 | do 23124 i=1,n |
---|
369 | xv = x(i) |
---|
370 | sz(i) = bvalue(knot,coef,nk,4,xv,0) |
---|
371 | 23124 continue |
---|
372 | if(.not.(icrit.eq.0))goto 23126 |
---|
373 | return |
---|
374 | 23126 continue |
---|
375 | call sinerp(abd,ld4,nk,p1ip,p2ip,ldnk,0) |
---|
376 | do 23128 i=1,n |
---|
377 | xv = x(i) |
---|
378 | call interv(knot(1),(nk+1),xv,ileft,mflag) |
---|
379 | if(.not.(mflag.eq.-1))goto 23130 |
---|
380 | ileft = 4 |
---|
381 | xv = knot(4)+eps |
---|
382 | 23130 continue |
---|
383 | if(.not.(mflag.eq.1))goto 23132 |
---|
384 | ileft = nk |
---|
385 | xv = knot(nk+1)-eps |
---|
386 | 23132 continue |
---|
387 | j=ileft-3 |
---|
388 | call bsplvd(knot,4,xv,ileft,work,vnikx,1) |
---|
389 | b0=vnikx(1,1) |
---|
390 | b1=vnikx(2,1) |
---|
391 | b2=vnikx(3,1) |
---|
392 | b3=vnikx(4,1) |
---|
393 | lev(i) = (p1ip(4,j)*b0**2 + 2.*p1ip(3,j)*b0*b1 +2.*p1ip(2,j)*b0* |
---|
394 | &b2 + 2.*p1ip(1,j)*b0*b3 +p1ip(4,j+1)*b1**2 + 2.*p1ip(3,j+1)*b1*b2 |
---|
395 | &+2.*p1ip(2,j+1)*b1*b3 +p1ip(4,j+2)*b2**2 + 2.*p1ip(3,j+2)*b2*b3 + |
---|
396 | &p1ip(4,j+3)*b3**2 )*w(i)**2 |
---|
397 | 23128 continue |
---|
398 | if(.not.(icrit.eq.1))goto 23134 |
---|
399 | rss = 0. |
---|
400 | df = 0. |
---|
401 | do 23136 i=1,n |
---|
402 | rss = rss + ((y(i)-sz(i))*w(i))**2 |
---|
403 | 23136 continue |
---|
404 | do 23138 i=1,n |
---|
405 | df = df + 1.-lev(i) |
---|
406 | 23138 continue |
---|
407 | crit = (rss/n)/((df/n)**2) |
---|
408 | goto 23135 |
---|
409 | 23134 continue |
---|
410 | crit = 0. |
---|
411 | do 23140 i=1,n |
---|
412 | crit = crit +(((y(i)-sz(i))*w(i))/(1-lev(i)))**2 |
---|
413 | 23140 continue |
---|
414 | 23135 continue |
---|
415 | return |
---|
416 | 23127 continue |
---|
417 | end |
---|
418 | subroutine stxwx(x,z,w,k,xknot,n,y,hs0,hs1,hs2,hs3) |
---|
419 | real z(k),w(k),x(k),xknot(n+4),y(n),hs0(n),hs1(n),hs2(n),hs3(n), |
---|
420 | &eps,vnikx(4,1),work(16) |
---|
421 | integer k,n,j,i,ilo,ileft,mflag |
---|
422 | do 23142 i=1,n |
---|
423 | y(i)=0. |
---|
424 | hs0(i)=0. |
---|
425 | hs1(i)=0. |
---|
426 | hs2(i)=0. |
---|
427 | hs3(i)=0. |
---|
428 | 23142 continue |
---|
429 | ilo=1 |
---|
430 | eps = .1e-9 |
---|
431 | do 23144 i=1,k |
---|
432 | call interv(xknot(1),(n+1),x(i),ileft,mflag) |
---|
433 | if(.not.(mflag.eq.-1))goto 23146 |
---|
434 | write(6,'("Error in hess ",i2)')mflag |
---|
435 | stop |
---|
436 | 23146 continue |
---|
437 | if(.not.(mflag.eq. 1))goto 23148 |
---|
438 | if(.not.(x(i).le.(xknot(ileft)+eps)))goto 23150 |
---|
439 | ileft=ileft-1 |
---|
440 | goto 23151 |
---|
441 | 23150 continue |
---|
442 | write(6,'("Error in hess ",i2)')mflag |
---|
443 | stop |
---|
444 | 23151 continue |
---|
445 | 23148 continue |
---|
446 | call bsplvd (xknot,4,x(i),ileft,work,vnikx,1) |
---|
447 | j= ileft-4+1 |
---|
448 | y(j) = y(j)+w(i)**2*z(i)*vnikx(1,1) |
---|
449 | hs0(j)=hs0(j)+w(i)**2*vnikx(1,1)**2 |
---|
450 | hs1(j)=hs1(j)+w(i)**2*vnikx(1,1)*vnikx(2,1) |
---|
451 | hs2(j)=hs2(j)+w(i)**2*vnikx(1,1)*vnikx(3,1) |
---|
452 | hs3(j)=hs3(j)+w(i)**2*vnikx(1,1)*vnikx(4,1) |
---|
453 | j= ileft-4+2 |
---|
454 | y(j) = y(j)+w(i)**2*z(i)*vnikx(2,1) |
---|
455 | hs0(j)=hs0(j)+w(i)**2*vnikx(2,1)**2 |
---|
456 | hs1(j)=hs1(j)+w(i)**2*vnikx(2,1)*vnikx(3,1) |
---|
457 | hs2(j)=hs2(j)+w(i)**2*vnikx(2,1)*vnikx(4,1) |
---|
458 | j= ileft-4+3 |
---|
459 | y(j) = y(j)+w(i)**2*z(i)*vnikx(3,1) |
---|
460 | hs0(j)=hs0(j)+w(i)**2*vnikx(3,1)**2 |
---|
461 | hs1(j)=hs1(j)+w(i)**2*vnikx(3,1)*vnikx(4,1) |
---|
462 | j= ileft-4+4 |
---|
463 | y(j) = y(j)+w(i)**2*z(i)*vnikx(4,1) |
---|
464 | hs0(j)=hs0(j)+w(i)**2*vnikx(4,1)**2 |
---|
465 | 23144 continue |
---|
466 | return |
---|
467 | end |
---|
468 | subroutine sknotl(x,n,knot,k) bind(C,name='sknotl') |
---|
469 | USE ISO_C_BINDING |
---|
470 | !DEC$ ATTRIBUTES DLLEXPORT::sknotl |
---|
471 | real x(n),knot(n+6),a1,a2,a3,a4 |
---|
472 | integer n,k,ndk,j |
---|
473 | a1 = log(50.)/log(2.) |
---|
474 | a2 = log(100.)/log(2.) |
---|
475 | a3 = log(140.)/log(2.) |
---|
476 | a4 = log(200.)/log(2.) |
---|
477 | if(.not.(n.lt.50))goto 23152 |
---|
478 | ndk = n |
---|
479 | goto 23153 |
---|
480 | 23152 continue |
---|
481 | if(.not.(n.ge.50 .and. n.lt.200))goto 23154 |
---|
482 | ndk = 2.**(a1+(a2-a1)*(n-50.)/150.) |
---|
483 | goto 23155 |
---|
484 | 23154 continue |
---|
485 | if(.not.(n.ge.200 .and. n.lt.800))goto 23156 |
---|
486 | ndk = 2.**(a2+(a3-a2)*(n-200.)/600.) |
---|
487 | goto 23157 |
---|
488 | 23156 continue |
---|
489 | if(.not.(n.ge.800 .and. n.lt.3200))goto 23158 |
---|
490 | ndk = 2.**(a3+(a4-a3)*(n-800.)/2400.) |
---|
491 | goto 23159 |
---|
492 | 23158 continue |
---|
493 | if(.not.(n.ge.3200))goto 23160 |
---|
494 | ndk = 200. + (n-3200)**.2 |
---|
495 | 23160 continue |
---|
496 | 23159 continue |
---|
497 | 23157 continue |
---|
498 | 23155 continue |
---|
499 | 23153 continue |
---|
500 | k = ndk + 6 |
---|
501 | do 23162 j=1,3 |
---|
502 | knot(j) = x(1) |
---|
503 | 23162 continue |
---|
504 | do 23164 j=1,ndk |
---|
505 | knot(j+3) = x( 1 + (j-1)*(n-1)/(ndk-1) ) |
---|
506 | 23164 continue |
---|
507 | do 23166 j=1,3 |
---|
508 | knot(ndk+3+j) = x(n) |
---|
509 | 23166 continue |
---|
510 | return |
---|
511 | end |
---|
512 | subroutine setreg(x,y,w,n,xw,nx,min,range,knot,nk) |
---|
513 | . bind(C,name='setreg') |
---|
514 | USE ISO_C_BINDING |
---|
515 | !DEC$ ATTRIBUTES DLLEXPORT::setreg |
---|
516 | real x(n),y(n),w(n),xw(n),min,range,knot(n+6) |
---|
517 | integer n,nk,nx |
---|
518 | real eps |
---|
519 | integer ycnt,i,k |
---|
520 | call scopy(n,x,1,xw,1) |
---|
521 | C call sortpr(x,n,w) |
---|
522 | C call sortpr(xw,n,y) |
---|
523 | call ssort(x,w,n,2) |
---|
524 | call ssort(xw,y,n,2) |
---|
525 | range = x(n)-x(1) |
---|
526 | min = x(1) |
---|
527 | eps = .1e-9 |
---|
528 | do 23168 i=1,n |
---|
529 | x(i) = (x(i)-min)/range |
---|
530 | 23168 continue |
---|
531 | call scopy(n,x,1,xw,1) |
---|
532 | nx = 1 |
---|
533 | x(nx) = x(1) |
---|
534 | w(nx) = w(1) |
---|
535 | y(nx) = y(1)*w(1) |
---|
536 | i=2 |
---|
537 | 23170 if(.not.(i.le.n))goto 23172 |
---|
538 | if(.not.(xw(i).gt.x(nx)+eps))goto 23173 |
---|
539 | if(.not.(w(nx).gt.0.0))goto 23175 |
---|
540 | y(nx) = y(nx)/w(nx) |
---|
541 | 23175 continue |
---|
542 | nx = nx + 1 |
---|
543 | x(nx) = x(i) |
---|
544 | y(nx) = y(i)*w(i) |
---|
545 | w(nx) = w(i) |
---|
546 | goto 23174 |
---|
547 | 23173 continue |
---|
548 | y(nx) = y(nx)+y(i)*w(i) |
---|
549 | w(nx) = w(nx) + w(i) |
---|
550 | 23174 continue |
---|
551 | i=i+1 |
---|
552 | goto 23170 |
---|
553 | 23172 continue |
---|
554 | if(.not.(w(nx).gt.0.0))goto 23177 |
---|
555 | y(nx) = y(nx)/w(nx) |
---|
556 | 23177 continue |
---|
557 | i=1 |
---|
558 | 23179 if(.not.(i.le.nx))goto 23181 |
---|
559 | if(.not.(w(i).gt.0))goto 23182 |
---|
560 | w(i)=sqrt(w(i)) |
---|
561 | 23182 continue |
---|
562 | i=i+1 |
---|
563 | goto 23179 |
---|
564 | 23181 continue |
---|
565 | call sknotl(x,nx,knot,k) |
---|
566 | nk = k-4 |
---|
567 | return |
---|
568 | end |
---|
569 | end |
---|