NLSL
lfun.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 beta 11/23/95
2 c----------------------------------------------------------------------
3 c> @file
4 c> This seems to be the central hub, which is called by the L-M
5 c> optimization, and which calls the functions that generate the
6 c> simulated spectrum
7 c =======================
8 c subroutine LFUN
9 c =======================
10 c
11 c> @brief Subroutine for interfacing EPRLL spectral calculations with the MINPACK
12 c> version of the Levenberg-Marquardt nonlinear least squares algorithm
13 c> for fitting experimental spectra. The function performed
14 c> function determined by the IFLAG argument.
15 c> MOMDLS subroutine to calculate a spectrum or series of spectra, each
16 c> of which consists of one or more spectral components, or "sites".
17 c>
18 c> @param iflag
19 c> @parblock
20 c> > IFLAG=0 (Printing of iterates)
21 c> >
22 c> > IFLAG=1 (Evaluation of residuals function)
23 c> >
24 c> > For a function evaluation, this routine takes the input vector X
25 c> > of search parameters, loads them into the fparm parameter arrays in
26 c> > /parcom/ (using additional indices available in /parcom/) and then
27 c> > performs spectral calculations by calling the MOMDLS subroutine.
28 c> >
29 c> > The experimental data to which the calculations are to be compared
30 c> > is contained in the common block /expdat/. LFUN optionally shifts
31 c> > the calculated spectra to minimize the squared differences between
32 c> > calculation and data, and then automatically determines a least-
33 c> > squares scale factor (or set of scale factors for multiple sites/spectra)
34 c> > which are returned in the vector sfac in common /mspctr/. The
35 c> > The unscaled spectra are returned in the matrix spectr in common
36 c> > /mspctr/. Each column of spectr corresponds to a single site;
37 c> > multiple spectra are stored sequentially in the columns of spectr.
38 c> >
39 c> > The residuals (differences between calculation and data after shifting
40 c> > and scaling) are scaled by the standard deviation of the estimated
41 c> > noise in each experimental spectrum, and returned to MINPACK through
42 c> > the fvec argument.
43 c> >
44 c> > IFLAG=2 (Jacobian evaluation)
45 c> >
46 c> > For evaluation of the Jacobian, the routine uses the forward-differences
47 c> > approximation to calculate the partial derivative of the spectrum with
48 c> > respect to each of the parameters being varied. The routine assumes
49 c> > that the individual spectra calculated at the point in parameter space
50 c> > where the Jacobian is being evaluated are contained in in the spectr
51 c> > matrix, and the appropriate scaling factors in the sfac array (both
52 c> > in common /mspctr/ ). The Jacobian matrix is returned in the fjac
53 c> > argument.
54 c> >
55 c> > IFLAG=3
56 c> > Fit is terminating. Output most recent set of parameter values and
57 c> > exit.
58 c> @endparblock
59 c Uses:
60 c setspc
61 c sshift
62 c sscale
63 c momdls
64 c fstplt (X-Windows interface)
65 c
66 c----------------------------------------------------------------------
67  subroutine lfun( m,n,x,fvec,fjac,ldfjac,iflag )
68 c
69  use nlsdim
70  use eprprm
71  use expdat
72  use tridag
73  use mspctr
74  use ftwork
75  use parcom
76  use basis
77  use errmsg
78  use pidef
79  use stdio
80 c use prmeqv
81  use iterat
82  use rnddbl
83 c
84  implicit none
85  integer m,n,iflag,ldfjac,ixbp
86 c
87  double precision x(n),fvec(m),fjac(ldfjac,n+nsite+nspc)
88  integer i,icalc,ierr,ise,isi,isp,ix,ixb,ixs,ixt,j,k,ld,lastsp,nj
89  double precision cgerr,fieldi,shift,snm,wdth,xtemp,dummy
90  character dashes*132,tmpnam*10,shftnd*1
91  logical shiftOK,glbscal,sppar
92 c
93  integer FULL,CFONLY
94  double precision ZERO
95  parameter(full=111,cfonly=0,zero=0.0d0)
96 c
97  integer itrim
98  double precision sshift,wtdres
99  logical spcpar,hltchk
100  external sshift,itrim,spcpar,hltchk,wtdres
101 c
102 c######################################################################
103 c
104  ierr=0
105  ld=0
106  glbscal=nsite.gt.1 .and. nspc.gt.1
107  if (iter.gt.1) warn=.false.
108 c
109 c**********************************************************************
110 c**********************************************************************
111 c IFLAG=0: Output iteration information
112 c**********************************************************************
113 c**********************************************************************
114 c
115 c
116 c ------------------------------------------
117 c Draw a line of dashes for first iteration
118 c ------------------------------------------
119  if (iflag.eq.0) then
120  if (iter.eq.1) then
121  ld=12*n+16
122  if (ld.gt.132) ld=132
123  do i=1,ld
124  dashes(i:i)='-'
125  end do
126 c
127 c -----------------------------------
128 c Output names of current parameters
129 c -----------------------------------
130 c
131  if (luout.ne.luttyo) then
132  if (ld.gt.0)write (luttyo,1006) dashes(:ld)
133  if (iwflag.ne.0) then
134  write (luttyo,1007) (tag(i)(:itrim(tag(i))),i=1,n)
135  else
136  write (luttyo,1010) (tag(i)(:itrim(tag(i))),i=1,n)
137  end if
138  if (ld.gt.0)write (luttyo,1008) dashes(:ld)
139  end if
140 c
141  if (ld.gt.0)write (luout,1006) dashes(:ld)
142  if (iwflag.ne.0) then
143  write (luout,1007) (tag(i)(:itrim(tag(i))),i=1,n)
144  else
145  write (luout,1010) (tag(i)(:itrim(tag(i))),i=1,n)
146  end if
147  if (ld.gt.0)write (luout,1008) dashes(:ld)
148  end if
149 c
150  if (ishglb.ne.0 .and. (nshift.eq.0 .or.iter.le.nshift) ) then
151  shftnd='S'
152  else
153  shftnd=' '
154  end if
155 c
156  if (iwflag.ne.0) then
157  rdchsq=fnorm*fnorm/float(m-n)
158  write (luout,1009) shftnd,iter,rdchsq,(x(i),i=1,n)
159  if (luout.ne.luttyo)
160  # write (luttyo,1009) shftnd,iter,rdchsq,(x(i),i=1,n)
161  else
162  write (luout,1009) shftnd,iter,fnorm,(x(i),i=1,n)
163  if (luout.ne.luttyo)
164  # write (luttyo,1009) shftnd,iter,fnorm,(x(i),i=1,n)
165  end if
166 c
167 c ------------------------------------------------------------
168 c Output any other requested information regarding iterates
169 c ------------------------------------------------------------
170  if (iitrfl.ne.0) call wrfun(iter)
171  if (itridg.ne.0) call wrtrdg(iter)
172  if (jacobi.ne.0) call wrjac(iter)
173 c
174 c**********************************************************************
175 c**********************************************************************
176 c IFLAG=1: Function evaluation
177 c**********************************************************************
178 c**********************************************************************
179 c
180  else if (iflag.eq.1) then
181 c
182  xreset=.false.
183  written=0
184 c
185  do i=1,n
186  call setprm(ixpr(i),ixst(i),x(i))
187  end do
188 c
189  shiftok = (nshift.eq.0 .or. iter.le.nshift)
190  call tdsqz()
191 c
192 c -----------------------------------
193 c Loop over all spectra in a series
194 c -----------------------------------
195 c
196  ixs=1
197  do ise=1,nser
198  ixsp(ise)=ixs
199  tmpshft(ise)=zero
200 c
201 c --------------------------------------
202 c Loop over all sites for the spectrum
203 c --------------------------------------
204 c
205  do isi=1,nsite
206  call tdchek(isi,ise,ierr)
207  call setspc(isi,ise)
208  if (hltchk(ierr,isi,ise,iflag)) return
209 c
210  ixt=ixtd(isi,ise)
211  if(basno(isi,ise).gt.0) then
212  ixb=ixbas( basno(isi,ise) )
213  ixbp=ixb
214  else
215  ixb=0
216  ixbp=1
217  endif
218  icalc=cfonly
219  if (modtd(isi,ise).ne.0) icalc=full
220 c
221 c ---------------------------------
222 c Calculate an individual spectrum
223 c ---------------------------------
224 c
225  call momdls( fparm(1,isi),iparm(1,isi),icalc,
226  # alpha(ixt),beta(ixt),ibasis(1,ixbp),ixb,
227  # spectr(ixs,isi),wspec,nft(ise),ltd(isi,ise),
228  # ierr )
229 c
230  modtd(isi,ise)=0
231  if (hltchk(ierr,isi,ise,iflag)) return
232  end do
233 c
234 c
235 c --------------------------------------------------------
236 c If data are available for this calculation, determine
237 c optimal shift and scaling and calculate fvec accordingly
238 c (This feature permits calculations w/out data)
239 c --------------------------------------------------------
240  if (ise.le.nspc) then
241 c
242 c ----------------------------------------------------------
243 c If shifting is enabled, find the shift that will give
244 c optimal overlap between a calculated spectrum (or set of
245 c spectra) and the data for each element of a series
246 c -----------------------------------------------------------
247  if (shiftok .and. ishft(ise).ne.0) then
248 c
249  shift=sshift( data(ixs),spectr(ixs,1),mxpt,nsite,
250  # npts(ise),nft(ise),idrv(ise),srange,ctol,
251  # noneg,tmpdat,tmpclc,wspec,work,sfac(1,ise) )
252  tmpshft(ise)=shift*sdb(ise)
253 c
254 c --------------------------------------------------------
255 c Re-calculate shifted spectrum with continued-fractions
256 c (this will keep the frwrd.-diff. derivative reasonable)
257 c --------------------------------------------------------
258  if (dabs(shift).gt.1.0d-3*sdb(ise)) then
259  do isi=1,nsite
260  icalc=cfonly
261  call setspc(isi,ise)
262  if (hltchk(ierr,isi,ise,iflag)) return
263 c
264  ixt=ixtd(isi,ise)
265  if (basno(isi,ise).gt.0) then
266  ixb=ixbas( basno(isi,ise) )
267  ixbp=ixb
268  else
269  ixb=0
270  ixbp=1
271  endif
272 c
273  call momdls( fparm(1,isi),iparm(1,isi),icalc,
274  # alpha(ixt),beta(ixt),ibasis(1,ixbp),ixb,
275  # spectr(ixs,isi),wspec,nft(ise),
276  # ltd(isi,ise),ierr )
277 c
278  if (hltchk(ierr,isi,ise,iflag)) return
279  end do
280  end if
281 c
282 c ----------------------------------------------------------
283 c For the current experimental spectrum calculate residuals
284 c from shifted spectra and scale factors (unless sites
285 c scales are to be determined globally)
286 c ----------------------------------------------------------
287  if (.not. glbscal) then
288  do j=1,npts(ise)
289  k=ixs+j-1
290  fvec(k)=data(k)
291  do isi=1,nsite
292  fvec(k)=fvec(k)-sfac(isi,ise)*spectr(k,isi)
293  end do
294  end do
295  end if
296 c
297 c ---------------------------------
298 c No shifting: use sscale routine
299 c ---------------------------------
300  else
301 c
302 c -----------------------------------------------------------
303 c Find least-squares site scaling factors and return
304 c residuals in fvec; however, multiple site/multiple spectra
305 c problems must be scaled globally, outside loop over spectra
306 c -----------------------------------------------------------
307  if (.not. glbscal)
308  # call sscale( data(ixs),spectr(ixs,1),wspec,mxpt,
309  # work,nsite,npts(ise),ctol,noneg,
310  # iscal,sfac(1,ise),fvec(ixs) )
311 c
312 c *** End of "if (shifting allowed)...else..."
313  end if
314 c
315 c *** End of "if (ise.le.nspc)"
316  else
317  do j=1,npts(ise)
318  k=ixs+j-1
319  fvec(k)=zero
320  do isi=1,nsite
321  fvec(k)=fvec(k)-sfac(isi,ise)*spectr(k,isi)
322  end do
323  end do
324 
325  end if
326 
327  ixs=ixs+npts(ise)
328 
329 c *** End of loop over spectra in a series
330  end do
331 c
332 c --------------------------------------------------------------
333 c For multiple-site multiple spectral fits, the scaling factors
334 c are determined by a global least-squares fit to all spectra
335 c in the series
336 c --------------------------------------------------------------
337  if (glbscal) then
338 c
339  call sscale( data,spectr,wspec,mxpt,work,nsite,m,ctol,
340  # noneg,iscal,sfac,fvec )
341 c
342 c --------------------------------------------
343 c Copy the global scale factors to all spectra
344 c --------------------------------------------
345  do ise=2,nser
346  do isi=1,nsite
347  sfac(isi,ise)=sfac(isi,1)
348  end do
349  end do
350 c
351  end if
352 c
353 c ------------------------------
354 c Plot the current calculation
355 c ------------------------------
356  do i=1,nser
357 
358  call fstplt( data(ixsp(i)),fvec(ixsp(i)),
359  # sbi(i)-shft(i)-tmpshft(i),sdb(i),npts(i),i )
360  end do
361 c
362 c%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
363 c New calling sequence when pltd is implemented:
364 c
365 c do i=1,nser
366 c ixs=ixsp(i)
367 c call pltwin( data(ixs),fvec(ixs),spectr(ixs,1),sfac(1,i),
368 c * MXPT,npts(i),nsite,i)
369 c end do
370 c%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
371 c
372 c ----------------------------------------------------------
373 c For weighted residual fit, normalize differences in fvec
374 c to the spectral noise (then |fvec|^2 will approximate
375 c the true chi^2)
376 c ----------------------------------------------------------
377  if (iwflag.ne.0) then
378  do isp=1,nspc
379  ixs=ixsp(isp)
380  do j=1,npts(isp)
381  k=ixs+j-1
382  fvec(k)=fvec(k)/rmsn(isp)
383  end do
384  end do
385  end if
386 c
387  if (xreset) call xpack(x,n)
388 c
389 c**********************************************************************
390 c**********************************************************************
391 c IFLAG = 2: Jacobian evaluation by forward-difference approximation:
392 c
393 c fjac = (f(x0+dx)-f(x0))/dx
394 c
395 c Note that f(x) = data-calc(x), so that f(x0+dx)-f(x0) = calc(x0)-calc(x0+dx)
396 c
397 c**********************************************************************
398 c**********************************************************************
399 c
400  else if (iflag.eq.2) then
401 c
402 c ------------------------------------------------
403 c Update shifts at the beginning of an iteration
404 c ------------------------------------------------
405  do isp=1,nspc
406  shft(isp)=shft(isp)+tmpshft(isp)
407  tmpshft(isp)=zero
408  end do
409 c
410 c -------------------------------------------------------------
411 c Loop over all parameters, introducing the forward-difference
412 c step into each parameter
413 c -------------------------------------------------------------
414  do ix=1,n
415  xtemp=x(ix)
416  x(ix)=x(ix)+xfdstp(ix)
417  call setprm(ixpr(ix),ixst(ix),x(ix))
418 c
419 c ------------------------------------
420 c Loop over all spectra in a series
421 c ------------------------------------
422  do isp=1,nspc
423  ixs=ixsp(isp)
424 c
425 c ---------------------
426 c Initialize Jacobian
427 c ---------------------
428  do j=1,npts(isp)
429  k=ixs+j-1
430  fjac(k,ix)=zero
431  end do
432  covarok=.false.
433 c
434 c ----------------------------------------
435 c Loop over all sites in the spectrum
436 c ----------------------------------------
437  do isi=1,nsite
438 c
439 c ---------------------------------------------------------
440 c Check whether this site/spectrum depends upon the
441 c given parameter before calculating the partial derivative
442 c ---------------------------------------------------------
443  sppar=spcpar( ixpr(ix) )
444  if ( (abs(sfac(isi,isp)).gt.rndoff) .and.
445  # (ixst(ix).le.0 .or.
446  # (ixst(ix).eq.isp .and. sppar) .or.
447  # (ixst(ix).eq.isi .and. .not. sppar) ) ) then
448 c
449  icalc=cfonly
450  if (modtd(isi,isp).ne.0) icalc=full
451 c
452  call setspc(isi,isp)
453  if (hltchk(ierr,isi,isp,iflag)) return
454 c
455  ixt=ixtd(isi,isp)
456  if(basno(isi,isp).gt.0) then
457  ixb=ixbas( basno(isi,isp) )
458  ixbp=ixb
459  else
460  ixb=0
461  ixbp=1
462  endif
463 c
464  call momdls( fparm(1,isi),iparm(1,isi),icalc,
465  # alpha(ixt),beta(ixt),ibasis(1,ixbp),ixb,
466  # work,wspec,nft(isp),ltd(isi,isp),ierr )
467 c
468  if (hltchk(ierr,isi,isp,iflag)) return
469 c
470 c ------------------------------------------
471 c Add in difference spectrum for this site
472 c ------------------------------------------
473  if (iwflag.ne.0) then
474 c --- weighted ---
475  do j=1,npts(isp)
476  k=ixs+j-1
477  fjac(k,ix)=fjac(k,ix)+(spectr(k,isi)-work(j))
478  # *sfac(isi,isp)/rmsn(isp)
479  end do
480 c --- unweighted ---
481  else
482  do j=1,npts(isp)
483  k=ixs+j-1
484  fjac(k,ix)=fjac(k,ix)+(spectr(k,isi)-work(j))
485  # *sfac(isi,isp)
486  end do
487 c
488  end if
489 c
490  endif
491 c
492 c *** end of loop over sites in spectrum
493  end do
494 c
495 c *** end of loop over spectra in series
496  end do
497 c
498 c -------------------------------------------------------
499 c Use the step size to calculate the partial derivative
500 c w.r.t. x(ix) for the current spectrum
501 c -------------------------------------------------------
502  do j=1,m
503  fjac(j,ix)=fjac(j,ix)/xfdstp(ix)
504  end do
505 c
506  x(ix)=xtemp
507  call setprm(ixpr(ix),ixst(ix),x(ix))
508 c
509 c *** end of loop over parameters
510  end do
511 c
512 c ------------------------------------------------------------
513 c Place unscaled spectra in last <nsite> columns of Jacobian
514 c matrix. (These are the partial derivative of the spectrum
515 c with respect to each of the scaling factors)
516 c These columns will not be allowed to pivot, but will be used
517 c to calculate the covariance matrix of the scale factors.
518 c ------------------------------------------------------------
519  if (nsite.gt.1) then
520  do isi=1,nsite
521  nj=n+isi
522  do isp=1,nspc
523  ixs=ixsp(isp)
524  do j=1,npts(isp)
525  k=ixs+j-1
526  fjac(k,nj)=spectr(k,isi)
527  if (iwflag.ne.0) fjac(k,nj)=fjac(k,nj)/rmsn(isp)
528  end do
529  end do
530  end do
531  else
532  do isp=1,nspc
533  nj=n+isp
534  lastsp=ixsp(isp)+npts(isp)-1
535  do j=1,ixsp(isp)-1
536  fjac(j,nj)=zero
537  end do
538  do j=ixsp(isp),lastsp
539  fjac(j,nj)=spectr(j,1)
540  if (iwflag.ne.0) fjac(j,nj)=fjac(j,nj)/rmsn(isp)
541  end do
542  do j=lastsp+1,m
543  fjac(j,nj)=zero
544  end do
545  end do
546  end if
547 
548 c
549 c**********************************************************************
550 c**********************************************************************
551 c IFLAG=3
552 c Fit has terminated: print result and save x-vector
553 c**********************************************************************
554 c**********************************************************************
555 c
556  else if (abs(iflag).eq.3) then
557 
558 c --------------------------------------------------------
559 c Print final parameters if termination is by convergence
560 c --------------------------------------------------------
561  if (iflag.lt.0) then
562 
563  if (iwflag.ne.0) then
564  rdchsq=fnorm*fnorm/float(m-n)
565  write (luout,1009) 'T',iter,rdchsq,(x(i),i=1,n)
566  if (luout.ne.luttyo)
567  # write (luttyo,1009) 'T',iter,rdchsq,(x(i),i=1,n)
568  else
569  write (luout,1009) 'T',iter,fnorm,(x(i),i=1,n)
570  if (luout.ne.luttyo)
571  # write (luttyo,1009) 'T',iter,fnorm,(x(i),i=1,n)
572  end if
573 
574  end if
575 c
576 c ----------------------------------------
577 c Store parameters currently in x-vector
578 c ----------------------------------------
579  do i=1,n
580  call setprm(ixpr(i),ixst(i),x(i))
581  end do
582 c
583  if (ld.gt.0) then
584  write (luttyo,1008) dashes(:ld)
585  if (luout.ne.luttyo) write (luout,1008) dashes(:ld)
586  end if
587  end if
588 
589 c
590  return
591 c
592 c ### format statements ##############################################
593 c
594  1006 format(/a)
595  1007 format(' Itr RedChSq ',10(1x,a9))
596  1008 format(a)
597  1009 format(a1,i3,2x,g10.5,1x,10(f10.5))
598  1010 format(' Itr RsdNorm ',10(1x,a9))
599  end
600 
601 
602 
603 c----------------------------------------------------------------------
604 c =========================
605 c function HLTCHK
606 c =========================
607 c
608 c> @brief Check whether a user halt (control-C) or other error
609 c> has occurred during the spectral calculation and set the
610 c> iflag error flag for LMNLS accordingly. Returns .true.
611 c> for fatal errors or user halts, and .false. otherwise.
612 c
613 c----------------------------------------------------------------------
614  function hltchk(ierr,isite,ispec,iflag)
615 c
616  use stdio
617  use errmsg
618 c
619  implicit none
620  integer ierr,iflag,isite,ispec
621  logical hltchk
622 c
623  hltchk=.false.
624  if (ierr.eq.0) return
625 c
626 c ------------------------------------
627 c Non-fatal errors: warn if requested
628 c ------------------------------------
629  if (ierr.lt.fatal) then
630  if (warn) then
631  write (luout,1000) 'Warning',isite,ispec,eprerr(ierr)
632  if (luout.ne.luttyo) write (luttyo,1000) 'Warning',
633  # isite,ispec,eprerr(ierr)
634  warn=.false.
635  end if
636 c
637 c ------------------------------
638 c Fatal errors (non-user related)
639 c ------------------------------
640  else if (ierr.ge.fatal .and. ierr.ne.mtxhlt
641  # .and. ierr.ne.cghlt) then
642 
643  write (luout,1000) 'Fatal err',isite,ispec,eprerr(ierr)
644  if (luout.ne.luttyo) then
645  write (luttyo,1000) 'Fatal err',isite,ispec,eprerr(ierr)
646  end if
647  hltchk=.true.
648  iflag=-1
649 c
650 c ------------------------------
651 c Other (user halt)
652 c ------------------------------
653  else
654  write (luout,1001) eprerr(ierr)
655  if (luout.ne.luttyo) write (luttyo,1001) eprerr(ierr)
656  hltchk=.true.
657  iflag=-1
658  end if
659 c
660  return
661 c
662 c----------------------------------------------------------------------
663 c
664  1000 format(/2x,'*** ',a,' site',i2,' spctrm',i2,': ',a)
665  1001 format(/2x,'*** ',a)
666  end
667 
subroutine tdchek(isite, ispc, ierr)
Definition: tdchek.f90:12
integer, save luout
Definition: stdio.f90:32
integer, dimension(mxvar), save ixst
Definition: parcom.f90:62
integer, parameter mxpt
Definition: nlsdim.f90:39
double precision, dimension(mxspc), save sbi
Definition: expdat.f90:40
double precision, dimension(mxpt), save tmpdat
Definition: ftwork.f90:26
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 noneg
Definition: parcom.f90:62
double precision, save srange
Definition: parcom.f90:56
double complex, dimension(mxtdg), save alpha
Definition: tridag.f90:32
double precision, dimension(mxspc), save sdb
Definition: expdat.f90:40
subroutine setspc(isite, ise)
Definition: setspc.f90:31
void FORTRAN() fstplt(double *y1, double *y2, double *xmin1, double *xstep1, int *indx1, int *wnum)
Definition: pltx.c:805
integer, dimension(5, mxdim), save ibasis
Definition: basis.f90:23
subroutine wrtrdg(iter)
Definition: writr.f90:212
integer, dimension(mxspc), save nft
Definition: expdat.f90:45
integer, dimension(mxspc), save ixsp
Definition: expdat.f90:45
integer, save iitrfl
Definition: parcom.f90:62
integer, save itridg
Definition: parcom.f90:62
integer, dimension(mxspc), save ishft
Definition: expdat.f90:45
integer, dimension(mxsite, mxspc), save basno
Definition: basis.f90:23
double precision, dimension(mxpt, mxsite), save wspec
Definition: mspctr.f90:33
double precision, dimension(mxspc), save tmpshft
Definition: expdat.f90:40
integer, dimension(mxsite, mxspc), save modtd
Definition: tridag.f90:28
integer, save iter
Definition: iterat.f90:14
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
subroutine wrjac(iter)
Definition: writr.f90:167
Definition: stdio.f90:26
double precision, dimension(nfprm, mxsite), target, save fparm
Definition: parcom.f90:54
subroutine momdls(fparmi, iparmi, icalc, al, be, bss, iprune, spectr, work, nft, ntotal, ierr)
Subroutine version of EPRLL family of programs. This routine is intended for use with nonlinear least...
Definition: momdls.f90:52
Definition: basis.f90:19
double precision, save fnorm
Definition: iterat.f90:15
integer, save iwflag
Definition: iterat.f90:14
subroutine tdsqz()
Definition: tdsqz.f90:13
integer, save jacobi
Definition: parcom.f90:62
double precision, dimension(mxvar), save xfdstp
Definition: parcom.f90:56
integer, dimension(mxvar), save ixpr
Definition: parcom.f90:62
integer, save written
Definition: expdat.f90:45
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 nshift
Definition: parcom.f90:62
integer, parameter cghlt
Definition: errmsg.f90:51
integer, dimension(niprm, mxsite), target, save iparm
Definition: parcom.f90:60
integer, dimension(mxsite), save iscal
Definition: mspctr.f90:32
integer, save nser
Definition: parcom.f90:62
subroutine wrfun(iter)
Definition: writr.f90:28
integer, dimension(mxspc), save idrv
Definition: expdat.f90:45
integer, dimension(mxsite, mxspc), save ixtd
Definition: tridag.f90:28
integer, parameter luttyo
Definition: stdio.f90:29
logical, save covarok
Definition: iterat.f90:18
double precision, save ctol
Definition: parcom.f90:56
logical function hltchk(ierr, isite, ispec, iflag)
Check whether a user halt (control-C) or other error has occurred during the spectral calculation and...
Definition: lfun.f90:615
double precision, dimension(mxspc), save rmsn
Definition: expdat.f90:40
logical, save warn
Definition: stdio.f90:33
integer, parameter fatal
Definition: errmsg.f90:49
integer, parameter mtxhlt
Definition: errmsg.f90:51
logical, save xreset
Definition: iterat.f90:18
subroutine xpack(x, n)
Definition: fitl.f90:222
character *50, dimension(neperr), save eprerr
Definition: errmsg.f90:59
integer, dimension(mxspc), save npts
Definition: expdat.f90:45
subroutine setprm(ixparm, ixsite, fval)
This file contains two routines that set a given parameter, specified by an index into the fparm or i...
Definition: setprm.f90:34
double complex, dimension(mxtdg), save beta
Definition: tridag.f90:32
integer, save ishglb
Definition: expdat.f90:45
double precision, dimension(mxpt), save tmpclc
Definition: ftwork.f90:26
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67
subroutine sscale(data, spct, wspct, ldspct, work, nsite, npt, ctol, noneg, iscal, sfac, resid)
Definition: sscale.f90:45
double precision, parameter rndoff
Definition: rnddbl.f90:86
integer, dimension(mxtdm), save ixbas
Definition: basis.f90:23
double precision, dimension(mxpt), save work
Definition: ftwork.f90:26
Definition: pidef.f90:12
double precision, dimension(mxpt, mxsite), save spectr
Definition: mspctr.f90:33