NLSL
setmts.f90
Go to the documentation of this file.
1 c Version 1.5.1 beta 3/21/96
2 c----------------------------------------------------------------------
3 c =========================
4 c function SETMTS
5 c =========================
6 c
7 c Given a set of slow-motional EPR parameters, check the MTS indices
8 c and adjust them accordingly.
9 c
10 c Inputs
11 c fparmi Array of floating-point EPR parameters appearing in the
12 c order defined in /eprprm/
13 c iparmi Array of integer EPR parameters appearing in the
14 c order defined in /eprprm/
15 c
16 c Output
17 c mts Array to receive the corrected MTS parameters
18 c Information is packed as documented in mtsdef.inc
19 c
20 c Returns 0 for successful completion
21 c 8 if lemx exceeded limit for calculations w/potential
22 c 9 if mts parameters were changed from those in iparmi array
23 c
24 c----------------------------------------------------------------------
25  function setmts(fparmi,iparmi,mts)
26 c
27  use nlsdim
28  use eprprm
29  use errmsg
30  use mtsdef
31  use maxl
32  use rnddbl
33 c
34  implicit none
35  double precision fparmi(nfprm)
36  integer setmts,iparmi(niprm),mts(mxmts)
37 c
38  integer i
39  logical knon0p,changed
40 c
41  integer ipar
42  logical isaxial
43  external ipar,isaxial
44 c
45 c----------------------------------------------------------------------
46 c If there is no associated list of basis set indices,
47 c check basis set parameters
48 c----------------------------------------------------------------------
49 c
50  setmts=0
51  knon0p = abs(fparmi(ic20+1)).gt.rndoff .or.
52  # abs(fparmi(ic20+2)).gt.rndoff .or.
53  # abs(fparmi(ic20+3)).gt.rndoff .or.
54  # abs(fparmi(ic20+4)).gt.rndoff
55 c
56  do i=1,ntrc
57  mts(i)=iparmi(ilemx+i-1)
58  end do
59 c
60  mts(nipsi0)=0
61  if ( (abs(fparmi(ipsi)).gt.rndoff .and.
62  # abs(fparmi(ipsi)-1.8d2).gt.rndoff)
63  # .or. iparmi(inort).gt.1) mts(nipsi0)=1
64  mts(nin2)=iparmi(iin2)
65 c
66 c --- Allow antisymmetric K combinations if any tensor has
67 c imaginary elements (nonzero alm, gam, ald, or gad)
68 c
69  mts(njkmn)=1
70  if (abs(fparmi(ialm)).ge.rndoff .or.
71  # abs(fparmi(igam)).ge.rndoff .or.
72  # abs(fparmi(iald)).ge.rndoff .or.
73  # abs(fparmi(igad)).ge.rndoff) mts(njkmn)=-1
74 c
75 c --- Allow antisymmetric M combinations if there is a nonzero
76 c nuclear Zeeman interaction
77 c
78  if (abs(fparmi(igaman)).ge.rndoff .and.
79  # iparmi(iin2).gt.0) then
80  mts(njmmn)=-1
81  else
82  mts(njmmn)=1
83  end if
84 
85 c
86  if((iparmi(ilemx).gt.mxlval).and.(ipt.ne.0)) then
87  setmts=lemxhi
88  iparmi(ilemx)=mxlval
89  endif
90 c
91  if (mts(nlemx).gt.mxlval) mts(nlemx)=mxlval
92  if (ipar(mts(nlemx)).ne.1) mts(nlemx)=mts(nlemx)-1
93  if (ipar(mts(nlomx)).ne.-1) mts(nlomx)=mts(nlomx)-1
94  if (mts(nlomx).gt.mts(nlemx)) mts(nlomx)=mts(nlemx)-1
95  if (mts(nkmx).gt.mts(nlemx)) mts(nkmx)=mts(nlemx)
96  if (mts(nmmx).gt.mts(nlemx)) mts(nmmx)=mts(nlemx)
97  if (ipar(mts(nkmx)).ne.1) mts(nkmx)=mts(nkmx)-1
98  if (mts(nipnmx).gt.in2) mts(nipnmx)=iparmi(iin2)
99  if (mts(nipsi0).eq.0.and.mts(nmmx).gt.mts(nipnmx) )
100  # mts(nmmx)=mts(nipnmx)
101  if (mts(njkmn).eq.1) mts(nkmn)=max(mts(nkmn),0)
102  if (mts(njmmn).eq.1) mts(nmmn)=max(mts(nmmn),0)
103 c
104  if(mts(nlemx).lt.0) mts(nlemx)=0
105  if(mts(nlomx).lt.0) mts(nlomx)=0
106  if(mts(nkmx).lt.0) mts(nkmx)=0
107  if(mts(nmmx).lt.0) mts(nmmx)=0
108  if(mts(nipnmx).lt.0) mts(nipnmx)=0
109 c
110  changed=.false.
111  do i=1,ntrc
112  changed = changed .or.(mts(i).ne.iparmi(ilemx+i-1))
113  end do
114 c
115  if (changed) then
116  setmts=bssadj
117  write (eprerr(bssadj)(24:),'(6(i3,'',''),i2)')
119  end if
120 c
121 c----------------------------------------------------------------------
122 c Determine basis set using rules in M,I,I,M & F
123 c (1) If there is no diffusion tilt, only even K values are needed
124 c----------------------------------------------------------------------
125  if(abs(ald).gt.rndoff .or. abs(bed).gt.rndoff .or.
126  # (abs(bed)-1.8d2).gt.rndoff .or.abs(gad).gt.rndoff) then
127  mts(nkdel)=1
128  else
129  mts(nkdel)=2
130  end if
131 c
132 c----------------------------------------------------------------------
133 c (2) In addition, if the magnetic and linewidth tensors are axial,
134 c and there is no diffusion tilt or K.ne.0 potential terms,
135 c only even L values and no K values are needed
136 c----------------------------------------------------------------------
137  if( isaxial(fparmi(igxx),iparmi(iigflg)) .and.
138  # isaxial(fparmi(iaxx),iparmi(iiaflg)) .and.
139  # isaxial(fparmi(iwxx),iparmi(iiwflg)) .and.
140  # (mts(nkdel).eq.2).and.(.not.knon0p)) then
141  mts(nldel)=2
142  mts(nkmx)=0
143  else
144  mts(nldel)=1
145  end if
146 c
147  return
148  end
double precision, pointer, save ald
Definition: eprprm.f90:58
integer, parameter nkmn
Definition: mtsdef.f90:29
integer, save ipt
Definition: eprprm.f90:82
integer, parameter igad
Definition: eprprm.f90:92
integer, parameter njkmn
Definition: mtsdef.f90:29
integer, parameter igam
Definition: eprprm.f90:92
integer, parameter iiaflg
Definition: eprprm.f90:101
integer, parameter iald
Definition: eprprm.f90:92
integer, parameter niprm
Definition: nlsdim.f90:57
Definition: maxl.f90:22
double precision, pointer, save gad
Definition: eprprm.f90:58
integer, parameter iiwflg
Definition: eprprm.f90:101
integer, parameter ntrc
Definition: mtsdef.f90:29
integer, parameter bssadj
Definition: errmsg.f90:51
integer, parameter iaxx
Definition: eprprm.f90:92
integer, parameter mxmts
Definition: nlsdim.f90:39
integer, pointer, save mmx
Definition: eprprm.f90:69
integer function setmts(fparmi, iparmi, mts)
Definition: setmts.f90:26
integer, parameter ipsi
Definition: eprprm.f90:92
integer, parameter nmmx
Definition: mtsdef.f90:29
integer, parameter nlemx
Definition: mtsdef.f90:29
integer, parameter nipsi0
Definition: mtsdef.f90:29
integer, parameter lemxhi
Definition: errmsg.f90:51
integer, pointer, save mmn
Definition: eprprm.f90:69
integer, pointer, save ipnmx
Definition: eprprm.f90:69
integer, parameter nkdel
Definition: mtsdef.f90:29
integer, parameter nin2
Definition: mtsdef.f90:29
integer, parameter iwxx
Definition: eprprm.f90:92
integer, parameter igaman
Definition: eprprm.f90:92
integer, parameter igxx
Definition: eprprm.f90:92
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 ic20
Definition: eprprm.f90:92
integer, pointer, save in2
Definition: eprprm.f90:69
integer, parameter nipnmx
Definition: mtsdef.f90:29
integer, parameter nmmn
Definition: mtsdef.f90:29
integer, pointer, save kmx
Definition: eprprm.f90:69
integer, parameter ilemx
Definition: eprprm.f90:101
integer, parameter inort
Definition: eprprm.f90:101
integer, parameter iigflg
Definition: eprprm.f90:101
integer, parameter mxlval
Definition: maxl.f90:25
integer, pointer, save lemx
Definition: eprprm.f90:69
integer, parameter njmmn
Definition: mtsdef.f90:29
character *50, dimension(neperr), save eprerr
Definition: errmsg.f90:59
integer, parameter nkmx
Definition: mtsdef.f90:29
integer, parameter nlomx
Definition: mtsdef.f90:29
double precision, parameter rndoff
Definition: rnddbl.f90:86
integer, parameter nfprm
Definition: nlsdim.f90:57
integer, parameter ialm
Definition: eprprm.f90:92
integer, parameter iin2
Definition: eprprm.f90:101
integer, parameter nldel
Definition: mtsdef.f90:29