NLSL
covrpt.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 beta 11/23/95
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine COVRPT
5 c =========================
6 c
7 c Produces a report of the covariance matrix and parameter
8 c uncertainties on the specified logical unit.
9 c----------------------------------------------------------------------
10 c
11  subroutine covrpt( lu )
12 c
13  use nlsdim
14  use eprprm
15  use expdat
16  use parcom
17  use lmcom
18  use lpnam
19  use mspctr
20  use iterat
21  use rnddbl
22  use stdio
23 c
24  implicit none
25  integer lu
26 c
27  integer i,isi,isp,j,k,mrgn
28  double precision denom,s(5,mxsite)
29  character*132 line
30  logical potflg
31 c
32  integer itrim
33  double precision chi,ts,fs
34  external itrim,chi,ts,fs
35 c
36 c ---------------------------------------------
37 c Append labels for scale factors to tag array
38 c ---------------------------------------------
39  if (nsite.gt.1) then
40  do isi=1,nsite
41  write(tag(nprm+isi),1010) isi
42  x(nprm+isi)=sfac(isi,1)
43  end do
44  else
45  do isp=1,nser
46  write(tag(nprm+isp),1011) isp
47  x(nprm+isp)=sfac(1,isp)
48  end do
49  end if
50 c
51 c-----------------------------------------------------------------------
52 c Calculate residuals and covariance matrix
53 c COVAR calculates covariance matrix from the R matrix of QR
54 c factorization, which lmder stored in upper triangular form in fjac
55 c ipvt is permutation array returned by lmder
56 c------------------------------------------------------------------------
57 c
58  if (.not.covarok) then
59  call covar( njcol,fjac,mxpt,ipvt,xtol*xtol,work1 )
60  covarok=.true.
61 c --- weighted residual fit
62  if (iwflag.ne.0) then
63  delchi1=chi(confid,1)
64  ch2bnd=chi(confid,2)
65  chnbnd=chi(confid,njcol)
66 c --- unweighted residual fit
67  else
68  tbound=ts( (1.0d0-confid)/2.0d0, ndatot-nprm )
69  f2bnd=2.0d0*fs( 1.0d0-confid, 2, ndatot-2 )*fnorm*fnorm/
70  # float(ndatot-2)
71  fnbnd=float(njcol)*fs( 1.0d0-confid, njcol, ndatot-njcol )*
72  # fnorm*fnorm/float(ndatot-njcol)
73  end if
74 c
75  end if
76 c
77 c ------------------------------------------------------------
78 c Calculate correlation matrix from covariance and output it
79 c ------------------------------------------------------------
80  line=' '
81  mrgn=35-5*njcol
82  if (mrgn.lt.1) mrgn=1
83  write (lu,2000)
84  write(line(mrgn:),2002) (tag(i),i=1,njcol)
85  write(lu,2001) line(:10*(njcol+1)+mrgn)
86 c
87  do i=1,njcol
88  line=' '
89  do j=i,njcol
90 c
91 c----------------------------------------------------------------------
92 c Estimate of error bounds:
93 c For weighted residuals, use chi-squared statistics
94 c For unweighted residuals use t-statistics
95 c----------------------------------------------------------------------
96  if (i.eq.j) then
97  if (iwflag.ne.0) then
98  xerr(i)=sqrt(delchi1*fjac(i,i))
99  else
100  xerr(i)=tbound*fnorm*sqrt(fjac(i,i))
101  # /sqrt(float(ndatot-njcol))
102  end if
103  end if
104 
105  denom = dsqrt( fjac(i,i) * fjac(j,j) )
106  if (denom .ne. 0.0d0) corr(i,j) = fjac(i,j)/denom
107  k=10*(j-1)+mrgn
108  write(line(k:),2003) corr(i,j)
109  end do
110  write (lu,2001) line(:10*(njcol+1)+mrgn)
111  end do
112 c
113 c --------------------------------
114 c output final fit of parameters
115 c --------------------------------
116  write (lu,1000)
117  if (iwflag.ne.0) then
118  write (lu,2004) 'Chi'
119  else
120  write(lu,2004) 'T'
121  end if
122 c
123  write (lu,2007) (tag(i)(:itrim(tag(i))),x(i),xerr(i),i=1,njcol)
124 c
125 c --- weighted residual fit ---
126  if (iwflag.ne.0) then
127  write (lu,2008) confid,delchi1,ch2bnd,njcol,chnbnd
128 c
129 c --- weighted residual fit ---
130  else
131  write (lu,2009) confid,tbound,f2bnd,njcol,fnbnd
132  end if
133 c
134 c --------------------------------------------------------
135 c Calculate and output order parameters for each site if
136 c appropriate
137 c --------------------------------------------------------
138  do isi=1,nsite
139  potflg=.false.
140  do j=0,4
141  potflg=potflg.or.(abs(fparm(ic20+j,isi)).gt.rndoff)
142  end do
143 c
144  if (potflg) then
145  call ordrpr( fparm(ic20,isi),s(1,isi) )
146  write (lu,1060) isi
147  do j=0,4
148  if (abs(fparm(ic20+j,isi)).gt.rndoff)
149  # write (lu,1070) parnam(ic20+j)(2:3),s(j+1,isi)
150  end do
151 
152  end if
153 
154  end do
155 c
156  return
157 c
158 c##### Formats #########################################################
159 c
160  1000 format(/,2x,70('-'),/)
161  1010 format('SITE',i1)
162  1011 format('SPCTRM',i1)
163  2000 format(/24x,'*** Correlation Matrix ***'/)
164  2001 format(a)
165  2002 format(12(2x,a8))
166  2003 format(f8.4)
167  2004 format(/9x,'*** Final Parameters ***',/,
168  1 7x,' Parameter Value',13x,'Uncertainty (',a,
169  2 ' statistics)'/7x,46('-')/)
170  2007 format(10x,a9,' = ',g13.7,' +/- ',g13.7)
171  2008 format(/7x,'Confidence= ',f5.3,3x,'Del-Chi-Sqr= ',g11.5/
172  # 7x,'Del-Chi-Sqr for confidence regions:'/
173  # 7x,'2 parameters: ',g11.5,4x,i2,' parameters: ',g11.5)
174  2009 format(/7x,'Confidence = ',f5.3,3x,'T-bound=',g11.5/
175  # 7x,'F-bound for confidence region:'/
176  # 7x,'2 parameters: ',g11.5,4x,i2,' parameters: ',g11.5)
177  1060 format(/10x,'Order parameters for site',i2,':')
178  1070 format(15x,'<D',a2,'> = ',f10.4,' +/- ')
179 c
180  end
181 
182 
integer, parameter mxpt
Definition: nlsdim.f90:39
subroutine covar(n, r, ldr, ipvt, tol, wa)
Definition: covar.f90:2
integer, save ndatot
Definition: expdat.f90:45
double precision, dimension(mxsite, mxspc), save sfac
Definition: mspctr.f90:33
double precision, save confid
Definition: iterat.f90:15
Definition: stdio.f90:26
double precision, dimension(nfprm, mxsite), target, save fparm
Definition: parcom.f90:54
double precision, save f2bnd
Definition: iterat.f90:15
double precision, save fnorm
Definition: iterat.f90:15
integer, save iwflag
Definition: iterat.f90:14
double precision, dimension(mxjcol, mxjcol), save corr
Definition: lmcom.f90:17
double precision, save chnbnd
Definition: iterat.f90:15
integer, save nprm
Definition: parcom.f90:62
integer, save nsite
Definition: parcom.f90:62
double precision, save tbound
Definition: iterat.f90:15
double precision, dimension(mxpt, mxjcol), save fjac
Definition: lmcom.f90:17
integer, parameter mxsite
Definition: nlsdim.f90:39
subroutine covrpt(lu)
Definition: covrpt.f90:12
integer, save nser
Definition: parcom.f90:62
double precision, dimension(mxjcol), save xerr
Definition: parcom.f90:56
double precision, save ch2bnd
Definition: iterat.f90:15
subroutine ordrpr(cf, d)
Definition: ordrpr.f90:41
double precision, save fnbnd
Definition: iterat.f90:15
double precision, pointer, save xtol
Definition: lmcom.f90:29
Definition: lpnam.f90:18
Definition: lmcom.f90:13
integer, parameter ic20
Definition: eprprm.f90:92
double precision, save delchi1
Definition: iterat.f90:15
logical, save covarok
Definition: iterat.f90:18
double precision, dimension(mxjcol), save work1
Definition: lmcom.f90:17
character *6, dimension(nfprm), save parnam
Definition: lpnam.f90:27
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67
double precision, parameter rndoff
Definition: rnddbl.f90:86
double precision, dimension(mxjcol), save x
Definition: lmcom.f90:17
integer, save njcol
Definition: parcom.f90:62