NLSL
stats.f90
Go to the documentation of this file.
1 c Version 1.5.1 beta 2/3/96
2 c----------------------------------------------------------------------
3 c This file contains various functions for making statistical
4 c inferences from the residuals and curvature matrix returned
5 c by the least-squares procedure. (Subroutines marked with *)
6 c
7 c chi Returns chi-squared for given confidence, degrees freedom
8 c pchi Chi-squared probability function
9 c residx Returns residual index of fit vs. data
10 c corrl Returns correlation coefficient of fit vs. data
11 c gammq Incomplete gamma function
12 c * gcf Continued-fraction approximation to incomplete gamma fn
13 c * gser Series approximation of incomplete gamma function
14 c gammln Returns ln( gamma )
15 c * zbrac Bracketing routine for Brent root-finder
16 c zbrent Brent's method for finding root of function pchi
17 c
18 c Auxiliary functions: betai, betacf, alphat, alphaf, ts, fs
19 c Auxiliary subroutine: pausex
20 c
21 c----------------------------------------------------------------------
22 c =========================
23 c function CHI
24 c =========================
25 c Finds delta-chi for a given confidence level p and number of
26 c degrees of freedom, nu.
27 c
28 c This is done by finding the value of the chi-squared probability
29 c function, gammq( nu/2, chi^2/2 ) corresponding to the desired
30 c confidence level for a given number of degrees of freedom. The
31 c value is obtained by finding the root of function pchi (coded
32 c below), which returns the confidence level p as a function of
33 c chi. The root is first bracketed using routine zbrac, and the
34 c root found by Brent's method (function zbrent).
35 c----------------------------------------------------------------------
36  function chi( p, nu )
37  implicit none
38  double precision ts, p, chi
39  integer nu
40 c
41  double precision x1,x2
42  logical success
43 c
44  double precision halfnu, conf
45  common /chicom/ halfnu, conf
46 c
47  double precision pchi,zbrent
48  external pchi,zbrent
49 c
50  halfnu=0.5d0*nu
51  conf=p
52  x1=1.0d0
53  x2=2.0d0
54  call zbrac( pchi, x1, x2, success )
55  if (success) chi=zbrent( pchi, x1, x2, 1.0d-6 )
56  return
57  end
58 
59 
60 c----------------------------------------------------------------------
61 c =========================
62 c function PCHI
63 c =========================
64 c Returns the chi-squared probability function given chi-squared
65 c and the number degrees of freedom, nu.
66 c----------------------------------------------------------------------
67  function pchi(chi)
68  implicit none
69  double precision pchi,chi
70 c
71  double precision halfnu, conf
72  common /chicom/ halfnu, conf
73 c
74  double precision gammq
75  external gammq
76 c
77  pchi = gammq( halfnu, 0.5d0*chi ) - 1.0d0 + conf
78  return
79  end
80 c----------------------------------------------------------------------
81 c =========================
82 c function RESIDX
83 c =========================
84 c Returns the residual index between experimental data
85 c and calculated fit function in common blocks /expdat/ and /lmcom/
86 c----------------------------------------------------------------------
87  function residx()
88 c
89  use nlsdim
90  use expdat
91  use lmcom
92 c
93  implicit none
94  integer i,j,k
95  double precision osum,residx,rsum,yc,yo
96 c
97  osum=0.0d0
98  rsum=0.0d0
99  do i=1,nspc
100  do j=1,npts(i)
101  k=ixsp(i)+j-1
102  yo=data(k)
103  yc=data(k)-rmsn(i)*fvec(k)
104  rsum=rsum+abs(yo-yc)
105  osum=osum+abs(yo)
106  end do
107  end do
108  residx=rsum/osum
109  return
110  end
111 
112 c----------------------------------------------------------------------
113 c =========================
114 c function CORRL
115 c =========================
116 c Returns the correlation function between experimental data
117 c and calculated fit function in common blocks /expdat/ and /lmcom/
118 c----------------------------------------------------------------------
119  function corrl()
120 c
121  use nlsdim
122  use expdat
123  use lmcom
124 c
125  implicit none
126  integer i,j,k
127  double precision cbar,cdsum,c2sum,dbar,dc,dd,d2sum,corrl
128 c
129  dbar=0.0d0
130  cbar=0.0d0
131  do i=1,nspc
132  do j=1,npts(i)
133  k=ixsp(i)+j-1
134  dbar=dbar+data(k)
135  cbar=cbar+(data(k)-rmsn(i)*fvec(k))
136  end do
137  end do
138  dbar=dbar/float(ndatot)
139  cbar=cbar/float(ndatot)
140 c
141  cdsum=0.0d0
142  d2sum=0.0d0
143  c2sum=0.0d0
144  do i=1,nspc
145  do j=1,npts(i)
146  k=ixsp(i)+j-1
147  dd=(data(k)-dbar)
148  dc=(data(k)-rmsn(i)*fvec(k)-cbar)
149  cdsum=cdsum+dd*dc
150  d2sum=d2sum+dd*dd
151  c2sum=c2sum+dc*dc
152  end do
153  end do
154  corrl=cdsum/sqrt(d2sum*c2sum)
155 c
156  return
157  end
158 c----------------------------------------------------------------------
159 c =========================
160 c function BETAI
161 c =========================
162 c Incomplete beta function, used in calculation of t-distribution
163 c and F-distribution used in statistical analysis of parameters
164 c in NLS fits.
165 c
166 c Function returns the function I_x(A,B)
167 c
168 c From Numerical Recipes by W.H. Press et al.
169 c----------------------------------------------------------------------
170 c
171  function betai(a,b,x)
172  implicit none
173  double precision betai,a,b,x,bt
174 c
175  double precision ONE,TWO,ZERO
176  parameter(one=1.0d0,two=2.0d0,zero=0.0d0)
177 c
178  double precision betacf,gammln
179  external betacf,gammln
180 c
181  if (x.lt.zero .or. x.gt.one) then
182  call pausex('bad argument X in BETAI')
183  end if
184  if (x.eq.zero .or. x.eq.one) then
185  bt=zero
186  else
187 c
188 c ----------------------------------------------------
189 c Calculate factors in front of the continued fraction
190 c ----------------------------------------------------
191  bt=dexp( gammln(a+b)-gammln(a)-gammln(b)
192  # +a*dlog(x)+b*dlog(one-x) )
193  end if
194 c
195 c --------------------------------
196 c Use continued fraction directly
197 c --------------------------------
198  if (x.lt.(a+one)/(a+b+two)) then
199  betai=bt*betacf(a,b,x)/a
200  return
201 c
202 c ----------------------------------------------------------------
203 c Use continued fraction after making the symmetry transformation
204 c ----------------------------------------------------------------
205  else
206  betai=one-bt*betacf(b,a,one-x)/b
207  return
208  end if
209  end
210 
211 
212  function betacf(a,b,x)
213  implicit none
214  double precision betacf,a,b,x
215 c
216  integer m
217  double precision am,aold,ap,app,az,bp,bm,bpp,bz,d,em,qab,qap,
218  # qam,tem
219 c
220  integer ITMAX
221  double precision ONE,TWO,ZERO,EPS
222  parameter(itmax=100,eps=3.0d-7,one=1.0d0,two=2.0d0,zero=0.0d0)
223 c
224  am=one
225  bm=one
226  az=one
227  qab=a+b
228  qap=a+one
229  qam=a-one
230  bz=one-qab*x/qap
231  do m=1,itmax
232  em=m
233  tem=em+em
234 c
235 c ----------------------------
236 c Even step of the recurrence
237 c ----------------------------
238  d=em*(b-m)*x/((qam+tem)*(a+tem))
239  ap=az+d*am
240  bp=bz+d*bm
241 c
242 c ----------------------------
243 c Odd step of the recurrence
244 c ----------------------------
245  d=-(a+em)*(qab+em)*x/((a+tem)*(qap+tem))
246  app=ap+d*az
247  bpp=bp+d*bz
248 c
249 c ----------------------------------------------------
250 c Save old answer and renormalize to prevent overflow
251 c ----------------------------------------------------
252  aold=az
253  am=ap/bpp
254  bm=bp/bpp
255  az=app/bpp
256  bz=one
257 c
258 c ----------------------
259 c Check for convergence
260 c ----------------------
261  if (abs(az-aold).lt.eps*abs(az)) goto 1
262  end do
263 c
264  write(*,1000) abs((az-aold)/az)
265  1000 format('BETACF: ITMAX exceeded, final error was ',g11.3)
266 c
267  1 betacf=az
268  return
269  end
270 
271  function alphat(t)
272  implicit none
273  double precision alphat,t,xnu
274 c
275  double precision a,b,alpha
276  common /bcom/ a,b,alpha
277 c
278  double precision betai
279  external betai
280 c
281  alphat = alpha - 0.5*betai( a, b, a/(a+0.5d0*t*t) )
282  return
283  end
284 
285 
286  function alphaf(f)
287  implicit none
288  double precision alphaf,f
289 c
290  double precision a,b,alpha
291  common /bcom/ a,b,alpha
292 c
293  double precision betai
294  external betai
295 c
296  alphaf = betai( b, a, b/(b+a*f) ) - alpha
297  return
298  end
299 
300  function ts( al, nu )
301  implicit none
302  double precision ts, al
303  integer nu
304 c
305  double precision x1,x2
306  logical success
307 c
308  double precision a, b, alpha
309  common /bcom/ a, b, alpha
310 c
311  double precision alphat,zbrent
312  external alphat,zbrent
313 c
314  alpha=al
315  a=0.5d0*nu
316  b=0.5d0
317  x1=0.0d0
318  x2=1.0d0
319  call zbrac( alphat, x1, x2, success )
320  if (success) ts=zbrent( alphat, x1, x2, 1.0d-6 )
321  return
322  end
323 
324 
325  function fs( al, nu1, nu2 )
326  implicit none
327  double precision fs, al
328  integer nu1,nu2
329 c
330  double precision x1,x2
331  logical success
332 c
333  double precision a, b, alpha
334  common /bcom/ a, b, alpha
335 c
336  double precision alphaf,zbrent
337  external alphaf,zbrent
338 c
339  alpha=al
340  a=0.5d0*nu1
341  b=0.5d0*nu2
342  x1=0.0d0
343  x2=1.0d0
344  call zbrac( alphaf, x1, x2, success )
345  if (success) fs=zbrent( alphaf, x1, x2, 1.0d-6 )
346  return
347  end
348 
349 c----------------------------------------------------------------------
350 c =========================
351 c function GAMMQ
352 c =========================
353 c
354 c Incomplete gamma function. From Numerical Recipes by Press, et al.
355 c
356 c----------------------------------------------------------------------
357 c
358  FUNCTION gammq(a,x)
359  IMPLICIT none
360  DOUBLE PRECISION a,gammq,x
361 CU USES gcf,gser
362  DOUBLE PRECISION gammcf,gamser,gln
363  if(x.lt.0..or.a.le.0.) call pausex('bad arguments in gammq')
364  if(x.lt.a+1.)then
365  call gser(gamser,a,x,gln)
366  gammq=1.-gamser
367  else
368  call gcf(gammcf,a,x,gln)
369  gammq=gammcf
370  endif
371  return
372  END
373 c
374 c----------------------------------------------------------------------
375 c =========================
376 c function GCF
377 c =========================
378 c
379 c Continued-fraction calculation of incomplete gamma function.
380 c From Numerical Recipes by Press, et al.
381 c
382 c----------------------------------------------------------------------
383  SUBROUTINE gcf(gammcf,a,x,gln)
384  IMPLICIT none
385  INTEGER ITMAX
386  DOUBLE PRECISION a,gammcf,gln,x,EPS,FPMIN
387  parameter(itmax=100,eps=3.e-7,fpmin=1.e-30)
388 C USES gammln
389  INTEGER i
390  DOUBLE PRECISION an,b,c,d,del,h,gammln
391  gln=gammln(a)
392  b=x+1.-a
393  c=1./fpmin
394  d=1./b
395  h=d
396  do 11 i=1,itmax
397  an=-i*(i-a)
398  b=b+2.
399  d=an*d+b
400  if(abs(d).lt.fpmin)d=fpmin
401  c=b+an/c
402  if(abs(c).lt.fpmin)c=fpmin
403  d=1./d
404  del=d*c
405  h=h*del
406  if(abs(del-1.).lt.eps)goto 1
407 11 continue
408  call pausex('a too large, ITMAX too small in gcf')
409 1 gammcf=exp(-x+a*log(x)-gln)*h
410  return
411  END
412 
413 c----------------------------------------------------------------------
414 c =========================
415 c function GSER
416 c =========================
417 c
418 c Series approximation of incomplete gamma function.
419 c From Numerical Recipes by Press, et al.
420 c
421 c----------------------------------------------------------------------
422  SUBROUTINE gser(gamser,a,x,gln)
423  IMPLICIT none
424  INTEGER ITMAX
425  DOUBLE PRECISION a,gamser,gln,x,EPS
426  parameter(itmax=100,eps=3.e-7)
427 C USES gammln
428  INTEGER n
429  DOUBLE PRECISION ap,del,sum,gammln
430  gln=gammln(a)
431  if(x.le.0.)then
432  if(x.lt.0.) call pausex('x < 0 in gser')
433  gamser=0.
434  return
435  endif
436  ap=a
437  sum=1./a
438  del=sum
439  do 11 n=1,itmax
440  ap=ap+1.
441  del=del*x/ap
442  sum=sum+del
443  if(abs(del).lt.abs(sum)*eps)goto 1
444 11 continue
445  call pausex('a too large, ITMAX too small in gser')
446 1 gamser=sum*exp(-x+a*log(x)-gln)
447  return
448  END
449 c
450 c----------------------------------------------------------------------
451 c =========================
452 c function GAMMLN
453 c =========================
454 c
455 c Returns the value ln(gamma(xx) for xx>1, with full accuracy for xx>1
456 c From Numerical Recipes by Press, et al.
457 c
458 c----------------------------------------------------------------------
459  FUNCTION gammln(xx)
460  IMPLICIT none
461  DOUBLE PRECISION gammln,xx
462  INTEGER j
463  DOUBLE PRECISION ser,stp,tmp,x,y,cof(6)
464  SAVE cof,stp
465  DATA cof,stp/76.18009172947146d0,-86.50532032941677d0,
466  *24.01409824083091d0,-1.231739572450155d0,.1208650973866179d-2,
467  *-.5395239384953d-5,2.5066282746310005d0/
468  x=xx
469  y=x
470  tmp=x+5.5d0
471  tmp=(x+0.5d0)*log(tmp)-tmp
472  ser=1.000000000190015d0
473  do 11 j=1,6
474  y=y+1.d0
475  ser=ser+cof(j)/y
476 11 continue
477  gammln=tmp+log(stp*ser/x)
478  return
479  END
480 
481 c
482 c######################################################################
483 c
484  SUBROUTINE zbrac(func,x1,x2,succes)
485  IMPLICIT none
486  INTEGER NTRY
487  DOUBLE PRECISION x1,x2,func,FACTOR
488  EXTERNAL func
489  parameter(factor=1.6,ntry=50)
490  INTEGER j
491  DOUBLE PRECISION f1,f2
492  LOGICAL succes
493  if (x1.eq.x2)
494  * call pausex('you have to guess an initial range in zbrac')
495  f1=func(x1)
496  f2=func(x2)
497  succes=.true.
498  do 11 j=1,ntry
499  if(f1*f2.lt.0.)return
500  if(abs(f1).lt.abs(f2))then
501  x1=x1+factor*(x1-x2)
502  f1=func(x1)
503  else
504  x2=x2+factor*(x2-x1)
505  f2=func(x2)
506  endif
507 11 continue
508  succes=.false.
509  return
510  END
511  FUNCTION zbrent(func,x1,x2,tol)
512  IMPLICIT none
513  INTEGER ITMAX
514  DOUBLE PRECISION zbrent,tol,x1,x2,func,EPS
515  EXTERNAL func
516  parameter(itmax=100,eps=3.e-8)
517  INTEGER iter
518  DOUBLE PRECISION a,b,c,d,e,fa,fb,fc,p,q,r,s,tol1,xm
519  a=x1
520  b=x2
521  fa=func(a)
522  fb=func(b)
523  if((fa.gt.0..and.fb.gt.0.).or.(fa.lt.0..and.fb.lt.0.))
524  * call pausex('root must be bracketed for zbrent')
525  c=b
526  fc=fb
527  do 11 iter=1,itmax
528  if((fb.gt.0..and.fc.gt.0.).or.(fb.lt.0..and.fc.lt.0.))then
529  c=a
530  fc=fa
531  d=b-a
532  e=d
533  endif
534  if(abs(fc).lt.abs(fb)) then
535  a=b
536  b=c
537  c=a
538  fa=fb
539  fb=fc
540  fc=fa
541  endif
542  tol1=2.*eps*abs(b)+0.5*tol
543  xm=.5*(c-b)
544  if(abs(xm).le.tol1 .or. fb.eq.0.)then
545  zbrent=b
546  return
547  endif
548  if(abs(e).ge.tol1 .and. abs(fa).gt.abs(fb)) then
549  s=fb/fa
550  if(a.eq.c) then
551  p=2.*xm*s
552  q=1.-s
553  else
554  q=fa/fc
555  r=fb/fc
556  p=s*(2.*xm*q*(q-r)-(b-a)*(r-1.))
557  q=(q-1.)*(r-1.)*(s-1.)
558  endif
559  if(p.gt.0.) q=-q
560  p=abs(p)
561  if(2.*p .lt. min(3.*xm*q-abs(tol1*q),abs(e*q))) then
562  e=d
563  d=p/q
564  else
565  d=xm
566  e=d
567  endif
568  else
569  d=xm
570  e=d
571  endif
572  a=b
573  fa=fb
574  if(abs(d) .gt. tol1) then
575  b=b+d
576  else
577  b=b+sign(tol1,xm)
578  endif
579  fb=func(b)
580 11 continue
581  call pausex('zbrent exceeding maximum iterations')
582  zbrent=b
583  return
584  END
585 
586  subroutine pausex(mesg)
587  character*(*) mesg
588  print *, mesg
589  print *, '[execution paused, press enter to continue]'
590  read (*,*)
591  return
592  end
double precision function corrl()
Definition: stats.f90:120
double precision function betai(a, b, x)
Definition: stats.f90:172
double precision function zbrent(func, x1, x2, tol)
Definition: stats.f90:512
integer, save ndatot
Definition: expdat.f90:45
subroutine gser(gamser, a, x, gln)
Definition: stats.f90:423
double precision function gammln(xx)
Definition: stats.f90:460
integer, dimension(mxspc), save ixsp
Definition: expdat.f90:45
double precision function chi(p, nu)
Definition: stats.f90:37
double precision function alphat(t)
Definition: stats.f90:272
double precision function betacf(a, b, x)
Definition: stats.f90:213
double precision function alphaf(f)
Definition: stats.f90:287
integer, save nspc
Definition: expdat.f90:45
double precision function residx()
Definition: stats.f90:88
double precision function ts(al, nu)
Definition: stats.f90:301
double precision, dimension(mxpt), save fvec
Definition: lmcom.f90:17
subroutine zbrac(func, x1, x2, succes)
Definition: stats.f90:485
Definition: lmcom.f90:13
double precision function fs(al, nu1, nu2)
Definition: stats.f90:326
double precision function gammq(a, x)
Definition: stats.f90:359
double precision, dimension(mxspc), save rmsn
Definition: expdat.f90:40
integer, dimension(mxspc), save npts
Definition: expdat.f90:45
subroutine pausex(mesg)
Definition: stats.f90:587
subroutine gcf(gammcf, a, x, gln)
Definition: stats.f90:384
double precision function pchi(chi)
Definition: stats.f90:68