NLSL
momdls.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 beta 11/23/95
2 
3 c*********************************************************************
4 c
5 c ==================
6 c subroutine:MOMDLS
7 c ==================
8 c
9 c>@brief Subroutine version of EPRLL family of programs.
10 c> This routine is intended for use with nonlinear least-squares
11 c> applications. The routine calculates a MOMD spectrum by
12 c> averaging spectra for a specified number of orientation using
13 c> the parameters passed in the fparm and iparm arrays. Spectra
14 c> are calculated using the Conjugate Gradients version of the Lanczos
15 c> algorithm.
16 c>
17 c> If the parameter specifying the number of MOMD orientations
18 c> is not greater than 1, or if there is no potential defined,
19 c> a single spectral calculation is carried out.
20 c
21 c Subroutine arguments are exactly output as follows:
22 c
23 c>@param al(mxdim) Diagonal of tridiagonal matrix for spectrum
24 c>@param be(mxdim) Off-diagonal of tridiagonal matrix for spectrum
25 c>
26 c>@param cgerr Residual error in CG solution vector
27 c>
28 c>@param ntotal Number of CG steps taken (total number for
29 c> all orientations in MOMD calculations)
30 c>
31 c>@param ierr Error flag
32 c Return values are defined in 'errmsg.f'
33 c
34 c
35 c>@author Written by DEB, 26-Sep-91 based on programs from DJS
36 c
37 c Includes :
38 c nlsdim.inc
39 c eprmat.inc
40 c eprprm.inc
41 c pidef.inc
42 c stdio.inc
43 c rndoff.inc
44 c
45 c Uses :
46 c eprls.f
47 c
48 c*********************************************************************
49 c
50  subroutine momdls( fparmi,iparmi,icalc,al,be,bss,iprune,
51  # spectr,work,nft,ntotal,ierr )
52 c
53  use eprprm
54  use nlsdim
55  use eprmat
56  use tridag
57 c use prmeqv
58  use errmsg
59  use pidef
60  use stdio
61  use rnddbl
62  use dfunc
63 c
64  implicit none
65  integer bss(5,mxdim),iparmi(niprm),ntotal,icalc,ierr,iprune,
66  # nft,isi
67 c double precision fparmi(NFPRM),spectr(iparmi(INFLD)),
68 c # work(iparmi(INFLD)),cgerr
69  double precision fparmi(nfprm),spectr(*),
70  # work(*),cgerr
71  double complex al(ntotal),be(ntotal)
72 c
73  integer i,id,ixt,j,nptr,nstp,nevl
74  double precision acc,cospsi,cosxi,dcos,dnorm,dwt,onorm,sinxi
75  logical init
76 c
77  double precision ZERO,ONE,SMALL,D180
78  parameter(zero=0.0d0,one=1.0d0,small=1.0d-16,d180=180.0d0)
79 c
80  double precision fu20,fu20phi
81  external fu20,fu20phi
82 c
83 c#####################################################################
84 c
85  init = mod(icalc,1000)/100 .ne. 0
86 c
87 c ----------------------------------------------
88 c Load values from parameter array into /eprprm/
89 c ----------------------------------------------
90 c
91 c *** In the future this could be done by select_site(isi) ***
92 c call select_site(isi)
93 c
94  do i=1,nfprm
95  fepr(i)=fparmi(i)
96  end do
97 
98  do i=1,niprm
99  iepr(i)=iparmi(i)
100  end do
101 c
102  if (init) then
103 c *** Fatal error in parameters
104  call lcheck(ierr)
105  if (ierr.ge.fatal) return
106 c
107 c -------------------------------------------------------
108 c If lcheck was not called, we still need to:
109 c (1) convert W tensor to Cartesian form and find w0
110 c (2) count potential coefficients
111 c (3) set director tilt flag
112 c -------------------------------------------------------
113  else
114  call tocart( fepr(iwxx), iwflg )
115  w0=(wxx+wyy+wzz)/3.0d0
116  ipt=0
117  if (dabs(c20).gt.rndoff) ipt=ipt+1
118  if (dabs(c22).gt.rndoff) ipt=ipt+1
119  if (dabs(c40).gt.rndoff) ipt=ipt+1
120  if (dabs(c42).gt.rndoff) ipt=ipt+1
121  if (dabs(c44).gt.rndoff) ipt=ipt+1
122  ipsi0=0
123  if( ((abs(psi).gt.rndoff).and.(abs(psi)-d180.gt.rndoff))
124  # .or. nort.gt.1) ipsi0=1
125  end if
126 c
127  do j=1,nfld
128  spectr(j)=zero
129  end do
130 c
131  if (nort.le.1 .or. ipt.le.0) then
132 c
133 c --------------------------
134 c --- Single orientation ---
135 c --------------------------
136 c
137  call eprls( icalc,al,be,bss,iprune,spectr,nft,ntotal,ierr)
138 c
139 c --------------------------
140 c ---- MOMD calculation ----
141 c --------------------------
142 c
143  else
144  ixt=1
145  nptr=0
146  onorm=one/dfloat(nort-1)
147 c
148 c ------------------------------------------------------
149 c Partial MOMD
150 c If there is ordering of the director, calculate
151 c normalization factor for the pseudopotential function
152 c In the case of partial director ordering, the parameter
153 c psi actually represents the angle between the director
154 c ordering axis and the field, whereas the calculation
155 c subroutine always interprets psi as the director tilt.
156 c Save the given psi value as the angle xi.
157 c ------------------------------------------------------
158 c
159  if (abs(dc20).gt.rndoff) then
160  xi=abs(mod(psi,d180))
161  cosxi=cos(radian*xi)
162  sinxi=sin(radian*xi)
163 c
164  acc=1.0d-4
165  dlam=dc20
166  call ccrint(zero,one,acc,small,dnorm,nevl,fu20,id)
167  else
168  dlam=zero
169  dnorm=one
170  end if
171 c
172 c -----------------------
173 c Loop over orientations
174 c -----------------------
175 c
176  cospsi=zero
177  dcos=one/(nort-1)
178 
179  do i=1,nort
180  cospsi=min(cospsi,one)
181  fepr(ipsi)=acos(cospsi)/radian
182 c
183  if (icalc.ne.0) then
184  nstp=nstep
185  if (nstep.le.0) nstp=min0(ndim,ntotal/nort)
186 c
187 c --------------------------------------------------
188 c Repeat continued-fractions for this orientation
189 c with the tridiagonal matrix: find number of steps
190 c by locating next zero in beta array
191 c --------------------------------------------------
192 c
193  else
194  nstp=0
195  3 nstp=nstp+1
196  nptr=nptr+1
197  if (nptr.gt.ntotal) then
198  ierr=tdgerr
199  return
200  end if
201  if ( abs(be(nptr)).gt.rndoff ) go to 3
202  end if
203 c
204  call eprls( icalc,al(ixt),be(ixt),bss,iprune,work,
205  # nft,nstp,ierr )
206  if (ierr.ge.fatal) return
207  ixt=ixt+nstp
208 c
209 c --------------------------------------------------
210 c If there is ordering of the director, calculate the
211 c weighting coefficient for this tilt angle and
212 c weight the spectrum
213 c --------------------------------------------------
214 c
215  if (abs(dc20).gt.rndoff) then
216  if (xi.gt.rndoff) then
217  ss=sinxi*sqrt(one-cospsi*cospsi)
218  cc=cosxi*cospsi
219  call ccrint(zero,pi,acc,small,dwt,nevl,fu20phi,id)
220  dwt=dwt/(dnorm*pi)
221  else
222  dwt=fu20(cospsi)/dnorm
223  end if
224 c
225  do j=1,nfld
226  work(j)=work(j)*dwt
227  end do
228  end if
229 c
230 c --------------------------------------------------------
231 c Trapezoidal rule: first and last pts. get factor of 1/2
232 c --------------------------------------------------------
233 c
234  if (i.eq.1 .or. i.eq.nort) then
235  do j=1,nfld
236  work(j)=work(j)*0.5d0
237  end do
238  end if
239 c
240 c ------------------------------------------------
241 c Add spectrum for current orientation into total
242 c ------------------------------------------------
243 c
244  do j=1,nfld
245  spectr(j)=spectr(j)+work(j)
246  end do
247 c
248  cospsi=cospsi+dcos
249  end do
250 c
251 c -------------------
252 c Normalize spectrum
253 c -------------------
254 c
255  do j=1,nfld
256  spectr(j)=spectr(j)*onorm
257  end do
258 c
259 c -------------------------------------------
260 c Return total number of steps actually taken
261 c -------------------------------------------
262 c
263  if (icalc.ne.0) ntotal=ixt-1
264 c
265  end if
266 c
267  return
268  end
269 
integer, save ipt
Definition: eprprm.f90:82
subroutine lcheck(ierr)
Definition: lcheck.f90:40
double precision, save w0
Definition: eprprm.f90:78
integer, pointer, save nfld
Definition: eprprm.f90:69
double precision, pointer, save dc20
Definition: eprprm.f90:58
integer, parameter niprm
Definition: nlsdim.f90:57
double precision, pointer, save c22
Definition: eprprm.f90:58
integer, pointer, save ndim
Definition: eprprm.f90:69
integer, dimension(niprm), target, save iepr
Definition: eprprm.f90:56
subroutine tocart(t, iflg)
Definition: tensym.f90:40
Definition: stdio.f90:26
integer, parameter tdgerr
Definition: errmsg.f90:51
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
integer, parameter ipsi
Definition: eprprm.f90:92
integer, pointer, save nort
Definition: eprprm.f90:69
double precision, save cc
Definition: dfunc.f90:19
double precision, save ss
Definition: dfunc.f90:19
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, dimension(nfprm), target, save fepr
Definition: eprprm.f90:55
double precision, pointer, save c44
Definition: eprprm.f90:58
double precision, pointer, save c40
Definition: eprprm.f90:58
double precision, pointer, save psi
Definition: eprprm.f90:58
double precision, pointer, save c20
Definition: eprprm.f90:58
integer, pointer, save nstep
Definition: eprprm.f90:69
integer, parameter iwxx
Definition: eprprm.f90:92
double precision, save xi
Definition: dfunc.f90:19
subroutine ccrint(bndlow, bndhi, epsiln, small, sum, neval, f, id)
Definition: ccrints.f90:32
integer, pointer, save iwflg
Definition: eprprm.f90:69
double precision, pointer, save c42
Definition: eprprm.f90:58
double precision, parameter radian
Definition: pidef.f90:15
double precision, parameter pi
Definition: pidef.f90:15
integer, parameter fatal
Definition: errmsg.f90:49
integer, parameter mxdim
Definition: nlsdim.f90:39
double precision, save dlam
Definition: dfunc.f90:19
Definition: dfunc.f90:16
double precision, pointer, save wyy
Definition: eprprm.f90:58
integer, save ipsi0
Definition: eprprm.f90:82
double precision, pointer, save wzz
Definition: eprprm.f90:58
double precision, pointer, save wxx
Definition: eprprm.f90:58
double precision, parameter rndoff
Definition: rnddbl.f90:86
integer, parameter nfprm
Definition: nlsdim.f90:57
Definition: pidef.f90:12