NLSL
statc.f90
Go to the documentation of this file.
1 c Version 1.5.1 beta 2/3/96
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine STATC
5 c =========================
6 c
7 c Prints out the current status of the fitting procedure to the given
8 c logical unit number. Includes information on fitting variables and
9 c user-specified convergence tolerances, etc.
10 c
11 c Includes:
12 c nlsdim.inc
13 c eprprm.inc
14 c expdat.inc
15 c parcom.inc
16 c lmcom.inc
17 c
18 c----------------------------------------------------------------------
19  subroutine statc( line )
20 c
21  use nlsdim
22  use eprprm
23  use expdat
24  use parcom
25  use lmcom
26  use stdio
27 c
28  implicit none
29  character*80 line
30 c
31  integer NKEYWD
32  parameter(nkeywd=7)
33 c
34  integer i,lth
35  logical fitflg,datflg,basflg,tdgflg,sclflg
36  character*8 keywrd(nkeywd)
37  character*30 token
38 c
39  integer itrim
40  external itrim
41 c
42  fitflg=.false.
43  datflg=.false.
44  basflg=.false.
45  tdgflg=.false.
46  sclflg=.false.
47 c
48  data keywrd /'FIT','DATA','BASIS','TRIDIAG','SHIFT','SCALE','ALL'/
49 c
50 c----------------------------------------------------------------------
51 c Look for a keyword
52 c----------------------------------------------------------------------
53  5 call gettkn(line,token,lth)
54  if (lth.ne.0) then
55  lth=min(lth,8)
56  call touppr(token,lth)
57  do i=1,nkeywd
58  if (token(:lth).eq.keywrd(i)(:lth)) go to 7
59  end do
60  write (luout,1000) token(:itrim(token))
61  if (luout.ne.luttyo) write (luttyo,1000) token(:itrim(token))
62  go to 5
63 c
64  7 if (i.eq.1) then
65  fitflg=.true.
66  else if (i.eq.2) then
67  datflg=.true.
68  else if (i.eq.3) then
69  basflg=.true.
70  else if (i.eq.4) then
71  tdgflg=.true.
72  else if (i.eq.5 .or. i.eq.6) then
73  sclflg=.true.
74  else if (i.eq.7) then
75  fitflg=.true.
76  datflg=.true.
77  basflg=.true.
78  tdgflg=.true.
79  sclflg=.true.
80  end if
81  go to 5
82  end if
83 c
84  fitflg = .not.(datflg.or.basflg.or.tdgflg.or.sclflg)
85  if (fitflg) then
86  call fitstt( luout )
87  call sclstt( luout )
88  if (luout.ne.luttyo) then
89  call fitstt( luttyo )
90  call sclstt( luttyo )
91  end if
92  end if
93 c
94  if (sclflg.and..not.fitflg) then
95  call sclstt( luout )
96  if (luout.ne.luttyo) call sclstt( luttyo )
97  end if
98 c
99  if (datflg) then
100  call datstt( luout )
101  if (luout.ne.luttyo) call datstt( luttyo )
102  end if
103 c
104  if (basflg) then
105  call basstt( luout )
106  if (luout.ne.luttyo) call basstt( luttyo )
107  end if
108 c
109  if (tdgflg) then
110  call tdgstt( luout )
111  if (luout.ne.luttyo) call tdgstt( luttyo )
112  end if
113 c
114  return
115 c
116 c######################################################################
117 c
118  1000 format(' *** unrecognized STATUS keyword: ''',a,''' ***')
119  end
120 
121 c----------------------------------------------------------------------
122 c =========================
123 c subroutine FITSTT
124 c =========================
125 c
126 c Output status report for NLS fitting parameters on logical unit lu
127 c----------------------------------------------------------------------
128 
129  subroutine fitstt( lu )
130 c
131  use nlsdim
132  use eprprm
133  use expdat
134  use parcom
135  use lmcom
136  use errmsg
137  use iterat
138  use stdio
139 c
140  implicit none
141  integer lu
142 c
143  integer i
144 c
145  double precision wtdres
146  external wtdres
147 c
148 c----------------------------------------------------------------------
149 c parameters for fit command
150 c----------------------------------------------------------------------
151  write(lu,1000)
152  write(lu,1001) xtol,ftol,gtol
153  write(lu,1002) maxitr,maxev,factor
154 c
155 c----------------------------------------------------------------------
156 c Options for fit command
157 c----------------------------------------------------------------------
158  if (nshift.lt.0) write (lu,1005)
159  if (nshift.eq.0) write (lu,1006) 100.0d0*srange
160  if (nshift.gt.0) write (lu,1007) nshift,nshift,100.0d0*srange
161  if (noneg.eq.0) write (lu,1008)
162  if (noneg.ne.0) write (lu,1009)
163  if (output.eq.1) write (lu,1014)
164  if (itrace.eq.1) write (lu,1015)
165 c
166  if (iwflag.ne.0) then
167  write (lu,1013) 'WEIGHTED'
168  else
169  write (lu,1013) 'UNWEIGHTED'
170  end if
171 c
172 c----------------------------------------------------------------------
173 c variable information
174 c----------------------------------------------------------------------
175  if (nprm.ge.1) then
176  write(lu,1010) nprm
177  write(lu,1011)
178 c
179  call xpack(x,nprm)
180  do i=1,nprm
181  write(lu,1012) tag(i),x(i),prscl(i),xfdstp(i)
182  end do
183  else
184  write(lu,1020)
185  end if
186 c
187  if (info.lt.0) info=0
188  if (info.gt.10) info=10
189  if (lmflag.ne.0) then
190  if (iwflag.ne.0) then
192  rdchsq=chisqr/dfloat(ndatot-nprm)
193  else
194  chisqr=wtdres( fvec,ndatot,nspc,ixsp,npts,rmsn )
195  rdchsq=chisqr/dfloat(ndatot-nprm)
196  end if
197 
198  write(lu,1030) fnorm,chisqr,rdchsq,minerr(info)
199  else
200  write(lu,1031)
201  end if
202  return
203 c
204 c# formats ##############################################################
205 c
206  1000 format(/,2x,20('#'),' FIT status ',20('#'))
207  1001 format(5x,'xtol = ',1p,e8.1,3x,'ftol = ',1p,e8.1,3x,'gtol = ',
208  # 1p,e8.1)
209  1002 format(5x,'maxitr=',i3,' maxfun=',i4,' bound=',f8.1//
210  #2x,'Options in effect:')
211  1005 format(5x,'NOSHIFT: No spectral shifting' )
212  1006 format(5x,'SHIFT: Spectral shifting (srange= +/-',f5.1,'%)')
213  1007 format(5x,'SHIFT',i4,': Spectral shifting first',i4,
214  #' iterations (srange= +/-',f5.1,'%)')
215  1008 format(5x,'NEG: negative scaling coefficients allowed')
216  1009 format(5x,'NONEG: negative scaling coefficients not allowed')
217  1010 format(/,2x,'VARY parameters : ',i1,' variables')
218  1011 format(2x,56('-'),/,6x,'Parameter',8x,'Value',10x,'Scale',5x,
219  # 'FD step',/,2x,56('-'))
220  1012 format(7x,a9,5x,g14.7,3x,f6.1,5x,1p,e8.1)
221  1013 format(5x,a,' residuals used for minimization')
222  1014 format(5x,'WRITE: always write .spc file')
223  1015 format(5x,'TRACE: write .trc file')
224  1020 format(/,2x,'No parameters are being varied.')
225  1030 format(/,2x,'Rsd. norm=',g14.7,2x,'Chi-sqr= ',g14.7,
226  # 2x,'Red. chi-sqr=',g14.7,
227  # /,2x,'Status: ',a/)
228  1031 format(/,2x,'No fit has been completed')
229  end
230 
231 c----------------------------------------------------------------------
232 c =========================
233 c subroutine DATSTT
234 c =========================
235 c
236 c Output status report for NLS datafile space on logical unit lu
237 c----------------------------------------------------------------------
238 
239  subroutine datstt( lu )
240 c
241  use nlsdim
242  use expdat
243  use parcom
244 c
245  implicit none
246  integer lu
247 c
248  character*2 shftopt,normopt
249  integer i
250 c
251  character*6 formopt(2)
252  data formopt/'ASCII','BINARY'/
253 c
254  write (lu,1000) nspc,ndatot,mxpt-ndatot
255  if (nspc.lt.nser) write (lu,1004) nser
256  if (nspc.gt.0) then
257  write (lu,1001)
258  do i=1,nspc
259  shftopt='N'
260  normopt='N'
261  if (ishft(i).ne.0) shftopt='Y'
262  if (nrmlz(i).ne.0) normopt='Y'
263  write (lu,1002) i,dataid(i),npts(i),ibase(i),idrv(i),
264  # shftopt,normopt,sb0(i),shft(i),rmsn(i)
265  end do
266  end if
267 c
268  shftopt=' '
269  normopt=' '
270  if (shftflg.eq.0) shftopt='NO'
271  if (normflg.eq.0) normopt='NO'
272  write (lu,1003) bcmode,nspline,shftopt,normopt,formopt(inform+1)
273  return
274 c
275 c######################################################################
276 c
277  1000 format(/2x,25('#'),' Datafile status ',25('#')//
278  # i2,' data files and',i4,' data points in buffer'/
279  # 'Space available for',i5,' points ')
280  1001 format(/2x,'File Filename # Pts Bsln Drv',
281  #' Shft Nrm B0 Offset RMS noise'/70('-'))
282  1002 format(2x,i3,4x,a,t20,2(2x,i4),3x,i1,2(4x,a1),f9.1,f9.2,g9.3)
283  1003 format(/'Options in effect:'/5x,'bcmode=',i3,', nspline=',i4,
284  #', ',a2,'SHIFT, ',a2,'NORM, ',a/)
285  1004 format(i2,' data files are required for the series fit')
286  end
287 
288 
289 
290 
291 
292 
293 c----------------------------------------------------------------------
294 c =========================
295 c subroutine BASSTT
296 c =========================
297 c
298 c Output status report for NLS datafile space on logical unit lu
299 c----------------------------------------------------------------------
300 
301  subroutine basstt( lu )
302 c
303  use nlsdim
304  use expdat
305  use parcom
306  use basis
307  use mtsdef
308 c
309  implicit none
310  integer lu
311 c
312  integer i,j
313 c
314  write (lu,1000) nbas,mxdim-nextbs+1
315  if (nbas.gt.0) then
316  write (lu,1001)
317  do i=1,nbas
318  write (lu,1002) i,basisid(i),ltbas(i),(mts(j,i),j=1,ntrc)
319  #
320  end do
321 c
322  write (lu,1003) (i,i=1,nsite)
323  do j=1,nser
324  write(lu,1004) dataid(j),(basno(i,j),i=1,nsite)
325  end do
326 c
327  end if
328 c
329  write(lu,'(/)')
330  return
331 c
332 c######################################################################
333 c
334  1000 format(/2x,i2,' BASIS SETS in buffer; Space for ',i6,' elements'/)
335  1001 format(2x,'Set Identification ndim Trunc. indices'/
336  #2x,60('-'))
337  1002 format(2x,i3,3x,a,t29,i5,4x,6(i3,','),i2)
338  1003 format(/2x,'Basis set usage:'/t24,5('Site',i2,2x:))
339  1004 format(2x,a,t24,5(i3,5x))
340  end
341 
342 
343 
344 c----------------------------------------------------------------------
345 c =========================
346 c subroutine TDGSTT
347 c =========================
348 c
349 c
350 c Output status report for NLS datafile space on logical unit lu
351 c----------------------------------------------------------------------
352 
353  subroutine tdgstt( lu )
354 c
355  use nlsdim
356  use expdat
357  use parcom
358  use tridag
359 c
360  implicit none
361  integer lu
362 c
363  integer i,isi,isp
364 c
365  write (lu,1000) ntd,mxtdg-nexttd+1
366  if (ntd.gt.0) then
367  write (lu,1001)
368  do i=1,ntd
369  isi=tdsite(i)
370  isp=tdspec(i)
371  if (modtd(isi,isp).ne.0) then
372  write (lu,1002) i,isi,dataid(isp)
373  else
374  write (lu,1003) i,isi,dataid(isp),ltd(isi,isp)
375  end if
376  end do
377  end if
378 c
379  write(lu,'(/)')
380  return
381 c
382 c######################################################################
383 c
384  1000 format(/2x,i2,' TRIDIAGONAL MATRICES in buffer; Space for',i6,
385  # ' elements'/)
386  1001 format(2x,'Matrix Site Spectrum Lth'/2x,40('-'))
387  1002 format(2x,i3,4x,i3,4x,a,t34,'(modified)')
388  1003 format(2x,i3,4x,i3,4x,a,t34,i5)
389  end
390 
391 
392 
393 c----------------------------------------------------------------------
394 c =========================
395 c subroutine SCLSTT
396 c =========================
397 c
398 c Output shifting and scaling information to specified logical unit
399 c----------------------------------------------------------------------
400  subroutine sclstt(lu)
401 c
402  use nlsdim
403  use parcom
404  use expdat
405  use mspctr
406  use nlsnam
407 c
408  implicit none
409  integer lu
410 c
411  integer i,j,k,nsp
412  double precision stot(mxspc)
413  character*1 shflg(mxspc),scflg(mxsite)
414 c
415 c ----------------------------------------------------------------
416 c Calculate total of scaling factor in order to report percentages
417 c Set flags showing which spectra are enabled for shifting
418 c ----------------------------------------------------------------
419  do j=1,nspc
420  shflg(j)='N'
421  if (ishft(j).ne.0) shflg(j)='Y'
422 c
423  stot(j)=0.0d0
424  do i=1,nsite
425  stot(j)=stot(j)+sfac(i,j)
426  end do
427  stot(j)=stot(j)/1.0d2
428  end do
429 c
430 c ---------------------------------------------
431 c Set flags showing which sites are autoscaled
432 c ---------------------------------------------
433  do i=1,nsite
434  scflg(i)=' '
435  if (iscal(i).ne.0) scflg(i)='*'
436  end do
437 c
438  if (nsite.gt.1.and.nspc.gt.1) then
439  nsp=1
440  else
441  nsp=nspc
442  end if
443 c
444 c -----------------------------------------------------
445 c Print out scaling information for each spectrum/site
446 c -----------------------------------------------------
447  if (nspc.gt.0) then
448  if (nsite.gt.1) then
449  write(lu,1047) (i,i=1,nsite)
450  do j=1,nsp
451  write(lu,1048) (sfac(i,j),sfac(i,j)/stot(j), i=1,nsite)
452  end do
453  end if
454 c
455 c --------------------------------------------------
456 c Print out shifting information for each spectrum
457 c --------------------------------------------------
458  k=lthdnm-4
459  write (lu,1049) (dataid(j)(:k),shft(j),shflg(j),sb0(j),
460  # sb0(j)+shft(j),j=1,nspc)
461  endif
462  return
463 c
464 c ##### formats ######################################################
465 c
466  1047 format(/2x,'Scales:',t10,5(11x,'site',i2,2x))
467  1048 format(12x,t16,5(g12.4,'(',f5.1,'%) '))
468  1049 format(/2x,'Shifts: file',t33,'shift',2x,'auto',2x,4x,'B0',
469  # 6x,'B0(eff)'/(12x,a,t30,f8.2,3x,a1,3x,f9.2,1x,f9.2))
470  end
471 
subroutine fitstt(lu)
Definition: statc.f90:130
integer, save luout
Definition: stdio.f90:32
subroutine basstt(lu)
Definition: statc.f90:302
integer, dimension(mxspc), save nrmlz
Definition: expdat.f90:45
integer, parameter mxpt
Definition: nlsdim.f90:39
char info[81]
Definition: genio.c:45
integer, save noneg
Definition: parcom.f90:62
double precision, save srange
Definition: parcom.f90:56
integer, save ndatot
Definition: expdat.f90:45
integer, dimension(mxtdm), save tdspec
Definition: tridag.f90:28
double precision, pointer, save ftol
Definition: lmcom.f90:29
integer, dimension(mxspc), save ixsp
Definition: expdat.f90:45
integer, save nexttd
Definition: tridag.f90:28
character *30, dimension(mxtdm), save basisid
Definition: basis.f90:27
integer, save inform
Definition: expdat.f90:45
integer, dimension(mxspc), save ishft
Definition: expdat.f90:45
integer, parameter ntrc
Definition: mtsdef.f90:29
integer, dimension(mxsite, mxspc), save basno
Definition: basis.f90:23
integer, dimension(mxspc), save ibase
Definition: expdat.f90:45
integer, dimension(mxsite, mxspc), save modtd
Definition: tridag.f90:28
double precision, dimension(mxsite, mxspc), save sfac
Definition: mspctr.f90:33
double precision, dimension(mxspc), save shft
Definition: expdat.f90:40
double precision, save rdchsq
Definition: iterat.f90:15
Definition: stdio.f90:26
character *30, dimension(mxspc), save dataid
Definition: expdat.f90:51
integer, save lthdnm
Definition: nlsnam.f90:27
double precision, pointer, save gtol
Definition: lmcom.f90:29
subroutine touppr(string, lth)
Definition: strutl2.f90:22
integer, save normflg
Definition: expdat.f90:45
integer, pointer, save maxev
Definition: lmcom.f90:32
Definition: basis.f90:19
double precision, save fnorm
Definition: iterat.f90:15
integer, save output
Definition: parcom.f90:62
integer, save iwflag
Definition: iterat.f90:14
subroutine gettkn(line, token, lth)
Written for free-form input of parameters for slow-motional calculations. Returns a token consisting ...
Definition: strutl1.f90:75
integer, dimension(9, mxmts), save mts
Definition: basis.f90:23
integer, dimension(mxtdm), save tdsite
Definition: tridag.f90:28
double precision, dimension(mxvar), save xfdstp
Definition: parcom.f90:56
integer, parameter mxspc
Definition: nlsdim.f90:39
integer, save bcmode
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
integer, dimension(mxsite, mxspc), save ltd
Definition: tridag.f90:28
integer, save nextbs
Definition: basis.f90:23
integer, save nshift
Definition: parcom.f90:62
character *32, dimension(0:nmnerr-1), save minerr
Definition: errmsg.f90:88
integer, dimension(mxsite), save iscal
Definition: mspctr.f90:32
integer, parameter mxsite
Definition: nlsdim.f90:39
double precision, dimension(mxpt), save fvec
Definition: lmcom.f90:17
integer, save nser
Definition: parcom.f90:62
double factor
Definition: genio.c:44
integer, save lmflag
Definition: lmcom.f90:24
integer, save nspline
Definition: expdat.f90:45
subroutine sclstt(lu)
Definition: statc.f90:401
subroutine tdgstt(lu)
Definition: statc.f90:354
subroutine statc(line)
Definition: statc.f90:20
integer, dimension(mxspc), save idrv
Definition: expdat.f90:45
double precision, pointer, save xtol
Definition: lmcom.f90:29
Definition: lmcom.f90:13
integer, parameter luttyo
Definition: stdio.f90:29
double precision, dimension(mxspc), save sb0
Definition: expdat.f90:40
double precision, dimension(mxspc), save rmsn
Definition: expdat.f90:40
subroutine xpack(x, n)
Definition: fitl.f90:222
integer, save shftflg
Definition: expdat.f90:45
integer, dimension(mxspc), save npts
Definition: expdat.f90:45
integer, pointer, save maxitr
Definition: lmcom.f90:32
integer, parameter mxdim
Definition: nlsdim.f90:39
integer, save itrace
Definition: stdio.f90:32
integer, parameter mxtdg
Definition: nlsdim.f90:39
subroutine datstt(lu)
Definition: statc.f90:240
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67
integer, save ntd
Definition: tridag.f90:28
integer, save nbas
Definition: basis.f90:23
integer, dimension(mxtdm), save ltbas
Definition: basis.f90:23
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