NLSL
eprls.f90
Go to the documentation of this file.
1 c NLSL Version 1.5.1 beta 11/25/95
2 c*********************************************************************
3 c ==================
4 c subroutine:EPRLS
5 c ==================
6 c
7 c>@brief Subroutine version of EPRLL family of programs by D. Schneider.
8 c> This routine is intended for use with nonlinear least-squares
9 c> applications. The routine calculates the tridiagonal matrix
10 c> for the parameters passed in the fparm and iparm arrays using
11 c> the Conjugate Gradients version of the Lanczos algorithm.
12 c>
13 c> Calculation parameters are input in the arrays fparm and iparm.
14 c> They are expected in the same order as that found in common
15 c> /eprprm/. The input flag icalc is specified as a three-digit
16 c> number ijk, where
17 c>
18 c> j.ne.0 => Perform matrix and CG calculations (matrll/cscg)
19 c> k.ne.0 => Perform starting vector and CG calculations (stvect/cscg)
20 c>
21 c> The continued-fractions calculation is always performed.
22 c
23 c Subroutine arguments are output as follows:
24 c
25 c> @param al(ndone) Diagonal of tridiagonal matrix for spectrum
26 c> @param be(ndone) Off-diagonal of tridiagonal matrix for spectrum
27 c> @param ndone Number of CG steps taken. Zero for Lanczos
28 c> calculations, and set negative if CG did
29 c> not converge
30 c> @param ierr Error flag. Meaning of values assigned to ierr
31 c> are defined in 'errmsg.inc'.
32 c>
33 c> Written by DEB, 26-Sep-91 based on programs from DJS
34 c
35 c Uses :
36 c pmatrl.f Liouville matrix calculation (pruning version)
37 c pstvec.f Starting vector calculation (pruning version)
38 c cscg.f Conjugate gradients tridiagonalization
39 c cfs.f Continued-fraction spectral calculation
40 c
41 c*********************************************************************
42 c
43  subroutine eprls( icalc,al,be,bss,iprune,spectr,nft,ndone,ierr)
44 c
45  use nlsdim
46  use eprprm
47  use eprmat
48  use tridag
49  use ftwork
50 c use prmeqv
51  use errmsg
52  use pidef
53  use stdio
54 c
55 c The parameter TOLEXP is used to determine the smallest number for
56 c which the Gaussian convolution function needs to be calculated,
57 c exp(-2*TOLEXP).
58 c
59  implicit none
60  double precision EIGHT,ONE,THIRD,TWO,ZERO,EPS
61  double complex CZERO
62  parameter(zero=0.0d0,one=1.0d0,third=0.333333333333333d0,
63  # two=2.0d0,eight=8.0d0,eps=1.0d-3,czero=(0.0d0,0.0d0))
64 c
65  integer bss(5,mxdim),ndone,icalc,ierr,iprune,nft
66 c double precision spectr(iepr(INFLD))
67  double precision spectr(*)
68  double complex al(mxdim),be(mxdim)
69  logical matrix,start
70 c
71  double precision cgerr,wline,df,f,g,gib,gnorm
72  integer i,m,n,no2,nstp
73 c
74  double precision dblint
75  external dblint
76 c
77 c#####################################################################
78 c
79  matrix = mod(icalc,100)/10 .ne. 0
80  start = mod(icalc,10) .ne. 0
81 c
82  if (matrix) then
83  if (iprune.ne.0) then
84  call pmatrl(bss,ierr)
85  else
86  call matrll(ierr)
87  end if
88 c
89  call wpoll
90  if (hltcmd.ne.0 .or. hltfit.ne.0 .or. ierr.gt.fatal) return
91  end if
92 c
93 c----------------------------------------------------------------------
94 c Generate starting vector
95 c----------------------------------------------------------------------
96  if (start) then
97  if (iprune.ne.0) then
98  call pstvec(bss,stv,ierr)
99  else
100  call stvect(stv,ierr)
101  endif
102  end if
103 c
104  if (start .or. matrix) then
105 c
106 c---------------------------------------------------------------------
107 c Call conjugate gradients tridiagonalization routine
108 c---------------------------------------------------------------------
109 c
110  do i=1,ndim
111  y(i)=czero
112  end do
113 c
114 c -----------------------------------------------------
115 c Set the number of steps for CG tridiagonalization
116 c If nstep is unspecified, use the matrix dimension or
117 c the number of steps previously taken, whichever is
118 c smaller
119 c -----------------------------------------------------
120  nstp=nstep
121  if (nstep.le.0) nstp=min(ndim,ndone)
122 c
123  call cscg( stv,ndim,nstp,cgtol,dcmplx(shiftr,shifti),
124  # y,al,be,ndone,cgerr )
125 c
126 c Check for user input or halt before proceeding
127 c
128  call wpoll
129  if (hltcmd.ne.0 .or. hltfit.ne.0) then
130  ierr=cghlt
131  return
132  end if
133 c
134 c *** CG did not converge
135  if (ndone.lt.0) then
136  ndone=-ndone
137  ierr=noconvrg
138  write(eprerr(ierr)(27:31),'(i5)') ndone
139  end if
140 c
141  end if
142 c
143 c ---------------------------------------------
144 c Continued-fraction calculation of spectrum
145 c ---------------------------------------------
146  call cfs( al,be,ndone,b0-fldi,-dfld,nfld,w0+lb,ideriv,phase,
147  # spectr )
148 c
149 c ------------------------------------
150 c Convolution with Gaussian lineshape
151 c ------------------------------------
152  gib=gib0+gib2*dsin( radian*psi )**2
153  wline=dsqrt(gib*gib+lb*lb)
154  call gconvl(spectr,wline,dfld,nfld,nft)
155 c
156 c ------------------
157 c Normalize spectrum (PI comes from Lorentzian shape)
158 c ------------------
159  do i=1,nfld
160  spectr(i)=spectr(i)/pi
161  end do
162 c
163  return
164 c
165  end
double complex, dimension(mxdim), save stv
Definition: tridag.f90:32
double precision, save w0
Definition: eprprm.f90:78
double precision, pointer, save dfld
Definition: eprprm.f90:58
void FORTRAN() wpoll()
Definition: pltx.c:860
integer, pointer, save nfld
Definition: eprprm.f90:69
double precision, pointer, save lb
Definition: eprprm.f90:58
integer, pointer, save ndim
Definition: eprprm.f90:69
integer, save hltcmd
Definition: stdio.f90:32
double complex, dimension(mxdim), save y
Definition: tridag.f90:32
double precision, pointer, save shifti
Definition: eprprm.f90:58
double precision, pointer, save cgtol
Definition: eprprm.f90:58
integer, save hltfit
Definition: stdio.f90:32
subroutine cscg(b, ndim, mxcgs, cgtol, shift, x, al, bl, ndone, error)
Definition: cscg.f90:68
Definition: stdio.f90:26
subroutine pmatrl(bss, ierr)
Definition: pmatrl.f90:79
double precision, pointer, save fldi
Definition: eprprm.f90:58
subroutine eprls(icalc, al, be, bss, iprune, spectr, nft, ndone, ierr)
Subroutine version of EPRLL family of programs by D. Schneider. This routine is intended for use with...
Definition: eprls.f90:44
double precision, pointer, save psi
Definition: eprprm.f90:58
integer, parameter cghlt
Definition: errmsg.f90:51
double precision, pointer, save gib0
Definition: eprprm.f90:58
double precision, pointer, save gib2
Definition: eprprm.f90:58
integer, pointer, save nstep
Definition: eprprm.f90:69
integer, pointer, save ideriv
Definition: eprprm.f90:69
subroutine pstvec(bss, v, ierr)
Definition: pstvec.f90:32
double precision, pointer, save phase
Definition: eprprm.f90:58
subroutine cfs(a, b, n, z0, dz, nz, w, ideriv, phs, val)
Definition: cfs.f90:43
subroutine matrll(ierr)
Definition: matrll.f90:63
double precision, parameter radian
Definition: pidef.f90:15
double precision, parameter pi
Definition: pidef.f90:15
integer, parameter noconvrg
Definition: errmsg.f90:51
subroutine stvect(v, ierr)
Definition: stvect.f90:33
integer, parameter fatal
Definition: errmsg.f90:49
character *50, dimension(neperr), save eprerr
Definition: errmsg.f90:59
double precision, pointer, save b0
Definition: eprprm.f90:58
integer, parameter mxdim
Definition: nlsdim.f90:39
double precision, pointer, save shiftr
Definition: eprprm.f90:58
subroutine gconvl(spectr, wline, dfld, nfld, nft)
Definition: gconvl.f90:15
Definition: pidef.f90:12