MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
Numeric_subroutines.F90
1 
5 
6 
7 
8 
9 
10 SUBROUTINE bicub_derivative ( ya, x1a, x2a, y1a, y2a, y12a, i_max, k_max )
11 !
12  USE dft_module, ONLY: ndft, r_grid
13  IMPLICIT NONE
14 !
15 ! ----------------------------------------------------------------------
16  REAL, INTENT(IN) :: ya(r_grid,ndft)
17  REAL, INTENT(IN) :: x1a(r_grid)
18  REAL, INTENT(IN) :: x2a(ndft)
19  REAL, INTENT(OUT) :: y1a(r_grid,ndft)
20  REAL, INTENT(OUT) :: y2a(r_grid,ndft)
21  REAL, INTENT(OUT) :: y12a(r_grid,ndft)
22  INTEGER, INTENT(IN) :: i_max
23  INTEGER, INTENT(IN) :: k_max
24 !
25 ! ----------------------------------------------------------------------
26  INTEGER :: i, k
27 ! ----------------------------------------------------------------------
28 
29 
30 DO i = 2, i_max-1
31  DO k = 2, k_max-1
32  y1a(i,k) = (ya(i+1,k)-ya(i-1,k)) / (x1a(i+1)-x1a(i-1))
33  y2a(i,k) = (ya(i,k+1)-ya(i,k-1)) / (x2a(k+1)-x2a(k-1))
34  y12a(i,k)= (ya(i+1,k+1)-ya(i+1,k-1)-ya(i-1,k+1)+ya(i-1,k-1)) &
35  /((x1a(i+1)-x1a(i-1))*(x2a(k+1)-x2a(k-1)))
36  END DO
37 END DO
38 
39 i = 1
40 DO k = 1, k_max
41  y1a(i,k) = 0.0
42  y2a(i,k) = 0.0
43  y12a(i,k)= 0.0
44 END DO
45 
46 k = 1
47 DO i = 1, i_max
48  y1a(i,k) = 0.0
49  y2a(i,k) = 0.0
50  y12a(i,k)= 0.0
51 END DO
52 
53 
54 i = i_max
55 DO k = 2, k_max-1
56  y1a(i,k) = (ya(i,k)-ya(i-1,k)) / (x1a(i)-x1a(i-1))
57  y2a(i,k) = (ya(i,k+1)-ya(i,k-1)) / (x2a(k+1)-x2a(k-1))
58  y12a(i,k)= (ya(i,k+1)-ya(i,k-1)-ya(i-1,k+1)+ya(i-1,k-1)) &
59  /((x1a(i)-x1a(i-1))*(x2a(k+1)-x2a(k-1)))
60 END DO
61 
62 
63 k = k_max
64 DO i = 2, i_max-1
65  y1a(i,k) = (ya(i+1,k)-ya(i-1,k)) / (x1a(i+1)-x1a(i-1))
66  y2a(i,k) = (ya(i,k)-ya(i,k-1)) / (x2a(k)-x2a(k-1))
67  y12a(i,k)= (ya(i+1,k)-ya(i+1,k-1)-ya(i-1,k)+ya(i-1,k-1)) &
68  /((x1a(i+1)-x1a(i-1))*(x2a(k)-x2a(k-1)))
69 END DO
70 
71 k = k_max
72 i = i_max
73 y1a(i,k) = (ya(i,k)-ya(i-1,k)) / (x1a(i)-x1a(i-1))
74 y2a(i,k) = (ya(i,k)-ya(i,k-1)) / (x2a(k)-x2a(k-1))
75 y12a(i,k)= (ya(i,k)-ya(i,k-1)-ya(i-1,k)+ya(i-1,k-1)) &
76  /((x1a(i)-x1a(i-1))*(x2a(k)-x2a(k-1)))
77 
78 END SUBROUTINE bicub_derivative
79 
80 
81 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
82 !
83 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
84 !
85  SUBROUTINE bicub_c ( ya, x1a, x2a, y1a, y2a, y12a, c_bicub, i_max, k_max )
86 !
87  USE dft_module, ONLY: ndft, r_grid
88  IMPLICIT NONE
89 !
90 ! ----------------------------------------------------------------------
91  REAL, INTENT(IN) :: ya(r_grid,ndft)
92  REAL, INTENT(IN) :: x1a(r_grid)
93  REAL, INTENT(IN) :: x2a(ndft)
94  REAL, INTENT(IN) :: y1a(r_grid,ndft)
95  REAL, INTENT(IN) :: y2a(r_grid,ndft)
96  REAL, INTENT(IN) :: y12a(r_grid,ndft)
97  REAL, INTENT(OUT) :: c_bicub(r_grid,ndft,4,4)
98  INTEGER, INTENT(IN) :: i_max
99  INTEGER, INTENT(IN) :: k_max
100 !
101 ! ----------------------------------------------------------------------
102  INTEGER :: i, k, m, n
103  REAL :: y(4),y1(4),y2(4),y12(4),x1l,x1u,x2l,x2u
104  REAL :: c(4,4)
105 ! ----------------------------------------------------------------------
106 
107 DO i = 1, i_max-1
108  DO k = 1, k_max-1
109  y(1)=ya(i,k)
110  y(2)=ya(i+1,k)
111  y(3)=ya(i+1,k+1)
112  y(4)=ya(i,k+1)
113 
114  y1(1)=y1a(i,k)
115  y1(2)=y1a(i+1,k)
116  y1(3)=y1a(i+1,k+1)
117  y1(4)=y1a(i,k+1)
118 
119  y2(1)=y2a(i,k)
120  y2(2)=y2a(i+1,k)
121  y2(3)=y2a(i+1,k+1)
122  y2(4)=y2a(i,k+1)
123 
124  y12(1)=y12a(i,k)
125  y12(2)=y12a(i+1,k)
126  y12(3)=y12a(i+1,k+1)
127  y12(4)=y12a(i,k+1)
128 
129  x1l=x1a(i)
130  x1u=x1a(i+1)
131  x2l=x2a(k)
132  x2u=x2a(k+1)
133 
134  CALL bcucof(y,y1,y2,y12,x1u-x1l,x2u-x2l,c)
135  DO m=1,4
136  DO n=1,4
137  c_bicub(i,k,m,n)=c(m,n)
138  END DO
139  END DO
140 
141  END DO
142 END DO
143 
144 END SUBROUTINE bicub_c
145 
146 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
147 !
148 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
149 !
150  SUBROUTINE bcucof(y,y1,y2,y12,d1,d2,c)
151 !
152  IMPLICIT NONE
153 !
154 ! ----------------------------------------------------------------------
155  REAL, INTENT(IN) :: y(4)
156  REAL, INTENT(IN) :: y1(4)
157  REAL, INTENT(IN) :: y2(4)
158  REAL, INTENT(IN) :: y12(4)
159  REAL, INTENT(IN) :: d1
160  REAL, INTENT(IN) :: d2
161  REAL, INTENT(OUT) :: c(4,4)
162 !
163 ! ----------------------------------------------------------------------
164  INTEGER :: i,j,k,l
165  REAL :: d1d2,xx,cl(16),wt(16,16),x(16)
166  SAVE wt
167  DATA wt/1,0,-3,2,4*0,-3,0,9,-6,2,0,-6,4,8*0,3,0,-9,6,-2,0,6,-4,10* &
168  0,9,-6,2*0,-6,4,2*0,3,-2,6*0,-9,6,2*0,6,-4,4*0,1,0,-3,2,-2,0,6,-4, &
169  1,0,-3,2,8*0,-1,0,3,-2,1,0,-3,2,10*0,-3,2,2*0,3,-2,6*0,3,-2,2*0, &
170  -6,4,2*0,3,-2,0,1,-2,1,5*0,-3,6,-3,0,2,-4,2,9*0,3,-6,3,0,-2,4,-2, &
171  10*0,-3,3,2*0,2,-2,2*0,-1,1,6*0,3,-3,2*0,-2,2,5*0,1,-2,1,0,-2,4, &
172  -2,0,1,-2,1,9*0,-1,2,-1,0,1,-2,1,10*0,1,-1,2*0,-1,1,6*0,-1,1,2*0, &
173  2,-2,2*0,-1,1/
174 ! ----------------------------------------------------------------------
175 
176 d1d2 = d1 * d2
177 DO i = 1, 4
178  x(i) = y(i)
179  x(i+4) = y1(i)*d1
180  x(i+8) = y2(i)*d2
181  x(i+12) = y12(i)*d1d2
182 END DO
183 DO i = 1, 16
184  xx = 0.0
185  DO k = 1, 16
186  xx = xx + wt(i,k) * x(k)
187  END DO
188  cl(i) = xx
189 END DO
190 l = 0
191 DO i = 1, 4
192  DO j = 1, 4
193  l = l + 1
194  c(i,j) = cl(l)
195  END DO
196 END DO
197 
198 END SUBROUTINE bcucof
199 
200 
201 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
202 !
203 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
204 !
205  SUBROUTINE bi_cub_spline ( rho_rdf, xg, ya, x1a, x2a, y1a, y2a, y12a, &
206  c_bicub, rdf, dg_drho, dg_dr, i_max, ih, k )
207 !
208  USE dft_module, ONLY: ndft, r_grid
209  IMPLICIT NONE
210 !
211 ! ----------------------------------------------------------------------
212  REAL, INTENT(IN OUT) :: rho_rdf
213  REAL, INTENT(IN OUT) :: xg
214  REAL, INTENT(IN) :: ya(r_grid,ndft)
215  REAL, INTENT(IN) :: x1a(r_grid)
216  REAL, INTENT(IN) :: x2a(ndft)
217  REAL, INTENT(IN) :: y1a(r_grid,ndft)
218  REAL, INTENT(IN) :: y2a(r_grid,ndft)
219  REAL, INTENT(IN) :: y12a(r_grid,ndft)
220  REAL, INTENT(IN) :: c_bicub(r_grid,ndft,4,4)
221  REAL, INTENT(OUT) :: rdf
222  REAL, INTENT(OUT) :: dg_drho
223  REAL, INTENT(OUT) :: dg_dr
224  INTEGER, INTENT(IN OUT) :: i_max
225  !INTEGER, INTENT(IN OUT) :: k_max
226  INTEGER, INTENT(OUT) :: ih
227  INTEGER, INTENT(IN) :: k
228 !
229 ! ----------------------------------------------------------------------
230  INTEGER :: m, n
231 
232  REAL :: y(4),y1(4),y2(4),y12(4),x1l,x1u,x2l,x2u
233  REAL :: c(4,4)
234 ! ----------------------------------------------------------------------
235 
236 IF ( rho_rdf < x1a(1) ) THEN
237  dg_drho = 0.0
238  dg_dr = 0.0
239  rdf = 1.0
240  RETURN
241 END IF
242 IF ( x1a(ih) <= rho_rdf .AND. rho_rdf < x1a(ih+1) ) GO TO 10
243 IF ( ih > 2 ) THEN
244  IF ( x1a(ih-1) <= rho_rdf .AND. rho_rdf < x1a(ih) ) THEN
245  ih = ih - 1
246  GO TO 10
247  END IF
248 END IF
249 ! write (*,*) 'in ',ih
250 CALL hunt ( x1a, i_max, rho_rdf, ih )
251 ! write (*,*) 'out',ih
252 10 CONTINUE
253 IF ( x2a(k) /= xg ) THEN
254 ! write (*,*) 'error bi-cubic-spline',k,x2a(k),xg
255 ! DO k=1,NDFT
256 ! write (*,*) k,x2a(k)
257 ! ENDDO
258 ! stop
259 END IF
260 
261 
262 
263 y(1) = ya(ih,k)
264 y(2) = ya(ih+1,k)
265 y(3) = ya(ih+1,k+1)
266 y(4) = ya(ih,k+1)
267 
268 y1(1) = y1a(ih,k)
269 y1(2) = y1a(ih+1,k)
270 y1(3) = y1a(ih+1,k+1)
271 y1(4) = y1a(ih,k+1)
272 
273 y2(1) = y2a(ih,k)
274 y2(2) = y2a(ih+1,k)
275 y2(3) = y2a(ih+1,k+1)
276 y2(4) = y2a(ih,k+1)
277 
278 y12(1) = y12a(ih,k)
279 y12(2) = y12a(ih+1,k)
280 y12(3) = y12a(ih+1,k+1)
281 y12(4) = y12a(ih,k+1)
282 
283 x1l = x1a(ih)
284 x1u = x1a(ih+1)
285 x2l = x2a(k)
286 x2u = x2a(k+1)
287 
288 DO m = 1, 4
289  DO n = 1, 4
290  c(m,n) = c_bicub( ih, k, m, n )
291  END DO
292 END DO
293 CALL bcuint ( x1l, x1u, x2l, x2u, rho_rdf, xg, c, rdf, dg_drho, dg_dr )
294 
295 END SUBROUTINE bi_cub_spline
296 
297 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
298 ! SUBROUTINE hunt
299 !
300 ! Given an array xx(1:n), and given a value x, returns a value jlo
301 ! such that x is between xx(jlo) and xx(jlo+1). xx(1:n) must be
302 ! monotonic, either increasing or decreasing. jlo=0 or jlo=n is
303 ! returned to indicate that x is out of range. jlo on input is taken
304 ! as the initial guess for jlo on output.
305 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
306 !
307  SUBROUTINE hunt ( xx, n, x, jlo )
308 !
309  IMPLICIT NONE
310 !
311 ! ----------------------------------------------------------------------
312  INTEGER, INTENT(IN) :: n
313  INTEGER, INTENT(OUT) :: jlo
314  REAL, INTENT(IN) :: xx(n)
315  REAL :: x
316 !
317 ! ----------------------------------------------------------------------
318  INTEGER :: inc,jhi,jm
319  LOGICAL :: ascnd
320 ! ----------------------------------------------------------------------
321 
322 ascnd = xx(n) >= xx(1)
323 IF( jlo <= 0 .OR. jlo > n ) THEN
324  jlo = 0
325  jhi = n + 1
326  GO TO 3
327 END IF
328 inc = 1
329 IF( x >= xx(jlo) .EQV. ascnd ) THEN
330 1 jhi = jlo + inc
331  IF ( jhi > n ) THEN
332  jhi = n + 1
333  ELSE IF ( x >= xx(jhi) .EQV. ascnd ) THEN
334  jlo = jhi
335  inc = inc + inc
336  GO TO 1
337  END IF
338 ELSE
339  jhi = jlo
340 2 jlo = jhi - inc
341  IF ( jlo < 1 ) THEN
342  jlo = 0
343  ELSE IF ( x < xx(jlo) .EQV. ascnd ) THEN
344  jhi = jlo
345  inc = inc + inc
346  GO TO 2
347  END IF
348 END IF
349 3 IF (jhi-jlo == 1 ) THEN
350  IF ( x == xx(n)) jlo = n - 1
351  IF ( x == xx(1) ) jlo = 1
352  RETURN
353 END IF
354 jm = ( jhi + jlo ) / 2
355 IF ( x >= xx(jm) .EQV. ascnd ) THEN
356  jlo = jm
357 ELSE
358  jhi = jm
359 END IF
360 GO TO 3
361 END SUBROUTINE hunt
362 
363 
364 
365 !**********************************************************************
366 !
367 !**********************************************************************
368 !
369  !SUBROUTINE bcuint ( y, y1, y2, y12, x1l, x1u, x2l, x2u, x1, x2, c, ansy, ansy1, ansy2 )
370  SUBROUTINE bcuint ( x1l, x1u, x2l, x2u, x1, x2, c, ansy, ansy1, ansy2 )
371 !
372  IMPLICIT NONE
373 !
374 ! ----------------------------------------------------------------------
375  !REAL, INTENT(IN OUT) :: y(4)
376  !REAL, INTENT(IN OUT) :: y1(4)
377  !REAL, INTENT(IN OUT) :: y2(4)
378  !REAL, INTENT(IN OUT) :: y12(4)
379  REAL, INTENT(IN OUT) :: x1l
380  REAL, INTENT(IN OUT) :: x1u
381  REAL, INTENT(IN OUT) :: x2l
382  REAL, INTENT(IN OUT) :: x2u
383  REAL, INTENT(IN OUT) :: x1
384  REAL, INTENT(IN OUT) :: x2
385  REAL, INTENT(IN) :: c(4,4)
386  REAL, INTENT(OUT) :: ansy
387  REAL, INTENT(OUT) :: ansy1
388  REAL, INTENT(OUT) :: ansy2
389 !
390 ! ----------------------------------------------------------------------
391  !U USES bcucof
392  INTEGER :: i
393  REAL :: t, u
394 ! ----------------------------------------------------------------------
395 
396 ! call bcucof ( y, y1, y2, y12, x1u-x1l, x2u-x2l, c )
397 
398 IF ( x1u == x1l .OR. x2u == x2l ) pause 'bad input in bcuint'
399 t = (x1-x1l) / (x1u-x1l)
400 u = (x2-x2l) / (x2u-x2l)
401 ansy = 0.0
402 ansy2 = 0.0
403 ansy1 = 0.0
404 DO i = 4, 1, -1
405  ansy = t *ansy + ( (c(i,4)*u + c(i,3))*u+c(i,2) )*u + c(i,1)
406  ansy2 = t *ansy2 + ( 3.*c(i,4)*u+2.*c(i,3) )*u + c(i,2)
407  ansy1 = u *ansy1 + ( 3.*c(4,i)*t+2.*c(3,i) )*t + c(2,i)
408 END DO
409 ansy1 = ansy1 / (x1u-x1l)
410 ansy2 = ansy2 / (x2u-x2l)
411 
412 END SUBROUTINE bcuint
413 
414 
415 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
416 ! SUBROUTINE spline
417 !
418 ! Given arrays x(1:n) and y(1:n) containing a tabulated function,
419 ! i.e., yi = f(xi), with x1 < x2 < .. . < xN, and given values yp1 and
420 ! ypn for the first derivative of the interpolating function at points 1
421 ! and n, respectively, this routine returns an array y2(1:n) of length n
422 ! which contains the second derivatives of the interpolating function at
423 ! the tabulated points xi. If yp1 and/or ypn are equal to 1 1030 or
424 ! larger, the routine is signaled to set the corresponding boundary
425 ! condition for a natural spline, with zero second derivative on that
426 ! boundary.
427 ! Parameter: NMAX is the largest anticipated value of n.
428 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
429 !
430  SUBROUTINE spline ( x, y, n, yp1, ypn, y2 )
431 !
432  IMPLICIT NONE
433 !
434 ! ----------------------------------------------------------------------
435  INTEGER, INTENT(IN) :: n
436  REAL, INTENT(IN) :: x(n)
437  REAL, INTENT(IN) :: y(n)
438  REAL, INTENT(IN) :: yp1
439  REAL, INTENT(IN) :: ypn
440  REAL, INTENT(OUT) :: y2(n)
441 !
442 ! ----------------------------------------------------------------------
443  INTEGER, PARAMETER :: nmax = 1000
444  INTEGER :: i, k
445  REAL :: p, qn, sig, un, u(nmax)
446 ! ----------------------------------------------------------------------
447 
448 If(n > nmax) stop 'Increase NMAX in Spline!!'
449 
450 
451 
452  IF ( yp1 > 0.99e30 ) THEN
453  y2(1) = 0.0
454  u(1) = 0.0
455  ELSE
456  y2(1) = -0.5
457  u(1) = ( 3.0/(x(2)-x(1)) ) * ( (y(2)-y(1))/(x(2)-x(1))-yp1 )
458  END IF
459  DO i = 2, n-1
460  IF ( (x(i+1)-x(i)) == 0.0 .OR. (x(i)-x(i-1)) == 0.0 .OR. (x(i+1)-x(i-1)) == 0.0 ) THEN
461  write (*,*) 'Surface Tension Code: error in spline-interpolation'
462  stop 5
463  END IF
464  sig = (x(i)-x(i-1)) / (x(i+1)-x(i-1))
465  p = sig*y2(i-1) + 2.0
466  y2(i) = (sig-1.0) / p
467  u(i) = ( 6.0 * ((y(i+1)-y(i))/(x(i+1)-x(i))-(y(i)-y(i-1))/(x(i)-x(i-1))) / (x(i+1)-x(i-1)) &
468  - sig * u(i-1) ) / p
469  END DO
470  IF ( ypn > 0.99e30 ) THEN
471  qn = 0.0
472  un = 0.0
473  ELSE
474  qn = 0.5
475  un = (3.0/(x(n)-x(n-1))) * (ypn-(y(n)-y(n-1))/(x(n)-x(n-1)))
476  END IF
477  y2(n) = (un-qn*u(n-1)) / (qn*y2(n-1)+1.0)
478  DO k = n-1, 1, -1
479  y2(k) = y2(k) * y2(k+1) + u(k)
480  END DO
481 
482 END SUBROUTINE spline
483 
484 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
485 ! SUBROUTINE splint_integral
486 !
487 ! Given the arrays xa(1:n) and ya(1:n) of length n, which tabulate a
488 ! function (with the in order), and given the array y2a(1:n), which is
489 ! the output from spline above, and given a value of x, this routine
490 ! returns a cubic-spline interpolated value y.
491 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
492 !
493  SUBROUTINE splint_integral ( xa, ya, y2a, n, xlo, xhi, integral )
494 !
495  IMPLICIT NONE
496 !
497 ! ----------------------------------------------------------------------
498  INTEGER, INTENT(IN) :: n
499  REAL, INTENT(IN) :: xa(n)
500  REAL, INTENT(IN) :: ya(n)
501  REAL, INTENT(IN) :: y2a(n)
502  REAL, INTENT(IN) :: xlo
503  REAL, INTENT(IN) :: xhi
504  REAL, INTENT(OUT) :: integral
505 !
506 ! ----------------------------------------------------------------------
507  INTEGER :: k, khi_l, klo_l, khi_h, klo_h
508  REAL :: xl, xh, h, int, x0, x1, y0, y1, y20, y21
509 ! ----------------------------------------------------------------------
510 
511  integral = 0.0
512  klo_l = 1
513  khi_l = n
514 1 IF ( khi_l-klo_l > 1 ) THEN
515  k = ( khi_l + klo_l ) / 2
516  IF ( xa(k) > xlo ) THEN
517  khi_l = k
518  ELSE
519  klo_l = k
520  END IF
521  GO TO 1
522  END IF
523 
524  klo_h = 1
525  khi_h = n
526 2 IF ( khi_h-klo_h > 1 ) THEN
527  k = ( khi_h + klo_h ) / 2
528  IF ( xa(k) > xhi ) THEN
529  khi_h = k
530  ELSE
531  klo_h = k
532  END IF
533  GO TO 2
534  END IF
535 
536  ! integration in spline pieces, the lower interval, bracketed
537  ! by xa(klo_L) and xa(khi_L) is in steps shifted upward.
538 
539  ! first: determine upper integration bound
540  xl = xlo
541 3 CONTINUE
542  IF ( khi_h > khi_l ) THEN
543  xh = xa(khi_l)
544  ELSE IF ( khi_h == khi_l ) THEN
545  xh = xhi
546  ELSE
547  WRITE (*,*) 'error in spline-integration'
548  pause
549  END IF
550 
551  h = xa(khi_l) - xa(klo_l)
552  IF ( h == 0.0 ) pause 'bad xa input in splint'
553  x0 = xa(klo_l)
554  x1 = xa(khi_l)
555  y0 = ya(klo_l)
556  y1 = ya(khi_l)
557  y20= y2a(klo_l)
558  y21= y2a(khi_l)
559  ! int = -xL/h * ( (x1-.5*xL)*y0 + (0.5*xL-x0)*y1 &
560  ! +y20/6.*(x1**3-1.5*xL*x1*x1+xL*xL*x1-.25*xL**3) &
561  ! -y20/6.*h*h*(x1-.5*xL) &
562  ! +y21/6.*(.25*xL**3-xL*xL*x0+1.5*xL*x0*x0-x0**3) &
563  ! -y21/6.*h*h*(.5*xL-x0) )
564  ! int = int + xH/h * ( (x1-.5*xH)*y0 + (0.5*xH-x0)*y1 &
565  ! +y20/6.*(x1**3-1.5*xH*x1*x1+xH*xH*x1-.25*xH**3) &
566  ! -y20/6.*h*h*(x1-.5*xH) &
567  ! +y21/6.*(.25*xH**3-xH*xH*x0+1.5*xH*x0*x0-x0**3) &
568  ! -y21/6.*h*h*(.5*xH-x0) )
569  int = -1.0/h * ( xl*((x1-.5*xl)*y0 + (0.5*xl-x0)*y1) &
570  -y20/24.*(x1-xl)**4 + y20/6.*(0.5*xl*xl-x1*xl)*h*h &
571  +y21/24.*(xl-x0)**4 - y21/6.*(0.5*xl*xl-x0*xl)*h*h )
572  int = int + 1.0/h * ( xh*((x1-.5*xh)*y0 + (0.5*xh-x0)*y1) &
573  -y20/24.*(x1-xh)**4 + y20/6.*(0.5*xh*xh-x1*xh)*h*h &
574  +y21/24.*(xh-x0)**4 - y21/6.*(0.5*xh*xh-x0*xh)*h*h )
575 
576  integral = integral + int
577  ! write (*,*) integral,x1,xH
578  klo_l = klo_l + 1
579  khi_l = khi_l + 1
580  xl = x1
581  IF (khi_h /= (khi_l-1)) GO TO 3 ! the -1 in (khi_L-1) because khi_L was already counted up
582 
583 END SUBROUTINE splint_integral
584 
585 
586 
587 
588 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
589 !
590 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
591 !
592  REAL FUNCTION praxis( t0, machep, h0, n, prin, x, f, fmin )
593 
594  IMPLICIT NONE
595 
596 ! ----------------------------------------------------------------------
597  REAL, INTENT(IN OUT) :: t0
598  REAL, INTENT(IN) :: machep
599  REAL, INTENT(IN) :: h0
600  INTEGER :: n
601  INTEGER, INTENT(IN OUT) :: prin
602  REAL, INTENT(IN OUT) :: x(n)
603  REAL :: f
604  REAL, INTENT(IN OUT) :: fmin
605 ! ----------------------------------------------------------------------
606 
607 EXTERNAL f
608 
609 ! PRAXIS RETURNS THE MINIMUM OF THE FUNCTION F(X,N) OF N VARIABLES
610 ! USING THE PRINCIPAL AXIS METHOD. THE GRADIENT OF THE FUNCTION IS
611 ! NOT REQUIRED.
612 
613 ! FOR A DESCRIPTION OF THE ALGORITHM, SEE CHAPTER SEVEN OF
614 ! "ALGORITHMS FOR FINDING ZEROS AND EXTREMA OF FUNCTIONS WITHOUT
615 ! CALCULATING DERIVATIVES" BY RICHARD P BRENT.
616 
617 ! THE PARAMETERS ARE:
618 ! T0 IS A TOLERANCE. PRAXIS ATTEMPTS TO RETURN PRAXIS=F(X)
619 ! SUCH THAT IF X0 IS THE TRUE LOCAL MINIMUM NEAR X, THEN
620 ! NORM(X-X0) < T0 + SQUAREROOT(MACHEP)*NORM(X).
621 ! MACHEP IS THE MACHINE PRECISION, THE SMALLEST NUMBER SUCH THAT
622 ! 1 + MACHEP > 1. MACHEP SHOULD BE 16.**-13 (ABOUT
623 ! 2.22D-16) FOR REAL*8 ARITHMETIC ON THE IBM 360.
624 ! H0 IS THE MAXIMUM STEP SIZE. H0 SHOULD BE SET TO ABOUT THE
625 ! MAXIMUM DISTANCE FROM THE INITIAL GUESS TO THE MINIMUM.
626 ! (IF H0 IS SET TOO LARGE OR TOO SMALL, THE INITIAL RATE OF
627 ! CONVERGENCE MAY BE SLOW.)
628 ! N (AT LEAST TWO) IS THE NUMBER OF VARIABLES UPON WHICH
629 ! THE FUNCTION DEPENDS.
630 ! PRIN CONTROLS THE PRINTING OF INTERMEDIATE RESULTS.
631 ! IF PRIN=0, NOTHING IS PRINTED.
632 ! IF PRIN=1, F IS PRINTED AFTER EVERY N+1 OR N+2 LINEAR
633 ! MINIMIZATIONS. FINAL X IS PRINTED, BUT INTERMEDIATE X IS
634 ! PRINTED ONLY IF N IS AT MOST 4.
635 ! IF PRIN=2, THE SCALE FACTORS AND THE PRINCIPAL VALUES OF
636 ! THE APPROXIMATING QUADRATIC FORM ARE ALSO PRINTED.
637 ! IF PRIN=3, X IS ALSO PRINTED AFTER EVERY FEW LINEAR
638 ! MINIMIZATIONS.
639 ! IF PRIN=4, THE PRINCIPAL VECTORS OF THE APPROXIMATING
640 ! QUADRATIC FORM ARE ALSO PRINTED.
641 ! X IS AN ARRAY CONTAINING ON ENTRY A GUESS OF THE POINT OF
642 ! MINIMUM, ON RETURN THE ESTIMATED POINT OF MINIMUM.
643 ! F(X,N) IS THE FUNCTION TO BE MINIMIZED. F SHOULD BE A REAL*8
644 ! FUNCTION DECLARED EXTERNAL IN THE CALLING PROGRAM.
645 ! FMIN IS AN ESTIMATE OF THE MINIMUM, USED ONLY IN PRINTING
646 ! INTERMEDIATE RESULTS.
647 ! THE APPROXIMATING QUADRATIC FORM IS
648 ! Q(X') = F(X,N) + (1/2) * (X'-X)-TRANSPOSE * A * (X'-X)
649 ! WHERE X IS THE BEST ESTIMATE OF THE MINIMUM AND A IS
650 ! INVERSE(V-TRANSPOSE) * D * INVERSE(V)
651 ! (V(*,*) IS THE MATRIX OF SEARCH DIRECTIONS; D(*) IS THE ARRAY
652 ! OF SECOND DIFFERENCES). IF F HAS CONTINUOUS SECOND DERIVATIVES
653 ! NEAR X0, A WILL TEND TO THE HESSIAN OF F AT X0 AS X APPROACHES X0.
654 
655 ! IT IS ASSUMED THAT ON FLOATING-POINT UNDERFLOW THE RESULT IS SET
656 ! TO ZERO.
657 ! THE USER SHOULD OBSERVE THE COMMENT ON HEURISTIC NUMBERS AFTER
658 ! THE INITIALIZATION OF MACHINE DEPENDENT NUMBERS.
659 
660  LOGICAL :: illc
661  INTEGER :: nl,nf,kl,kt,ktm,idim,i,j,k,k2,km1,klmk,ii,im1
662  REAL :: s,sl,dn,dmin,fx,f1,lds,ldt,t,h,sf,df,qf1,qd0, qd1,qa,qb,qc
663  REAL :: m2,m4,small,vsmall,large,vlarge,scbd,ldfac,t2, dni,value
664  REAL :: random
665 
666 !.....IF N>20 OR IF N<20 AND YOU NEED MORE SPACE, CHANGE '20' TO THE
667 ! LARGEST VALUE OF N IN THE NEXT CARD, IN THE CARD 'IDIM=20', AND
668 ! IN THE DIMENSION STATEMENTS IN SUBROUTINES MINFIT,MIN,FLIN,QUAD.
669 
670  REAL :: d(20),y(20),z(20),q0(20),q1(20),v(20,20)
671  COMMON /global/ fx,ldt,dmin,nf,nl /q/ v,q0,q1,qa,qb,qc,qd0,qd1,qf1
672 
673 ! ---------------------------------
674 ! introduced by Joachim........
675  idim = n
676 ! ---------------------------------
677 
678 
679 
680 !.....INITIALIZATION.....
681 ! MACHINE DEPENDENT NUMBERS:
682 
683 small = machep*machep
684 vsmall = small*small
685 large = 1.d0/small
686 vlarge = 1.d0/vsmall
687 m2 = sqrt(machep)
688 m4 = sqrt(m2)
689 
690 ! HEURISTIC NUMBERS:
691 ! IF THE AXES MAY BE BADLY SCALED (WHICH IS TO BE AVOIDED IF
692 ! POSSIBLE), THEN SET SCBD=10. OTHERWISE SET SCBD=1.
693 ! IF THE PROBLEM IS KNOWN TO BE ILL-CONDITIONED, SET ILLC=TRUE.
694 ! OTHERWISE SET ILLC=FALSE.
695 ! KTM IS THE NUMBER OF ITERATIONS WITHOUT IMPROVEMENT BEFORE THE
696 ! ALGORITHM TERMINATES. KTM=4 IS VERY CAUTIOUS; USUALLY KTM=1
697 ! IS SATISFACTORY.
698 
699 scbd = 1.0
700 illc = .false.
701 ktm = 1
702 
703 ldfac = 0.01
704 IF (illc) ldfac = 0.1
705 kt = 0
706 nl = 0
707 nf = 1
708 fx = f(x,n)
709 qf1 = fx
710 t = small+abs(t0)
711 t2 = t
712 dmin = small
713 h = h0
714 IF (h < 100*t) h = 100*t
715 ldt = h
716 !.....THE FIRST SET OF SEARCH DIRECTIONS V IS THE IDENTITY MATRIX.....
717 DO i = 1,n
718  DO j = 1,n
719  v(i,j) = 0.0
720  END DO
721  v(i,i) = 1.0
722 END DO
723 d(1) = 0.0
724 qd0 = 0.0
725 DO i = 1,n
726  q0(i) = x(i)
727  q1(i) = x(i)
728 END DO
729 IF (prin > 0) CALL print(n,x,prin,fmin)
730 
731 !.....THE MAIN LOOP STARTS HERE.....
732 40 sf=d(1)
733 d(1)=0.d0
734 s=0.d0
735 
736 !.....MINIMIZE ALONG THE FIRST DIRECTION V(*,1).
737 ! FX MUST BE PASSED TO MIN BY VALUE.
738 value=fx
739 CALL min(n,1,2,d(1),s,value,.false.,f,x,t,machep,h)
740 IF (s > 0.d0) GO TO 50
741 DO i=1,n
742  v(i,1)=-v(i,1)
743 END DO
744 50 IF (sf > 0.9d0*d(1).AND.0.9d0*sf < d(1)) GO TO 70
745 DO i=2,n
746  d(i)=0.d0
747 END DO
748 
749 !.....THE INNER LOOP STARTS HERE.....
750 70 DO k=2,n
751  DO i=1,n
752  y(i)=x(i)
753  END DO
754  sf=fx
755  IF (kt > 0) illc=.true.
756  80 kl=k
757  df=0.d0
758 
759 !.....A RANDOM STEP FOLLOWS (TO AVOID RESOLUTION VALLEYS).
760 ! PRAXIS ASSUMES THAT RANDOM RETURNS A RANDOM NUMBER UNIFORMLY
761 ! DISTRIBUTED IN (0,1).
762 
763  IF(.NOT.illc) GO TO 95
764  DO i=1,n
765  s=(0.1d0*ldt+t2*(10**kt))*(random(n)-0.5d0)
766  z(i)=s
767  DO j=1,n
768  x(j)=x(j)+s*v(j,i)
769  END DO
770  END DO
771  fx=f(x,n)
772  nf=nf+1
773 
774 !.....MINIMIZE ALONG THE "NON-CONJUGATE" DIRECTIONS V(*,K),...,V(*,N)
775 
776  95 DO k2=k,n
777  sl=fx
778  s=0.d0
779  value=fx
780  CALL min(n,k2,2,d(k2),s,value,.false.,f,x,t,machep,h)
781  IF (illc) GO TO 97
782  s=sl-fx
783  GO TO 99
784  97 s=d(k2)*((s+z(k2))**2)
785  99 IF (df > s) cycle
786  df=s
787  kl=k2
788  END DO
789  IF (illc.OR.(df >= abs((100*machep)*fx))) GO TO 110
790 
791 !.....IF THERE WAS NOT MUCH IMPROVEMENT ON THE FIRST TRY, SET
792 ! ILLC=TRUE AND START THE INNER LOOP AGAIN.....
793 
794  illc=.true.
795  GO TO 80
796  110 IF (k == 2.AND.prin > 1) CALL vcprnt(1,d,n)
797 
798 !.....MINIMIZE ALONG THE "CONJUGATE" DIRECTIONS V(*,1),...,V(*,K-1)
799 
800  km1=k-1
801  DO k2=1,km1
802  s=0
803  value=fx
804  CALL min(n,k2,2,d(k2),s,value,.false.,f,x,t,machep,h)
805  END DO
806  f1=fx
807  fx=sf
808  lds=0
809  DO i=1,n
810  sl=x(i)
811  x(i)=y(i)
812  sl=sl-y(i)
813  y(i)=sl
814  lds=lds+sl*sl
815  END DO
816  lds=sqrt(lds)
817  IF (lds <= small) GO TO 160
818 
819 !.....DISCARD DIRECTION V(*,KL).
820 ! IF NO RANDOM STEP WAS TAKEN, V(*,KL) IS THE "NON-CONJUGATE"
821 ! DIRECTION ALONG WHICH THE GREATEST IMPROVEMENT WAS MADE.....
822 
823  klmk=kl-k
824  IF (klmk < 1) GO TO 141
825  DO ii=1,klmk
826  i=kl-ii
827  DO j=1,n
828  v(j,i+1)=v(j,i)
829  END DO
830  d(i+1)=d(i)
831  END DO
832  141 d(k)=0
833  DO i=1,n
834  v(i,k)=y(i)/lds
835  END DO
836 
837 !.....MINIMIZE ALONG THE NEW "CONJUGATE" DIRECTION V(*,K), WHICH IS
838 ! THE NORMALIZED VECTOR: (NEW X) - (0LD X).....
839 
840  value=f1
841  CALL min(n,k,4,d(k),lds,value,.true.,f,x,t,machep,h)
842  IF (lds > 0.d0) GO TO 160
843  lds=-lds
844  DO i=1,n
845  v(i,k)=-v(i,k)
846  END DO
847  160 ldt=ldfac*ldt
848  IF (ldt < lds) ldt=lds
849  IF (prin > 0) CALL print(n,x,prin,fmin)
850  t2=0.d0
851  DO i=1,n
852  t2=t2+x(i)**2
853  END DO
854  t2=m2*sqrt(t2)+t
855 
856 !.....SEE WHETHER THE LENGTH OF THE STEP TAKEN SINCE STARTING THE
857 ! INNER LOOP EXCEEDS HALF THE TOLERANCE.....
858 
859  IF (ldt > (0.5*t2)) kt=-1
860  kt=kt+1
861  IF (kt > ktm) GO TO 400
862 END DO
863 !.....THE INNER LOOP ENDS HERE.
864 
865 ! TRY QUADRATIC EXTRAPOLATION IN CASE WE ARE IN A CURVED VALLEY.
866 
867 CALL quad(n,f,x,t,machep,h)
868 dn=0.d0
869 DO i=1,n
870  d(i)=1.d0/sqrt(d(i))
871  IF (dn < d(i)) dn=d(i)
872 END DO
873 IF (prin > 3) CALL maprnt(1,v,idim,n)
874 DO j=1,n
875  s=d(j)/dn
876  DO i=1,n
877  v(i,j)=s*v(i,j)
878  END DO
879 END DO
880 
881 !.....SCALE THE AXES TO TRY TO REDUCE THE CONDITION NUMBER.....
882 
883 IF (scbd <= 1.d0) GO TO 200
884 s=vlarge
885 DO i=1,n
886  sl=0.d0
887  DO j=1,n
888  sl=sl+v(i,j)*v(i,j)
889  END DO
890  z(i)=sqrt(sl)
891  IF (z(i) < m4) z(i)=m4
892  IF (s > z(i)) s=z(i)
893 END DO
894 DO i=1,n
895  sl=s/z(i)
896  z(i)=1.d0/sl
897  IF (z(i) <= scbd) GO TO 189
898  sl=1.d0/scbd
899  z(i)=scbd
900  189 DO j=1,n
901  v(i,j)=sl*v(i,j)
902  END DO
903 END DO
904 
905 !.....CALCULATE A NEW SET OF ORTHOGONAL DIRECTIONS BEFORE REPEATING
906 ! THE MAIN LOOP.
907 ! FIRST TRANSPOSE V FOR MINFIT:
908 
909 200 DO i=2,n
910  im1=i-1
911  DO j=1,im1
912  s=v(i,j)
913  v(i,j)=v(j,i)
914  v(j,i)=s
915  END DO
916 END DO
917 
918 !.....CALL MINFIT TO FIND THE SINGULAR VALUE DECOMPOSITION OF V.
919 ! THIS GIVES THE PRINCIPAL VALUES AND PRINCIPAL DIRECTIONS OF THE
920 ! APPROXIMATING QUADRATIC FORM WITHOUT SQUARING THE CONDITION
921 ! NUMBER.....
922 
923 CALL minfit(idim,n,machep,vsmall,v,d)
924 
925 !.....UNSCALE THE AXES.....
926 
927 IF (scbd <= 1.d0) GO TO 250
928 DO i=1,n
929  s=z(i)
930  DO j=1,n
931  v(i,j)=s*v(i,j)
932  END DO
933 END DO
934 DO i=1,n
935  s=0.d0
936  DO j=1,n
937  s=s+v(j,i)**2
938  END DO
939  s=sqrt(s)
940  d(i)=s*d(i)
941  s=1/s
942  DO j=1,n
943  v(j,i)=s*v(j,i)
944  END DO
945 END DO
946 
947 250 DO i=1,n
948  dni=dn*d(i)
949  IF (dni > large) GO TO 265
950  IF (dni < small) GO TO 260
951  d(i)=1/(dni*dni)
952  cycle
953  260 d(i)=vlarge
954  cycle
955  265 d(i)=vsmall
956 END DO
957 
958 !.....SORT THE EIGENVALUES AND EIGENVECTORS.....
959 
960 CALL sort(idim,n,d,v)
961 dmin=d(n)
962 IF (dmin < small) dmin=small
963 illc=.false.
964 IF (m2*d(1) > dmin) illc=.true.
965 IF (prin > 1.AND.scbd > 1.d0) CALL vcprnt(2,z,n)
966 IF (prin > 1) CALL vcprnt(3,d,n)
967 IF (prin > 3) CALL maprnt(2,v,idim,n)
968 !.....THE MAIN LOOP ENDS HERE.....
969 
970 GO TO 40
971 
972 !.....RETURN.....
973 
974 400 IF (prin > 0) CALL vcprnt(4,x,n)
975 praxis=fx
976 
977 END FUNCTION praxis
978 
979 
980 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
981 !
982 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
983 !
984  SUBROUTINE minfit(m,n,machep,tol,ab,q)
985 
986  IMPLICIT NONE
987 
988  INTEGER, INTENT(IN OUT) :: m
989  INTEGER, INTENT(IN) :: n
990  REAL, INTENT(IN) :: machep
991  REAL, INTENT(IN OUT) :: tol
992  REAL, INTENT(IN OUT) :: ab(m,n)
993  REAL, INTENT(OUT) :: q(n)
994  INTEGER :: i,j,k,l, kk,kt,l2,ll2,ii,lp1
995 ! IMPLICIT REAL (A-H,O-Z)
996 
997 
998 REAL :: x,eps,e(20),g,s, f,h,y,c,z,temp
999 !...AN IMPROVED VERSION OF MINFIT (SEE GOLUB AND REINSCH, 1969)
1000 ! RESTRICTED TO M=N,P=0.
1001 ! THE SINGULAR VALUES OF THE ARRAY AB ARE RETURNED IN Q AND AB IS
1002 ! OVERWRITTEN WITH THE ORTHOGONAL MATRIX V SUCH THAT U.DIAG(Q) = AB.V,
1003 ! WHERE U IS ANOTHER ORTHOGONAL MATRIX.
1004 
1005 !...HOUSEHOLDER'S REDUCTION TO BIDIAGONAL FORM...
1006 IF (n == 1) GO TO 200
1007 eps = machep
1008 g = 0.d0
1009 x = 0.d0
1010 DO i=1,n
1011  e(i) = g
1012  s = 0.d0
1013  l = i + 1
1014  DO j=i,n
1015  s = s + ab(j,i)**2
1016  END DO
1017  g = 0.d0
1018  IF (s < tol) GO TO 4
1019  f = ab(i,i)
1020  g = sqrt(s)
1021  IF (f >= 0.d0) g = -g
1022  h = f*g - s
1023  ab(i,i)=f-g
1024  IF (l > n) GO TO 4
1025  DO j=l,n
1026  f = 0.d0
1027  DO k=i,n
1028  f = f + ab(k,i)*ab(k,j)
1029  END DO
1030  f = f/h
1031  DO k=i,n
1032  ab(k,j) = ab(k,j) + f*ab(k,i)
1033  END DO
1034  END DO
1035  4 q(i) = g
1036  s = 0.d0
1037  IF (i == n) GO TO 6
1038  DO j=l,n
1039  s = s + ab(i,j)*ab(i,j)
1040  END DO
1041  6 g = 0.d0
1042  IF (s < tol) GO TO 10
1043  IF (i == n) GO TO 16
1044  f = ab(i,i+1)
1045  16 g = sqrt(s)
1046  IF (f >= 0.d0) g = -g
1047  h = f*g - s
1048  IF (i == n) GO TO 10
1049  ab(i,i+1) = f - g
1050  DO j=l,n
1051  e(j) = ab(i,j)/h
1052  END DO
1053  DO j=l,n
1054  s = 0.d0
1055  DO k=l,n
1056  s = s + ab(j,k)*ab(i,k)
1057  END DO
1058  DO k=l,n
1059  ab(j,k) = ab(j,k) + s*e(k)
1060  END DO
1061  END DO
1062  10 y = abs(q(i)) + abs(e(i))
1063  IF (y > x) x = y
1064 END DO
1065 !...ACCUMULATION OF RIGHT-HAND TRANSFORMATIONS...
1066 ab(n,n) = 1.d0
1067 g = e(n)
1068 l = n
1069 DO ii=2,n
1070  i = n - ii + 1
1071  IF (g == 0.d0) GO TO 23
1072  h = ab(i,i+1)*g
1073  DO j=l,n
1074  ab(j,i) = ab(i,j)/h
1075  END DO
1076  DO j=l,n
1077  s = 0.d0
1078  DO k=l,n
1079  s = s + ab(i,k)*ab(k,j)
1080  END DO
1081  DO k=l,n
1082  ab(k,j) = ab(k,j) + s*ab(k,i)
1083  END DO
1084  END DO
1085  23 DO j=l,n
1086  ab(i,j) = 0.d0
1087  ab(j,i) = 0.d0
1088  END DO
1089  ab(i,i) = 1.d0
1090  g = e(i)
1091  l = i
1092 END DO
1093 !...DIAGONALIZATION OF THE BIDIAGONAL FORM...
1094 eps = eps*x
1095 DO kk=1,n
1096  k = n - kk + 1
1097  kt = 0
1098  101 kt = kt + 1
1099  IF (kt <= 30) GO TO 102
1100  e(k) = 0.d0
1101  WRITE (6,1000)
1102  1000 FORMAT (' QR FAILED')
1103  102 DO ll2=1,k
1104  l2 = k - ll2 + 1
1105  l = l2
1106  IF (abs(e(l)) <= eps) GO TO 120
1107  IF (l == 1) cycle
1108  IF (abs(q(l-1)) <= eps) EXIT
1109  END DO
1110 !...CANCELLATION OF E(L) IF L>1...
1111  c = 0.d0
1112  s = 1.d0
1113  DO i=l,k
1114  f = s*e(i)
1115  e(i) = c*e(i)
1116  IF (abs(f) <= eps) GO TO 120
1117  g = q(i)
1118 !...Q(I) = H = SQRT(G*G + F*F)...
1119  IF (abs(f) < abs(g)) GO TO 113
1120  IF (f == 0.0) THEN
1121  GO TO 111
1122  ELSE
1123  GO TO 112
1124  END IF
1125  111 h = 0.d0
1126  GO TO 114
1127  112 h = abs(f)*sqrt(1 + (g/f)**2)
1128  GO TO 114
1129  113 h = abs(g)*sqrt(1 + (f/g)**2)
1130  114 q(i) = h
1131  IF (h /= 0.d0) GO TO 115
1132  g = 1.d0
1133  h = 1.d0
1134  115 c = g/h
1135  s = -f/h
1136  END DO
1137 !...TEST FOR CONVERGENCE...
1138  120 z = q(k)
1139  IF (l == k) GO TO 140
1140 !...SHIFT FROM BOTTOM 2*2 MINOR...
1141  x = q(l)
1142  y = q(k-1)
1143  g = e(k-1)
1144  h = e(k)
1145  f = ((y - z)*(y + z) + (g - h)*(g + h))/(2*h*y)
1146  g = sqrt(f*f + 1.0d0)
1147  temp = f - g
1148  IF (f >= 0.d0) temp = f + g
1149  f = ((x - z)*(x + z) + h*(y/temp - h))/x
1150 !...NEXT QR TRANSFORMATION...
1151  c = 1.d0
1152  s = 1.d0
1153  lp1 = l + 1
1154  IF (lp1 > k) GO TO 133
1155  DO i=lp1,k
1156  g = e(i)
1157  y = q(i)
1158  h = s*g
1159  g = g*c
1160  IF (abs(f) < abs(h)) GO TO 123
1161  IF (f == 0.0) THEN
1162  GO TO 121
1163  ELSE
1164  GO TO 122
1165  END IF
1166  121 z = 0.d0
1167  GO TO 124
1168  122 z = abs(f)*sqrt(1 + (h/f)**2)
1169  GO TO 124
1170  123 z = abs(h)*sqrt(1 + (f/h)**2)
1171  124 e(i-1) = z
1172  IF (z /= 0.d0) GO TO 125
1173  f = 1.d0
1174  z = 1.d0
1175  125 c = f/z
1176  s = h/z
1177  f = x*c + g*s
1178  g = -x*s + g*c
1179  h = y*s
1180  y = y*c
1181  DO j=1,n
1182  x = ab(j,i-1)
1183  z = ab(j,i)
1184  ab(j,i-1) = x*c + z*s
1185  ab(j,i) = -x*s + z*c
1186  END DO
1187  IF (abs(f) < abs(h)) GO TO 129
1188  IF (f == 0.0) THEN
1189  GO TO 127
1190  ELSE
1191  GO TO 128
1192  END IF
1193  127 z = 0.d0
1194  GO TO 130
1195  128 z = abs(f)*sqrt(1 + (h/f)**2)
1196  GO TO 130
1197  129 z = abs(h)*sqrt(1 + (f/h)**2)
1198  130 q(i-1) = z
1199  IF (z /= 0.d0) GO TO 131
1200  f = 1.d0
1201  z = 1.d0
1202  131 c = f/z
1203  s = h/z
1204  f = c*g + s*y
1205  x = -s*g + c*y
1206  END DO
1207  133 e(l) = 0.d0
1208  e(k) = f
1209  q(k) = x
1210  GO TO 101
1211 !...CONVERGENCE: Q(K) IS MADE NON-NEGATIVE...
1212  140 IF (z >= 0.d0) cycle
1213  q(k) = -z
1214  DO j=1,n
1215  ab(j,k) = -ab(j,k)
1216  END DO
1217 END DO
1218 RETURN
1219 200 q(1) = ab(1,1)
1220 ab(1,1) = 1.d0
1221 
1222 END SUBROUTINE minfit
1223 
1224 
1225 
1226 
1227 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1228 !
1229 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1230 !
1231  SUBROUTINE min(n,j,nits,d2,x1,f1,fk,f,x,t,machep,h)
1232 
1233  IMPLICIT NONE
1234 
1235  INTEGER, INTENT(IN) :: n
1236  INTEGER :: j
1237  INTEGER :: nits
1238  REAL, INTENT(IN OUT) :: d2
1239  REAL, INTENT(IN OUT) :: x1
1240  REAL, INTENT(IN OUT) :: f1
1241  LOGICAL :: fk
1242  REAL :: f
1243  REAL, INTENT(IN OUT) :: x(n)
1244  REAL, INTENT(IN) :: t
1245  REAL, INTENT(IN) :: machep
1246  REAL, INTENT(IN) :: h
1247 
1248  INTEGER :: i,k
1249  EXTERNAL f
1250 
1251 
1252  REAL :: flin ! function
1253  REAL :: small,sf1,sx1,s,temp, xm,x2,f2,d1
1254  REAL :: fm,f0,t2
1255 !----------------------------------------------
1256  INTEGER :: nf,nl
1257  REAL :: fx,ldt,dmin
1258  REAL :: v(20,20),q0(20),q1(20),qa,qb,qc,qd0,qd1,qf1
1259  COMMON /global/ fx,ldt,dmin,nf,nl /q/ v,q0,q1,qa,qb,qc,qd0,qd1,qf1
1260 !----------------------------------------------
1261 !...THE SUBROUTINE MIN MINIMIZES F FROM X IN THE DIRECTION V(*,J) UNLESS
1262 ! J IS LESS THAN 1, WHEN A QUADRATIC SEARCH IS MADE IN THE PLANE
1263 ! DEFINED BY Q0,Q1,X.
1264 ! D2 IS EITHER ZERO OR AN APPROXIMATION TO HALF F".
1265 ! ON ENTRY, X1 IS AN ESTIMATE OF THE DISTANCE FROM X TO THE MINIMUM
1266 ! ALONG V(*,J) (OR, IF J=0, A CURVE). ON RETURN, X1 IS THE DISTANCE
1267 ! FOUND.
1268 ! IF FK=.TRUE., THEN F1 IS FLIN(X1). OTHERWISE X1 AND F1 ARE IGNORED
1269 ! ON ENTRY UNLESS FINAL FX IS GREATER THAN F1.
1270 ! NITS CONTROLS THE NUMBER OF TIMES AN ATTEMPT WILL BE MADE TO HALVE
1271 ! THE INTERVAL.
1272  LOGICAL :: dz
1273  REAL :: m2,m4
1274 
1275 small = machep**2
1276 m2 = sqrt(machep)
1277 m4 = sqrt(m2)
1278 sf1 = f1
1279 sx1 = x1
1280 k = 0
1281 xm = 0.d0
1282 fm = fx
1283 f0 = fx
1284 dz = d2 < machep
1285 !...FIND THE STEP SIZE...
1286 s = 0.d0
1287 DO i=1,n
1288  s = s + x(i)**2
1289 END DO
1290 s = sqrt(s)
1291 temp = d2
1292 IF (dz) temp = dmin
1293 t2 = m4*sqrt(abs(fx)/temp + s*ldt) + m2*ldt
1294 s = m4*s + t
1295 IF (dz.AND.t2 > s) t2 = s
1296 t2 = dmax1(t2,small)
1297 t2 = dmin1(t2,.01d0*h)
1298 IF (.NOT.fk.OR.f1 > fm) GO TO 2
1299 xm = x1
1300 fm = f1
1301 2 IF (fk.AND.abs(x1) >= t2) GO TO 3
1302 temp=1.d0
1303 IF (x1 < 0.d0) temp=-1.d0
1304 x1=temp*t2
1305 f1 = flin(n,j,x1,f,x,nf)
1306 3 IF (f1 > fm) GO TO 4
1307 xm = x1
1308 fm = f1
1309 4 IF (.NOT.dz) GO TO 6
1310 !...EVALUATE FLIN AT ANOTHER POINT AND ESTIMATE THE SECOND DERIVATIVE...
1311 x2 = -x1
1312 IF (f0 >= f1) x2 = 2.d0*x1
1313 f2 = flin(n,j,x2,f,x,nf)
1314 IF (f2 > fm) GO TO 5
1315 xm = x2
1316 fm = f2
1317 5 d2 = (x2*(f1 - f0)-x1*(f2 - f0))/((x1*x2)*(x1 - x2))
1318 !...ESTIMATE THE FIRST DERIVATIVE AT 0...
1319 6 d1 = (f1 - f0)/x1 - x1*d2
1320 dz = .true.
1321 !...PREDICT THE MINIMUM...
1322 IF (d2 > small) GO TO 7
1323 x2 = h
1324 IF (d1 >= 0.d0) x2 = -x2
1325 GO TO 8
1326 7 x2 = (-.5d0*d1)/d2
1327 8 IF (abs(x2) <= h) GO TO 11
1328 IF (x2 > 0.0) THEN
1329  GO TO 10
1330 END IF
1331 x2 = -h
1332 GO TO 11
1333 10 x2 = h
1334 !...EVALUATE F AT THE PREDICTED MINIMUM...
1335 11 f2 = flin(n,j,x2,f,x,nf)
1336 IF (k >= nits.OR.f2 <= f0) GO TO 12
1337 !...NO SUCCESS, SO TRY AGAIN...
1338 k = k + 1
1339 IF (f0 < f1.AND.(x1*x2) > 0.d0) GO TO 4
1340 x2 = 0.5d0*x2
1341 GO TO 11
1342 !...INCREMENT THE ONE-DIMENSIONAL SEARCH COUNTER...
1343 12 nl = nl + 1
1344 IF (f2 <= fm) GO TO 13
1345 x2 = xm
1346 GO TO 14
1347 13 fm = f2
1348 !...GET A NEW ESTIMATE OF THE SECOND DERIVATIVE...
1349 14 IF (abs(x2*(x2 - x1)) <= small) GO TO 15
1350 d2 = (x2*(f1-f0) - x1*(fm-f0))/((x1*x2)*(x1 - x2))
1351 GO TO 16
1352 15 IF (k > 0) d2 = 0.d0
1353 16 IF (d2 <= small) d2 = small
1354 x1 = x2
1355 fx = fm
1356 IF (sf1 >= fx) GO TO 17
1357 fx = sf1
1358 x1 = sx1
1359 !...UPDATE X FOR LINEAR BUT NOT PARABOLIC SEARCH...
1360 17 IF (j == 0) RETURN
1361 DO i=1,n
1362  x(i) = x(i) + x1*v(i,j)
1363 END DO
1364 
1365 END SUBROUTINE min
1366 
1367 
1368 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1369 !
1370 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1371 !
1372  SUBROUTINE vcprnt(option,v,n)
1373 
1374  IMPLICIT NONE
1375  INTEGER :: option
1376  REAL, INTENT(IN OUT) :: v(n)
1377  INTEGER :: n
1378 
1379  INTEGER :: i
1380 
1381 SELECT CASE ( option )
1382  CASE ( 1)
1383  GO TO 1
1384  CASE ( 2)
1385  GO TO 2
1386  CASE ( 3)
1387  GO TO 3
1388  CASE ( 4)
1389  GO TO 4
1390 END SELECT
1391 
1392 1 WRITE (6,101) (v(i),i=1,n)
1393 RETURN
1394 2 WRITE (6,102) (v(i),i=1,n)
1395 RETURN
1396 3 WRITE (6,103) (v(i),i=1,n)
1397 RETURN
1398 4 WRITE (6,104) (v(i),i=1,n)
1399 RETURN
1400 101 FORMAT (/' THE SECOND DIFFERENCE ARRAY D(*) IS:'/ (e32.14,4e25.14))
1401 102 FORMAT (/' THE SCALE FACTORS ARE:'/(e32.14,4e25.14))
1402 103 FORMAT (/' THE APPROXIMATING QUADR. FORM HAS PRINCIPAL VALUES:'/ &
1403  (e32.14,4e25.14))
1404 104 FORMAT (/' X IS:',e26.14/(e32.14))
1405 END SUBROUTINE vcprnt
1406 
1407 
1408 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1409 !
1410 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1411 !
1412  SUBROUTINE print(n,x,prin,fmin)
1413 
1414  IMPLICIT NONE
1415  INTEGER :: n
1416  REAL, INTENT(IN OUT) :: x(n)
1417  INTEGER, INTENT(IN OUT) :: prin
1418  REAL, INTENT(IN OUT) :: fmin
1419 
1420  INTEGER :: i
1421  REAL :: ln
1422 !----------------------------------------------
1423 INTEGER :: nf,nl
1424 REAL :: fx,ldt,dmin
1425 COMMON /global/ fx,ldt,dmin,nf,nl
1426 !----------------------------------------------
1427 WRITE (6,101) nl,nf,fx
1428 
1429 IF (fx <= fmin) GO TO 1
1430 ln = log10(fx-fmin)
1431 WRITE (6,102) fmin,ln
1432 GO TO 2
1433 1 WRITE (6,103) fmin
1434 2 IF (n > 4.AND.prin <= 2) RETURN
1435 WRITE (6,104) (x(i),i=1,n)
1436 RETURN
1437 101 FORMAT (/' AFTER',i6, &
1438  ' LINEAR SEARCHES, THE FUNCTION HAS BEEN EVALUATED',i6, &
1439  ' TIMES. THE SMALLEST VALUE FOUND IS F(X) = ',e21.14)
1440 102 FORMAT (' LOG (F(X)-',e21.14,') = ',e21.14)
1441 103 FORMAT (' LOG (F(X)-',e21.14,') IS UNDEFINED.')
1442 104 FORMAT (' X IS:',e26.14/(e32.14))
1443 END SUBROUTINE print
1444 
1445 
1446 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1447 !
1448 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1449 !
1450  SUBROUTINE maprnt(option,v,m,n)
1451 
1452  IMPLICIT NONE
1453  INTEGER :: option
1454  REAL, INTENT(IN OUT) :: v(m,n)
1455  INTEGER, INTENT(IN OUT) :: m
1456  INTEGER, INTENT(IN) :: n
1457  INTEGER :: i,j
1458 
1459  INTEGER :: low,upp
1460 !...THE SUBROUTINE MAPRNT PRINTS THE COLUMNS OF THE NXN MATRIX V
1461 ! WITH A HEADING AS SPECIFIED BY OPTION.
1462 ! M IS THE ROW DIMENSION OF V AS DECLARED IN THE CALLING PROGRAM...
1463 low = 1
1464 upp = 5
1465 SELECT CASE ( option )
1466  CASE ( 1)
1467  GO TO 1
1468  CASE ( 2)
1469  GO TO 2
1470 END SELECT
1471 1 WRITE (6,101)
1472 101 FORMAT (/' THE NEW DIRECTIONS ARE:')
1473 GO TO 3
1474 2 WRITE (6,102)
1475 102 FORMAT (' AND THE PRINCIPAL AXES:')
1476 3 IF (n < upp) upp = n
1477 DO i=1,n
1478  WRITE (6,104) (v(i,j),j=low,upp)
1479 END DO
1480 low = low + 5
1481 IF (n < low) RETURN
1482 upp = upp + 5
1483 WRITE (6,103)
1484 GO TO 3
1485 103 FORMAT (' ')
1486 104 FORMAT (e32.14,4e25.14)
1487 END SUBROUTINE maprnt
1488 
1489 
1490 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1491 !
1492 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1493 !
1494  REAL FUNCTION random(naught)
1495 
1496  IMPLICIT NONE
1497  INTEGER, INTENT(IN OUT) :: naught
1498 
1499  REAL :: ran1,ran3(127),half
1500  INTEGER :: i,j,ran2,q,r
1501  LOGICAL :: init
1502  DATA init/.false./
1503  SAVE init,ran2,ran1,ran3
1504 
1505 IF (init) GO TO 3
1506 r = mod(naught,8190) + 1
1507 ran2 = 128
1508 DO i=1,127
1509  ran2 = ran2 - 1
1510  ran1 = -2.d0**55
1511  DO j=1,7
1512  r = mod(1756*r,8191)
1513  q = r/32
1514  ran1 = (ran1 + q)*(1.0d0/256)
1515  END DO
1516  ran3(ran2) = ran1
1517 END DO
1518 init = .true.
1519 3 IF (ran2 == 1) ran2 = 128
1520 ran2 = ran2 - 1
1521 ran1 = ran1 + ran3(ran2)
1522 half = .5d0
1523 IF (ran1 >= 0.d0) half = -half
1524 ran1 = ran1 + half
1525 ran3(ran2) = ran1
1526 random = ran1 + .5d0
1527 
1528 END FUNCTION random
1529 
1530 
1531 
1532 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1533 !
1534 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1535 !
1536  REAL FUNCTION flin (n,j,l,f,x,nf)
1537  IMPLICIT NONE
1538 
1539  INTEGER, INTENT(IN) :: n
1540  INTEGER, INTENT(IN OUT) :: j
1541  REAL, INTENT(IN) :: l
1542  REAL :: f
1543  REAL, INTENT(IN) :: x(n)
1544  INTEGER, INTENT(OUT) :: nf
1545 
1546  INTEGER :: i
1547  REAL :: t(20)
1548 
1549  EXTERNAL f
1550 
1551 !...FLIN IS THE FUNCTION OF ONE REAL VARIABLE L THAT IS MINIMIZED
1552 ! BY THE SUBROUTINE MIN...
1553 !----------------------------------------------
1554  REAL :: v(20,20),q0(20),q1(20),qa,qb,qc,qd0,qd1,qf1
1555  COMMON /q/ v,q0,q1,qa,qb,qc,qd0,qd1,qf1
1556 !----------------------------------------------
1557 IF (j == 0) GO TO 2
1558 !...THE SEARCH IS LINEAR...
1559 DO i=1,n
1560  t(i) = x(i) + l*v(i,j)
1561 END DO
1562 GO TO 4
1563 !...THE SEARCH IS ALONG A PARABOLIC SPACE CURVE...
1564 2 qa = (l*(l - qd1))/(qd0*(qd0 + qd1))
1565 qb = ((l + qd0)*(qd1 - l))/(qd0*qd1)
1566 qc = (l*(l + qd0))/(qd1*(qd0 + qd1))
1567 DO i=1,n
1568  t(i) = (qa*q0(i) + qb*x(i)) + qc*q1(i)
1569 END DO
1570 !...THE FUNCTION EVALUATION COUNTER NF IS INCREMENTED...
1571 4 nf = nf + 1
1572 flin = f(t,n)
1573 
1574 END FUNCTION flin
1575 
1576 
1577 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1578 !
1579 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1580 !
1581  SUBROUTINE sort(m,n,d,v)
1582  IMPLICIT NONE
1583 !
1584  INTEGER, INTENT(IN OUT) :: m
1585  INTEGER, INTENT(IN) :: n
1586  REAL, INTENT(IN OUT) :: d(n)
1587  REAL, INTENT(IN OUT) :: v(m,n)
1588 
1589  INTEGER :: i,j,k,nm1,ip1
1590  REAL :: s
1591 !...SORTS THE ELEMENTS OF D(N) INTO DESCENDING ORDER AND MOVES THE
1592 ! CORRESPONDING COLUMNS OF V(N,N).
1593 ! M IS THE ROW DIMENSION OF V AS DECLARED IN THE CALLING PROGRAM.
1594 IF (n == 1) RETURN
1595 nm1 = n - 1
1596 DO i = 1,nm1
1597  k=i
1598  s = d(i)
1599  ip1 = i + 1
1600  DO j = ip1,n
1601  IF (d(j) <= s) cycle
1602  k = j
1603  s = d(j)
1604  END DO
1605  IF (k <= i) cycle
1606  d(k) = d(i)
1607  d(i) = s
1608  DO j = 1,n
1609  s = v(j,i)
1610  v(j,i) = v(j,k)
1611  v(j,k) = s
1612  END DO
1613 END DO
1614 END SUBROUTINE sort
1615 
1616 
1617 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1618 !
1619 !WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW
1620 !
1621  SUBROUTINE quad(n,f,x,t,machep,h)
1622  IMPLICIT NONE
1623 
1624  INTEGER, INTENT(IN) :: n
1625  REAL :: f
1626  REAL, INTENT(IN OUT) :: x(n)
1627  REAL, INTENT(IN OUT) :: t
1628  REAL :: machep
1629  REAL, INTENT(IN OUT) :: h
1630 ! IMPLICIT REAL (A-H,O-Z)
1631  EXTERNAL f
1632 
1633 !...QUAD LOOKS FOR THE MINIMUM OF F ALONG A CURVE DEFINED BY Q0,Q1,X...
1634  INTEGER :: i
1635  REAL :: l
1636  REAL :: s,value
1637 !----------------------------------------------
1638  INTEGER :: nf,nl
1639  REAL :: fx,ldt,dmin
1640 
1641 REAL :: v(20,20),q0(20),q1(20),qa,qb,qc,qd0,qd1,qf1
1642 COMMON /global/ fx,ldt,dmin,nf,nl /q/ v,q0,q1,qa,qb,qc,qd0,qd1,qf1
1643 !----------------------------------------------
1644 s = fx
1645 fx = qf1
1646 qf1 = s
1647 qd1 = 0.d0
1648 DO i=1,n
1649  s = x(i)
1650  l = q1(i)
1651  x(i) = l
1652  q1(i) = s
1653  qd1 = qd1 + (s-l)**2
1654 END DO
1655 qd1 = sqrt(qd1)
1656 l = qd1
1657 s = 0.d0
1658 IF (qd0 <= 0.d0 .OR. qd1 <= 0.d0 .OR. nl < 3*n*n) GO TO 2
1659 value=qf1
1660 CALL min(n,0,2,s,l,value,.true.,f,x,t,machep,h)
1661 qa = (l*(l-qd1))/(qd0*(qd0+qd1))
1662 qb = ((l+qd0)*(qd1-l))/(qd0*qd1)
1663 qc = (l*(l+qd0))/(qd1*(qd0+qd1))
1664 GO TO 3
1665 2 fx = qf1
1666 qa = 0.d0
1667 qb = qa
1668 qc = 1.d0
1669 3 qd0 = qd1
1670 DO i=1,n
1671  s = q0(i)
1672  q0(i) = x(i)
1673  x(i) = (qa*s + qb*x(i)) + qc*q1(i)
1674 END DO
1675 END SUBROUTINE quad
1676 
double sig
correlated to variance of initial distribution
WWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWWW Module DFT_MODULE This module...
Definition: modules.f90:272