NLSL
lcheck.f90
Go to the documentation of this file.
1 c Version 1.5.1 abeta 03/07/96
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine LCHECK
5 c =========================
6 c
7 c This procedure replaces the function of program LBLL for
8 c performing initial computations based on the input parameters
9 c for the EPR lineshape programs. It loads the floating point and
10 c integer variables in the fparm and iparm arrays into common block
11 c /eprprm/ and checks to make sure they are valid.
12 c
13 c Returns ierr=0 if parameters are valid; a negative number
14 c identifying a fatal error, and a positive number indicating
15 c other conditions.
16 c
17 c lumsg specifies a logical unit for output of any error messages.
18 c None are produced if lumsg=0.
19 c
20 c Note: Check of MTS indices that was originally coded here
21 c is now in file mtschk
22 c
23 c Fatal parameter errors
24 c 1) zero gxx, gyy, or gzz
25 c 2) zero B0
26 c
27 c Includes:
28 c nlsdim.inc
29 c maxl.inc
30 c rndoff.inc
31 c eprprm.inc
32 c
33 c Uses:
34 c ipar
35 c cd2km
36 c sp2car
37 c
38 c----------------------------------------------------------------------
39  subroutine lcheck(ierr)
40 c
41  use nlsdim
42  use maxl
43  use eprprm
44 c use prmeqv
45  use errmsg
46  use rnddbl
47 c
48  implicit none
49  integer ierr
50 c
51  integer i,itmp,j
52  integer inlemx,inlomx,inkmn,inkmx,inmmn,inmmx,inpnmx
53  double precision gmax,gmin,hmax,davg,dn
54  double precision d2km(2,5,5)
55  logical axiala,axialg,axialw,axialr
56 c
57  double precision ZERO,HALF,DSQ23,ONE,TEN
58  parameter(zero=0.0d0,half=0.5d0,one=1.0d0,ten=1.0d1)
59  parameter(dsq23=0.816496580927726d0)
60 c
61  integer ipar
62  external ipar
63 c
64 
65 c######################################################################
66 c
67  ierr=0
68 c
69 c----------------------------------------------------------------------
70 c initialize spherical tensors
71 c----------------------------------------------------------------------
72  do i=1,5
73  faa(i)=zero
74  fgm(i)=zero
75  fwm(i)=zero
76  do j=1,2
77  fam(j,i)=zero
78  fgd(j,i)=zero
79  fad(j,i)=zero
80  fwd(j,i)=zero
81  end do
82  end do
83 
84 c----------------------------------------------------------------------
85 c Check for zero field
86 c----------------------------------------------------------------------
87  if (b0 .lt. zero) then
88  ierr=posb0
89  b0=abs(b0)
90  elseif (b0 .eq. zero) then
91  ierr=zerob0
92  return
93  endif
94 c----------------------------------------------------------------------
95 c Check for nonzero number of points and field range
96 c----------------------------------------------------------------------
97  if (nfld.le.0) then
98  ierr=zeronfld
99  return
100  end if
101 c
102 c----------------------------------------------------------------------
103 c Convert g, A, and R tensors to Cartesian form
104 c----------------------------------------------------------------------
105  call tocart( fepr(igxx), igflg )
106  call tocart( fepr(iaxx), iaflg )
107  call tocart( fepr(iwxx), iwflg )
108  call tocart( fepr(idx), irflg )
109 c
110 c *** NOTE: in NLS version, diffusion rates are specified as
111 c a power of 10
112 c
113  dx=ten**dx
114  dy=ten**dy
115  dz=ten**dz
116  axialr=abs(dx-dy).le.rndoff
117 c
118 c----------------------------------------------------------------------
119 c Check for zero values in g-tensor and zero A tensor
120 c----------------------------------------------------------------------
121  if ( abs(gxx) .lt. rndoff .or.
122  # abs(gyy) .lt. rndoff .or.
123  # abs(gzz) .lt. rndoff) then
124  ierr=zerog
125  return
126  endif
127 
128  if ( (in2 .le. 0) .or.
129  # (abs(axx) .lt. rndoff .and.
130  # abs(ayy) .lt. rndoff .and.
131  # abs(azz) .lt. rndoff) ) then
132  in2=0
133  axx=zero
134  ayy=zero
135  azz=zero
136  ierr=zeroin2
137  endif
138 
139 c----------------------------------------------------------------------
140 c Calculate spherical components of F_g and F_A
141 c----------------------------------------------------------------------
142 c *** Electronic Zeeman ***
143  g0=(gxx+gyy+gzz)/3.0d0
144  fgm(1)=half*(gxx-gyy)*(b0/g0)
145  fgm(3)=dsq23*(gzz-half*(gxx+gyy))*(b0/g0)
146  fgm(5)=fgm(1)
147  axialg=abs(fgm(1)).lt.rndoff
148 c *** Hyperfine ***
149  a0=(axx+ayy+azz)/3.0d0
150  faa(1)=half*(axx-ayy)
151  faa(3)=dsq23*(azz-half*(axx+ayy))
152  faa(5)=faa(1)
153  axiala=abs(faa(1)).lt.rndoff
154 c *** Line-broadening ***
155  w0=(wxx+wyy+wzz)/3.0d0
156  fwm(1)=half*(wxx-wyy)
157  fwm(3)=dsq23*(wzz-half*(wxx+wyy))
158  fwm(5)=fwm(1)
159  axialw=abs(fwm(1)).lt.rndoff
160 c
161  gmin=dmin1( gxx, gyy, gzz )
162  gmax=dmax1( gxx, gyy, gzz )
163  hmax=dmax1( abs(axx), abs(ayy), abs(azz) )
164 c
165 c----------------------------------------------------------------------
166 c Issue warning if high-field approximation has been violated
167 c (this is a very rough criterion for the high-field approx.!)
168 c----------------------------------------------------------------------
169 
170  if (b0 .lt. 10.0d0*dmax1( hmax, (gmax-gmin)*b0/g0) ) then
171  ierr=badhifld
172  endif
173 
174 c----------------------------------------------------------------------
175 c Set ipt and cpot array according to the potential coefficients given
176 c----------------------------------------------------------------------
177  ipt=0
178  do j=1,5
179  do i=1,5
180  cpot(i,j)=zero
181  end do
182  end do
183 
184  if (abs(c20) .gt. rndoff) then
185  ipt=ipt+1
186  cpot(2,1)=c20
187  lptmx=2
188  endif
189  if (abs(c22) .gt. rndoff) then
190  ipt=ipt+1
191  cpot(2,2)=c22
192  lptmx=2
193  kptmx=2
194  endif
195  if (abs(c40) .gt. rndoff) then
196  ipt=ipt+1
197  cpot(3,1)=c40
198  lptmx=4
199  endif
200  if (abs(c42) .gt. rndoff) then
201  ipt=ipt+1
202  cpot(3,2)=c42
203  lptmx=4
204  kptmx=2
205  endif
206  if (abs(c44) .gt. rndoff) then
207  ipt=ipt+1
208  cpot(3,3)=c44
209  lptmx=4
210  kptmx=4
211  endif
212 c
213  if (lptmx.ge.2) then
214  lband=lptmx*2
215  else
216  lband=2
217  end if
218 c
219  kband=kptmx*2
220  if (axialr) kband=kband+2
221 
222 c----------------------------------------------------------------------
223 c Check for consistency of specified motion parameters with given model
224 c----------------------------------------------------------------------
225  if (ipdf.gt.2) ipdf=2
226  if (ipdf.lt.0) ipdf=0
227  if (ist.lt.0) ist=0
228 c
229 c Set exponents according to non-Brownian model
230 c
231  expl=one
232  expkxy=one
233  expkzz=one
234  if (ml.eq.1) expl=half
235  if (mxy.eq.1) expkxy=half
236  if (mzz.eq.1) expkzz=half
237 c
238 c
239 c *** Non-Brownian motion: must have at least one nonzero residence
240 c time, and no potential specified
241 c
242 c
243 c Set residence times to zero for Brownian model
244 c
245  if (ipdf .eq. 0) then
246  pml = zero
247  pmxy = zero
248  pmzz = zero
249 c
250 c NOTE: in NLS version, R*residence time products are specified as
251 c powers of ten
252 c
253  else
254  pml = ten**pml
255  pmxy = ten**pmxy
256  pmzz = ten**pmzz
257 c
258 c *** NOTE: in NLS version, djf, djfprp specified as powers of ten
259 c
260  djf=ten**djf
261  djfprp=ten**djfprp
262 c
263  if (ipt .gt. 0) then
264  ierr=badjmp
265  ipdf=0
266  endif
267  endif
268 c
269 c *** Anisotropic viscosity model: must have potential and
270 c no discrete jump motion specified
271 c
272  if (ipdf .eq. 2) then
273  if (ipt .eq. 0) then
274  ierr=badav
275  ipdf=0
276  endif
277 c
278  if (ist .gt. 0) then
279  ierr=baddj
280  ipdf=0
281  endif
282  endif
283 c
284  if (oss.ne.zero) oss=ten**oss
285 c
286  ipsi0=0
287  if( ((abs(psi).gt.rndoff).and.(abs(psi-180.0d0).gt.rndoff))
288  # .or. nort.gt.1) ipsi0=1
289 c
290 c----------------------------------------------------------------------
291 c Set magnetic tilt flag and A tensor in g frame
292 c----------------------------------------------------------------------
293 c
294 c --- mag. tilt angle gamma not needed if g,W-matrices are axial
295 c
296 c if (axialg .and. axialw .and. abs(gam).gt.RNDOFF) then
297 c gam=ZERO
298 cc ierr=NOALPHAM
299 c end if
300 c
301  if (abs(alm).gt.rndoff .or. bem.gt.rndoff
302  # .or. abs(gam).gt.rndoff) itm=1
303 c
304 c
305  if (itm .eq. 0) then
306 c *** no tilt: copy A tensor directly
307  do j=1,5
308  fam(1,j)=faa(j)
309  end do
310  else
311 c *** transform A tensor into g axis system
312 c
313  call cd2km(d2km,alm,bem,gam)
314  do i=1,5
315  do j=1,5
316  fam(1,i)=fam(1,i)+d2km(1,i,j)*faa(j)
317  fam(2,i)=fam(2,i)+d2km(2,i,j)*faa(j)
318  end do
319  end do
320  end if
321 c
322 c----------------------------------------------------------------------
323 c Set diffusion tilt flag and g,A,W tensors in diffusion frame
324 c----------------------------------------------------------------------
325 c
326 c if (axialr .and. abs(gad).gt.RNDOFF) then
327 c gad=ZERO
328 cc ierr=NOALPHAD
329 c end if
330 c
331  if ( abs(ald).gt.rndoff .or.
332  # abs(gad).gt.rndoff .or.
333  # bed.gt.rndoff) itd=1
334 c
335 c *** no tilt: copy tensors directly
336  if (itd .eq. 0) then
337  do i=1,5
338  fgd(1,i)=fgm(i)
339  fwd(1,i)=fwm(i)
340  fad(1,i)=fam(1,i)
341  fad(2,i)=fam(2,i)
342  end do
343  else
344 c *** Transform A, g, W tensors into the diffusion frame
345 c
346  call cd2km(d2km,ald,bed,gad)
347  do i=1,5
348  do j=1,5
349  fgd(1,i)=fgd(1,i)+d2km(1,i,j)*fgm(j)
350  fgd(2,i)=fgd(2,i)+d2km(2,i,j)*fgm(j)
351  fwd(1,i)=fwd(1,i)+d2km(1,i,j)*fwm(j)
352  fwd(2,i)=fwd(2,i)+d2km(2,i,j)*fwm(j)
353  fad(1,i)=fad(1,i)+d2km(1,i,j)*fam(1,j)
354  # -d2km(2,i,j)*fam(2,j)
355  fad(2,i)=fad(2,i)+d2km(1,i,j)*fam(2,j)
356  # +d2km(2,i,j)*fam(1,j)
357  end do
358  end do
359  end if
360 c
361 c --- Allow antisymmetric K combinations if any tensor has
362 c imaginary elements (nonzero alm, gam, ald, or gad)
363 c
364  jkmn=1
365  do i=1,5
366  if ( abs(fgd(2,i)).ge.rndoff
367  # .or. abs(fad(2,i)).ge.rndoff
368  # .or. abs(fwd(2,i)).ge.rndoff ) jkmn=-1
369  end do
370 c
371 c --- Allow antisymmetric M combinations if there is a nonzero
372 c nuclear Zeeman interaction
373 c
374  if (abs(gamman).ge.rndoff .and. in2.gt.0) then
375  jmmn=-1
376  else
377  jmmn=1
378  end if
379 
380 c----------------------------------------------------------------------
381 c Check basis set parameters
382 c----------------------------------------------------------------------
383  inlemx=lemx
384  inlomx=lomx
385  inkmn=kmn
386  inkmx =kmx
387  inmmn=mmn
388  inmmx =mmx
389  inpnmx=ipnmx
390 c
391  if((lemx.gt.mxlval).and.(ipt.ne.0)) then
392  ierr=lemxhi
393  lemx=mxlval
394  endif
395 c
396  if(ipar(lemx).ne.1) lemx=lemx-1
397  if(ipar(lomx).ne.-1) lomx=lomx-1
398  if(lomx.gt.lemx) lomx=lemx-1
399  if(kmx.gt.lemx) kmx=lemx
400  if(mmx.gt.lemx) mmx=lemx
401  if(ipar(kmx).ne.1) kmx=kmx-1
402  if(ipnmx.gt.in2) ipnmx=in2
403  if((ipsi0.eq.0).and.(mmx.gt.ipnmx)) mmx=ipnmx
404  if(-kmn.gt.kmx) kmn=-kmx
405  if(-mmn.gt.mmx) mmn=-mmx
406  if(jkmn.eq.1) kmn=max(kmn,0)
407  if(jmmn.eq.1) mmn=max(kmn,0)
408 c
409  if(lemx.lt.0) lemx=0
410  if(lomx.lt.0) lomx=0
411  if(kmx.lt.0) kmx=0
412  if(mmx.lt.0) mmx=0
413  if(ipnmx.lt.0) ipnmx=0
414 c
415  if (inlemx .ne. lemx .or.
416  # inlomx .ne. lomx .or.
417  # inkmx .ne. kmx .or.
418  # inkmn .ne. kmn .or.
419  # inmmx .ne. mmx .or.
420  # inmmn .ne. mmn .or.
421  # inpnmx .ne. ipnmx ) then
422 c
423  ierr=bssadj
424  write (eprerr(bssadj)(24:),'(6(i3,'',''),i2)')
426  endif
427 c
428 c----------------------------------------------------------------------
429 c Determine basis set using rules in M,I,I,M & F
430 c (1) If there is no diffusion or magnetic tilt, only even K values
431 c are needed
432 c----------------------------------------------------------------------
433  if(itm.eq.0 .and. itd.eq.0) then
434  kdelta=2
435  else
436  kdelta=1
437  end if
438 c
439 c----------------------------------------------------------------------
440 c (2) In addition, if the magnetic and linewidth tensors are axial
441 c (in the diffusion frame) only even L values and no K values
442 c are needed
443 c----------------------------------------------------------------------
444  if (axiala .and. axialg .and. axialw .and.
445  # kdelta.eq.2 .and. kptmx.eq.0) then
446  ldelta=2
447  kmx=0
448  else
449  ldelta=1
450  end if
451 c
452 c----------------------------------------------------------------------
453 c check number of Lanczos/CG steps to execute
454 c----------------------------------------------------------------------
455 c
456  if (nstep.gt.mxstep) then
457  ierr=nstephi
458  write (eprerr(nstephi)(30:),'(i5)') mxstep
459  nstep=mxstep
460  endif
461 c
462 c----------------------------------------------------------------------
463 c Check specified tolerances and shifts for CG calculation
464 c----------------------------------------------------------------------
465 c
466  if (abs(cgtol) .lt. ten*rndoff) then
467  cgtol=ten*rndoff
468  ierr=cgtolhi
469  write (eprerr(cgtolhi)(30:),'(g9.3)') shiftr
470  end if
471  if (abs(shiftr) .lt. ten*rndoff) then
472  shiftr=10.0d0*rndoff
473  ierr=shifthi
474  write (eprerr(shifthi)(31:),'(g9.3)') shiftr
475  end if
476 c
477  return
478 c
479  end
480 
double precision, pointer, save pml
Definition: eprprm.f90:58
double precision, pointer, save ald
Definition: eprprm.f90:58
integer, pointer, save irflg
Definition: eprprm.f90:69
integer, save ipt
Definition: eprprm.f90:82
double precision, pointer, save alm
Definition: eprprm.f90:58
double precision, dimension(5), save faa
Definition: eprprm.f90:78
double precision, pointer, save ayy
Definition: eprprm.f90:58
subroutine lcheck(ierr)
Definition: lcheck.f90:40
double precision, pointer, save gxx
Definition: eprprm.f90:58
double precision, pointer, save dz
Definition: eprprm.f90:58
integer, parameter badjmp
Definition: errmsg.f90:51
double precision, pointer, save dx
Definition: eprprm.f90:58
double precision, dimension(5), save fgm
Definition: eprprm.f90:78
double precision, save w0
Definition: eprprm.f90:78
integer, pointer, save nfld
Definition: eprprm.f90:69
integer, pointer, save iaflg
Definition: eprprm.f90:69
subroutine cd2km(d2km, alpha, beta, gamma)
Definition: cd2km.f90:33
integer, pointer, save mxy
Definition: eprprm.f90:69
double precision, save expkzz
Definition: eprprm.f90:78
integer, parameter mxstep
Definition: nlsdim.f90:39
Definition: maxl.f90:22
double precision, pointer, save gad
Definition: eprprm.f90:58
double precision, save expl
Definition: eprprm.f90:78
integer, save lband
Definition: eprprm.f90:82
double precision, pointer, save c22
Definition: eprprm.f90:58
integer, pointer, save ist
Definition: eprprm.f90:69
double precision, dimension(5, 5), save cpot
Definition: eprprm.f90:78
double precision, pointer, save azz
Definition: eprprm.f90:58
double precision, pointer, save gyy
Definition: eprprm.f90:58
double precision, dimension(2, 5), save fad
Definition: eprprm.f90:78
double precision, pointer, save cgtol
Definition: eprprm.f90:58
integer, parameter shifthi
Definition: errmsg.f90:51
integer, parameter bssadj
Definition: errmsg.f90:51
integer, parameter iaxx
Definition: eprprm.f90:92
subroutine tocart(t, iflg)
Definition: tensym.f90:40
integer, parameter zerob0
Definition: errmsg.f90:51
integer, parameter baddj
Definition: errmsg.f90:51
integer, pointer, save mmx
Definition: eprprm.f90:69
integer, pointer, save ipdf
Definition: eprprm.f90:69
double precision, pointer, save bem
Definition: eprprm.f90:58
double precision, pointer, save dy
Definition: eprprm.f90:58
integer, pointer, save jkmn
Definition: eprprm.f90:69
integer, parameter zerog
Definition: errmsg.f90:51
integer, pointer, save nort
Definition: eprprm.f90:69
integer, parameter nstephi
Definition: errmsg.f90:51
integer, parameter posb0
Definition: errmsg.f90:51
double precision, pointer, save gamman
Definition: eprprm.f90:58
double precision, dimension(2, 5), save fgd
Definition: eprprm.f90:78
integer, pointer, save mzz
Definition: eprprm.f90:69
integer, parameter lemxhi
Definition: errmsg.f90:51
integer, parameter badav
Definition: errmsg.f90:51
integer, save lptmx
Definition: eprprm.f90:82
integer, pointer, save mmn
Definition: eprprm.f90:69
integer, pointer, save ipnmx
Definition: eprprm.f90:69
integer, parameter badhifld
Definition: errmsg.f90:51
integer, pointer, save igflg
Definition: eprprm.f90:69
double precision, dimension(nfprm), target, save fepr
Definition: eprprm.f90:55
integer, parameter zeroin2
Definition: errmsg.f90:51
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
integer, parameter zeronfld
Definition: errmsg.f90:51
double precision, pointer, save gam
Definition: eprprm.f90:58
double precision, pointer, save c20
Definition: eprprm.f90:58
integer, parameter idx
Definition: eprprm.f90:92
double precision, dimension(2, 5), save fam
Definition: eprprm.f90:78
double precision, pointer, save pmxy
Definition: eprprm.f90:58
integer, save kptmx
Definition: eprprm.f90:82
integer, save kband
Definition: eprprm.f90:82
double precision, save expkxy
Definition: eprprm.f90:78
double precision, pointer, save djfprp
Definition: eprprm.f90:58
integer, pointer, save nstep
Definition: eprprm.f90:69
double precision, pointer, save gzz
Definition: eprprm.f90:58
integer, parameter iwxx
Definition: eprprm.f90:92
integer, parameter igxx
Definition: eprprm.f90:92
integer, save ldelta
Definition: eprprm.f90:82
integer, pointer, save lomx
Definition: eprprm.f90:69
double precision, pointer, save bed
Definition: eprprm.f90:58
integer, pointer, save kmn
Definition: eprprm.f90:69
integer, parameter cgtolhi
Definition: errmsg.f90:51
double precision, pointer, save pmzz
Definition: eprprm.f90:58
double precision, pointer, save axx
Definition: eprprm.f90:58
integer, pointer, save iwflg
Definition: eprprm.f90:69
double precision, dimension(5), save fwm
Definition: eprprm.f90:78
integer, pointer, save in2
Definition: eprprm.f90:69
double precision, pointer, save c42
Definition: eprprm.f90:58
integer, pointer, save kmx
Definition: eprprm.f90:69
double precision, save g0
Definition: eprprm.f90:78
integer, parameter mxlval
Definition: maxl.f90:25
integer, pointer, save jmmn
Definition: eprprm.f90:69
integer, pointer, save lemx
Definition: eprprm.f90:69
integer, pointer, save ml
Definition: eprprm.f90:69
character *50, dimension(neperr), save eprerr
Definition: errmsg.f90:59
double precision, pointer, save b0
Definition: eprprm.f90:58
integer, save itm
Definition: eprprm.f90:82
double precision, pointer, save wyy
Definition: eprprm.f90:58
integer, save ipsi0
Definition: eprprm.f90:82
double precision, pointer, save shiftr
Definition: eprprm.f90:58
double precision, dimension(2, 5), save fwd
Definition: eprprm.f90:78
double precision, pointer, save oss
Definition: eprprm.f90:58
double precision, pointer, save wzz
Definition: eprprm.f90:58
double precision, pointer, save wxx
Definition: eprprm.f90:58
double precision, save a0
Definition: eprprm.f90:78
double precision, parameter rndoff
Definition: rnddbl.f90:86
integer, save kdelta
Definition: eprprm.f90:82
integer, save itd
Definition: eprprm.f90:82
double precision, pointer, save djf
Definition: eprprm.f90:58