NLSL
fitl.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 beta 11/23/95
2 c----------------------------------------------------------------------
3 c =========================
4 c program FITL
5 c =========================
6 c
7 c Carry out nonlinear least-squares analysis of slow-motional spectra
8 c the EPRLL calculation. Uses a modification of the Levenberg-Marquardt
9 c trust-region algorithm available in MINPACK.
10 c
11 c Uses
12 c getdat
13 c setfil
14 c setdat
15 c xpack (coded below)
16 c wtdres (coded below)
17 c----------------------------------------------------------------------
18 c
19  subroutine fitl
20 c
21  use nlsdim
22  use nlsnam
23  use eprprm
24  use expdat
25  use parcom
26  use mspctr
27  use lpnam
28  use lmcom
29  use errmsg
30  use iterat
31  use stdio
32  use timer
33  use rnddbl
34 c
35  implicit none
36  integer i,ifile,iflag,iret,ix,nfev,njev
37  real cpu
38  logical snglclc,convrg
39 c
40  integer getdat,itrim
41  real dtime
42  double precision enorm,corrl,residx,wtdres
43  external getdat,lfun,enorm,dtime,itrim,corrl,residx,wtdres
44 c
45 c######################################################################
46 c
47  snglclc = maxitr.lt.1 .or. maxev.le.1
48  # .or. nprm.le.0 .or. nspc.le.0
49 c
50 c ------------------------------
51 c Open trace file (if desired)
52 c ------------------------------
53  if (itrace.ne.0) then
54  open(lutrc,file=trname(:lthfnm),status='unknown',
55  # access='sequential',form='formatted')
57  end if
58 c
59  call catchc( hltfit )
60  warn=.true.
61  call xpack(x,nprm)
62 c
63 c -----------------------------------------------------------
64 c Call lfun once if only a single calculation is specified
65 c -----------------------------------------------------------
66  if (snglclc) then
67  iflag=1
68 c cpu=dtime(tarray)
69  call lfun(ndatot,nprm,x,fvec,fjac,mxpt,iflag)
70 c cpu=dtime(tarray)
71  lmflag=1
72  info=11
73  if (itridg.ne.0) call wrtrdg(0)
74 
75 c
76  if (nspc.gt.0) then
77  do i=1,nspc
78  shft(i)=shft(i)+tmpshft(i)
79  tmpshft(i)=0.0d0
80  end do
81 c
82  fnorm=enorm(ndatot,fvec)
83 c --- weighted residual fit ---
84  if (iwflag.ne.0) then
86  rdchsq=chisqr/float(ndatot)
87 c --- unweighted residual fit ---
88  else
89  chisqr=wtdres( fvec,ndatot,nspc,ixsp,npts,rmsn )**2
90  rdchsq=chisqr/float(ndatot)
91  end if
92 c
93  write(luout,1046) cpu
94  write(luout,2036) fnorm,chisqr,rdchsq,corrl(),residx()
95 c
96  if (luout.ne.luttyo) then
97  write(luttyo,1046) cpu
98  write(luttyo,2036) fnorm,chisqr,rdchsq,corrl(),residx()
99  end if
100  else
101  write (luout,2037)
102  if (luout.ne.luttyo) write(luttyo,2037)
103  end if
104 c
105  else
106 c
107 c======================================================================
108 c Call Marquardt-Levenberg nonlinear least squares routine
109 c======================================================================
110 c
111  nprint=1
112  if (nsite.gt.1) then
114  else
115  njcol=nprm+nspc
116  end if
117 c cpu=dtime(tarray)
118 c
120  1 maxev,maxitr,diag,prscl,factor,nprint,info,nfev,njev,
121  2 ipvt,qtf,gnvec,gradf,work1,work2,work3,work4 )
122 c
123  lmflag=1
124 c cpu=dtime(tarray)
125 c
126 c *** MINPACK error return
127  convrg = info.ne.0 .and. info .le.3
128  if (.not. convrg) then
129  write (luout,2034) minerr(info)
130  write (luout,2037) nfev,cpu
131  if (luout.ne.luttyo) then
132  write (luttyo,2034) minerr(info)
133  write (luttyo,2037) nfev,cpu
134  end if
135  if (itrace.ne.0) write (lutrc,2034) minerr(info)
136 c
137 c *** Normal return
138  else
139 c
140  fnorm=enorm(ndatot,fvec)
141  if (iwflag.ne.0) then
143  rdchsq=chisqr/float(ndatot-nprm)
144  else
145  chisqr=wtdres( fvec,ndatot,nspc,ixsp,npts,rmsn )
146  rdchsq=chisqr/float(ndatot)
147  end if
148 c
149  write(luout,1000)
150  write(luout,2035) minerr(info),nprm,ndatot
151  write(luout,2037) nfev,cpu
152  write(luout,2036) fnorm,chisqr,rdchsq,corrl(),residx()
153  write(luout,1000)
154 c
155  if (luout.ne.luttyo) then
156  write(luttyo,2035) minerr(info),nprm,ndatot
157  write(luttyo,2037) nfev,cpu
158  write(luttyo,2036) fnorm,chisqr,rdchsq,corrl(),
159  # residx()
160  end if
161 c
162  if (itrace.ne.0) then
163  write(lutrc,2035) minerr(info),nprm,ndatot
164  write(lutrc,2037) nfev,cpu
165  write(lutrc,2036) fnorm,chisqr,rdchsq,corrl(),
166  # residx()
167  end if
168 c
169  call covrpt( luout )
170  if (luout.ne.luttyo) call covrpt( luttyo )
171 
172 c *** end 'if info.eq.0.or.info.gt.3 ... else'
173  end if
174 c
175 c *** end 'if single calculation...else'
176  end if
177 c
178 c ------------------------------------------------------------------
179 c Optionally output splined data and calculated spectrum/spectra
180 c for each datafile
181 c ------------------------------------------------------------------
182  if ((convrg .or. output.ne.0).and. written.eq.0) call wrspc()
183 c
184 c ----------------------------------------
185 c Report shifting and scaling information
186 c ----------------------------------------
187  call sclstt(luout)
188  if (luout.ne.luttyo) call sclstt(luttyo)
189 c
190 c ----------------------
191 c Close files and return
192 c ----------------------
193  if (itrace.ne.0) close (lutrc)
194  call uncatchc( hltfit )
195  return
196 c
197 c# format statements ##################################################
198 c
199  1000 format(/,2x,70('-'),/)
200  1046 format(//10x,'Single calculation: CPU time= ',f9.3,' sec')
201 c
202  2034 format(/2x,'MINPACK could not find a solution: ', a )
203  2035 format(/9x,'MINPACK completed: ', a/
204  # 12x,i2,' parameters;',i5,' data points')
205  2036 format(10x,'Residual norm= ',g13.6/
206  # 10x,'Chi-squared=',g13.6,2x,'Reduced Chi-sq=',g13.6/
207  # 10x,'Correlation= ',f8.5,7x,'Residual index=',f8.5/)
208  2037 format(12x,i4,' function evals (CPU time: ',f9.3,' sec)')
209  2038 format(10x,'No data files were specified')
210  end
211 c
212 c----------------------------------------------------------------------
213 c =========================
214 c subroutine XPACK
215 c =========================
216 c
217 c Pack parameters from fparm array (and spectral parameter arrays)
218 c into x vector.
219 c
220 c----------------------------------------------------------------------
221  subroutine xpack(x,n)
222 c
223  use nlsdim
224  use eprprm
225  use expdat
226  use parcom
227  use rnddbl
228 c
229  implicit none
230  integer i,ix2,n
231  double precision x(n)
232 c
233 c --- loop over all parameters in x-vector
234 c
235  do 10 i=1,n
236 c
237 c----------------------------------------------------------------------
238 c For certain parameters which cannot vary from site to site
239 c in the same spectrum, the index refers to a spectrum in the
240 c series rather than a site within a spectrum.
241 c These parameters include LB, PHASE, B0, and PSI.
242 c----------------------------------------------------------------------
243  ix2=ixst(i)
244  if (ix2.le.0) ix2=1
245  if (ixpr(i).eq.lb) then
246  x(i)=slb(ix2)
247 c
248 c --- B0: Spectrometer field
249 c
250  else if (ixpr(i).eq.ib0) then
251  x(i)=sb0(ix2)
252 c
253 c --- PHASE: Microwave phase
254 c
255  else if (ixpr(i).eq.iphase) then
256  x(i)=sphs(ix2)
257 c
258 c --- PSI: Director tilt
259 c
260  else if (ixpr(i).eq.ipsi) then
261  x(i)=spsi(ix2)
262 c
263 c----------------------------------------------------------------------
264 c All other parameters: index refers to a site in the spectrum.
265 c----------------------------------------------------------------------
266  else
267  x(i)=fparm(ixpr(i),ix2)
268  end if
269 
270  10 continue
271 c
272  return
273  end
274 
275 c----------------------------------------------------------------------
276 c =========================
277 c function WTDRES
278 c =========================
279 c
280 c Returns weighted residuals in fvec (weighted according to
281 c rms noise for each spectrum given in rmsn array)
282 c----------------------------------------------------------------------
283 
284  function wtdres( fvec, m, nspc, ixsp, npts, rmsn )
285  implicit none
286  integer m,nspc,ixsp(nspc),npts(nspc)
287  double precision fvec(m),rmsn(nspc),wtdres
288 c
289  integer isp,ixs,j,k
290 c
291  wtdres=0.0d0
292  do isp=1,nspc
293  ixs=ixsp(isp)
294  do j=1,npts(isp)
295  k=ixs+j-1
296  wtdres=wtdres+(fvec(k)/rmsn(isp))**2
297  end do
298  end do
299  return
300  end
301 
302 
double precision function wtdres(fvec, m, nspc, ixsp, npts, rmsn)
Definition: fitl.f90:285
integer, save luout
Definition: stdio.f90:32
integer, dimension(mxvar), save ixst
Definition: parcom.f90:62
integer, parameter mxpt
Definition: nlsdim.f90:39
integer, save lthfnm
Definition: nlsnam.f90:27
char info[81]
Definition: genio.c:45
subroutine lfun(m, n, x, fvec, fjac, ldfjac, iflag)
Subroutine for interfacing EPRLL spectral calculations with the MINPACK version of the Levenberg-Marq...
Definition: lfun.f90:68
integer, save ndatot
Definition: expdat.f90:45
subroutine wrtrdg(iter)
Definition: writr.f90:212
double precision, pointer, save ftol
Definition: lmcom.f90:29
double precision, pointer, save lb
Definition: eprprm.f90:58
integer, dimension(mxspc), save ixsp
Definition: expdat.f90:45
void FORTRAN() uncatchc(int *flag)
Definition: catch.c:49
integer, save itridg
Definition: parcom.f90:62
subroutine fitl
Definition: fitl.f90:20
Definition: timer.f90:10
double precision, dimension(mxspc), save tmpshft
Definition: expdat.f90:40
double precision, dimension(mxjcol), save work2
Definition: lmcom.f90:17
double precision, dimension(mxspc), save shft
Definition: expdat.f90:40
integer, save hltfit
Definition: stdio.f90:32
double precision, save rdchsq
Definition: iterat.f90:15
Definition: stdio.f90:26
double precision, dimension(mxspc), save spsi
Definition: expdat.f90:40
subroutine wrspc()
Definition: writr.f90:93
double precision, dimension(nfprm, mxsite), target, save fparm
Definition: parcom.f90:54
double precision, pointer, save gtol
Definition: lmcom.f90:29
integer, pointer, save maxev
Definition: lmcom.f90:32
double precision, dimension(mxpt), save work4
Definition: lmcom.f90:17
double precision, save fnorm
Definition: iterat.f90:15
integer, parameter ipsi
Definition: eprprm.f90:92
integer, save output
Definition: parcom.f90:62
integer, save iwflag
Definition: iterat.f90:14
double precision, dimension(mxjcol), save qtf
Definition: lmcom.f90:17
integer, dimension(mxvar), save ixpr
Definition: parcom.f90:62
integer, save written
Definition: expdat.f90:45
integer, save nprm
Definition: parcom.f90:62
integer, save nspc
Definition: expdat.f90:45
integer, save nsite
Definition: parcom.f90:62
double precision, dimension(mxspc), save sphs
Definition: expdat.f90:40
character *32, dimension(0:nmnerr-1), save minerr
Definition: errmsg.f90:88
double precision, dimension(mxpt, mxjcol), save fjac
Definition: lmcom.f90:17
subroutine covrpt(lu)
Definition: covrpt.f90:12
double precision, dimension(mxpt), save fvec
Definition: lmcom.f90:17
double factor
Definition: genio.c:44
integer, parameter ib0
Definition: eprprm.f90:92
integer, save lmflag
Definition: lmcom.f90:24
double precision, dimension(mxjcol), save gradf
Definition: lmcom.f90:17
integer, parameter lutrc
Definition: stdio.f90:29
integer, parameter iphase
Definition: eprprm.f90:92
subroutine sclstt(lu)
Definition: statc.f90:401
integer, save nprint
Definition: lmcom.f90:24
subroutine lmnls(fcn, m, n, x, fvec, fjac, ldfjac, ftol, xtol, gtol, maxfev, maxitr, diag, scale, factor, nprint, info, nfev, njev, ipvt, qtf, gnvec, gradf, wa1, wa2, wa3, wa4)
(L)evenberg-(M)arquardt (N)onlinear (L)east (S)quares This is a modification of the original lmder su...
Definition: lmnls.f90:23
void FORTRAN() catchc(int *flag)
Definition: catch.c:33
double precision, pointer, save xtol
Definition: lmcom.f90:29
Definition: lpnam.f90:18
Definition: lmcom.f90:13
character *30, save trname
Definition: nlsnam.f90:28
double precision, dimension(mxjcol), save diag
Definition: lmcom.f90:17
integer, parameter luttyo
Definition: stdio.f90:29
double precision, dimension(mxjcol), save work1
Definition: lmcom.f90:17
double precision, dimension(mxspc), save sb0
Definition: expdat.f90:40
double precision, dimension(mxspc), save rmsn
Definition: expdat.f90:40
logical, save warn
Definition: stdio.f90:33
double precision, dimension(mxjcol), save gnvec
Definition: lmcom.f90:17
subroutine xpack(x, n)
Definition: fitl.f90:222
integer, dimension(mxspc), save npts
Definition: expdat.f90:45
integer, pointer, save maxitr
Definition: lmcom.f90:32
double precision, dimension(mxspc), save slb
Definition: expdat.f90:40
integer, save itrace
Definition: stdio.f90:32
double precision, dimension(mxjcol), save work3
Definition: lmcom.f90:17
double precision, dimension(mxvar), save prscl
Definition: parcom.f90:56
double precision, save chisqr
Definition: iterat.f90:15
double precision, dimension(mxjcol), save x
Definition: lmcom.f90:17
integer, save njcol
Definition: parcom.f90:62