NLSL
setprm.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 beta 11/23/95
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine SETPRM
5 c =========================
6 c>@brief This file contains two routines that set a given parameter,
7 c> specified by an index into the fparm or iparm array, given the
8 c> parameter value and a site/spectrum index. The secondary index
9 c> is interpreted according to whether or not the parameter belongs to
10 c> a set of parameters that must be the same for all components of a
11 c> given spectrum. This set of "spectral parameters" includes
12 c>
13 c> PHASE (spectral phase angle)
14 c> PSI (director tilt angle)
15 c> B0 (spectrometer field)
16 c> LB (spectral line broadening)
17 c> RANGE (spectral range)
18 c>
19 c> and is identified by the logical function spcpar(index) coded
20 c> in this file. These parameters are kept arrays separate from the
21 c> fparm arrays in which the site parameters are kept.
22 c>
23 c> The index ixsite is interpreted as referring to a spectrum for
24 c> spectral parameters; otherwise, the index refers to a component
25 c> within a given spectrum.
26 c>
27 c> Whenever a parameter is changed, if recalculation of a spectrum
28 c> with the new parameter value would require recalculation of the
29 c> tridiagonal matrix, a flag is set for the given sites/spectra
30 c> affected by that parameter.
31 c
32 c----------------------------------------------------------------------
33  subroutine setprm(ixparm,ixsite,fval)
34 c
35  use nlsdim
36  use eprprm
37  use expdat
38  use tridag
39  use parcom
40  use symdef
41  use stdio
42 c
43  implicit none
44  integer ixparm,ixsite
45  double precision fval
46 c
47  integer i,ixmax,jx,jx1,jx2
48 c
49  logical spcpar,spcprm
50  external spcpar
51 c
52  spcprm=spcpar(ixparm)
53  if (spcprm) then
54  ixmax=mxspc
55  else
56  ixmax=mxsite
57  endif
58 c
59 c --- Index between 0 and maximum: single site or spectrum
60 c
61  if (ixsite.gt.0 .and. ixsite.le.ixmax) then
62  jx1=ixsite
63  jx2=ixsite
64 c
65 c --- Zero index: all sites or spectra
66 c
67  else if (ixsite.eq.0) then
68  jx1=1
69  jx2=ixmax
70 c
71 c --- Illegal index
72 c
73  else
74  write(luout,1000) jx
75  if (luout.ne.luttyo) write (luttyo,1000) jx
76  end if
77 c
78  do jx=jx1,jx2
79 c
80 c --- Spectral parameters are stored separately:
81 c
82  if (spcprm) then
83  if (ixparm.eq.ib0) then
84  sb0(jx)=fval
85  else if (ixparm.eq.iphase) then
86  sphs(jx)=fval
87  else if (ixparm.eq.ipsi) then
88  spsi(jx)=fval
89  else if (ixparm.eq.ilb) then
90  slb(jx)=fval
91  else if (ixparm.eq.irange) then
92  srng(jx)=fval
93  sbi(jx)=sb0(jx)-0.5d0*fval
94  if (npts(jx).gt.1) sdb(jx)=fval/float(npts(jx)-1)
95  end if
96 c
97 c Mark tridiagonal matrix for recalculation only if tilt angle
98 c has been changed (spectra can be calculated from an existing
99 c tridiagonal matrix if the other parameters above are changed)
100 c
101  if (ixparm.eq.ipsi) then
102  do i=1,nsite
103  modtd(i,jx)=1
104  end do
105  end if
106 c
107 c --- Site parameter: mark tridiagonal matrix for recalculation
108 c (except for Gaussian inhomogeneous broadening
109 c and isotropic Lorentzian broadening parameters)
110 c
111  else
112  fparm(abs(mod(ixparm,100)),jx)=fval
113  if (ixparm.ne.igib0 .and. ixparm.ne.igib2 .and. .not.
114  # (ixparm.eq.iwxx .and.iparm(iiwflg,ixsite).eq.spherical) )
115  # then
116  do i=1,nser
117  modtd(jx,i)=1
118  end do
119  end if
120 c
121  end if
122 c
123  end do
124 c
125  return
126  1000 format(' *** Illegal index: ',i2,' ***')
127  end
128 
129 
130 c----------------------------------------------------------------------
131 c =========================
132 c subroutine SETIPR
133 c =========================
134 c
135 c>@brief Analogous routine to setprm for integer parameters
136 c> There are only two user-settable integer spectrum parameters: nfield,
137 c> and ideriv. These are normally determined by the input data file, but
138 c> are needed when calculations are to performed without data or fitting.
139 c
140 c----------------------------------------------------------------------
141  subroutine setipr(ixparm,ixsite,ival)
142 c
143  use nlsdim
144  use eprprm
145  use expdat
146  use parcom
147  use tridag
148  use basis
149  use stdio
150 c
151  implicit none
152  integer ixparm,ixsite,ival
153 c
154  integer i,ixmax,ixp,j,jx,jx1,jx2
155 c
156  integer itrim
157  logical spcpar,spcprm
158  external itrim,spcpar
159 c
160  ixp=abs(mod(ixparm,100))
161  spcprm=spcpar(ixparm)
162 c
163  if (spcprm) then
164  ixmax=mxspc
165  else
166  ixmax=mxsite
167  endif
168 c
169 c Single site/spectrum
170 c
171  if (ixsite.gt.0 .and. ixsite.le.mxsite) then
172  jx1=ixsite
173  jx2=ixsite
174 c
175 c All sites/spectra
176 c
177  else if (ixsite.eq.0) then
178  jx1=1
179  jx2=ixmax
180 c
181 c Illegal index
182 c
183  else
184  write(luout,1000) jx
185  if (luout.ne.luttyo) write (luttyo,1000) jx
186  end if
187 c
188  do jx=jx1,jx2
189 c
190 c Spectrum parameters
191 c
192  if (spcprm) then
193  if (ixp.eq.infld) then
194 c
195 c Cannot change number of points in an existing datafile
196 c
197  if (jx.le.nspc) then
198  write(luout,1001) dataid(jx)(:itrim(dataid(jx)))
199  if (luout.ne.luttyo) write (luttyo,1001)
200  # dataid(jx)(:itrim(dataid(jx)))
201  else
202  npts(jx)=ival
203  if (ival.gt.1) sdb(jx)=srng(jx)/float(ival-1)
204 c
205  nft(jx)=1
206  9 nft(jx)=nft(jx)*2
207  if (nft(jx).lt.npts(jx)) go to 9
208  end if
209 c
210  else if (ixparm.eq.iiderv) then
211  idrv(jx)=ival
212  end if
213 c
214 c Site parameters
215 c Note that changes in site-related integer quantities always
216 c require recalculation of the tridiagonal matrix
217 c
218  else
219  iparm(ixp,jx)=ival
220  do i=1,mxspc
221  modtd(jx,i)=1
222  end do
223  end if
224  end do
225 c
226  return
227  1000 format('*** Illegal index: ',i2,' ***')
228  1001 format('*** Number of data points for file ',a,
229  # ' cannot be changed ***')
230  end
231 
232 
233 
234 c----------------------------------------------------------------------
235 c =========================
236 c function SPCPAR
237 c =========================
238 c
239 c>@brief Returns .true. if the index argument corresponds to a floating
240 c> point or integer parameter that cannot change for an individual spectrum,
241 c> regardless of the number of sites for which the calculation is made.
242 c> These include:
243 c>
244 c> PHASE (spectral phase angle)
245 c> PSI (director tilt angle)
246 c> LB (spectral line broadening)
247 c> B0 (spectrometer field)
248 c> FIELDI (initial field of spectrum)
249 c> DFLD (field step per point in spectrum)
250 c> RANGE (field range of spectrum)
251 c>
252 c> NFIELD (number of points in spectrum)
253 c> IDERIV (0th/1st derivative flag for spectrum)
254 c
255 c
256 c Includes:
257 c nlsdim.inc
258 c eprprm.inc
259 c parcom.inc
260 c
261 c----------------------------------------------------------------------
262  function spcpar( ix )
263 c
264  use nlsdim
265  use eprprm
266  use parcom
267 c
268  implicit none
269  integer ityp,ix,ixp
270  logical spcpar
271 c
272  ixp = mod(ix,100)
273  ityp = ix/100
274  spcpar = (ityp.eq.0 .and.
275  # (ixp.eq.iphase .or. ixp.eq.ipsi .or. ixp.eq.ilb
276  # .or. ixp.eq.ib0 .or. ixp.eq.ifldi .or. ixp.eq.idfld
277  # .or. ixp.eq.irange) )
278  # .or.(ityp.eq.1 .and.
279  # (ixp.eq.infld .or. ixp.eq.iiderv) )
280 c
281  return
282  end
283 
284  function ptype( ix )
285  implicit none
286  integer ix
287  character*7 ptype
288 c
289  logical spcpar
290  external spcpar
291  if (spcpar(ix)) then
292  ptype='spectra'
293  else
294  ptype='sites'
295  end if
296  return
297  end
298 
299  function ixlim( ix )
300 c
301  use nlsdim
302  use expdat
303  use parcom
304 c
305  implicit none
306  integer ixlim, ix
307  logical spcpar
308  external spcpar
309 c
310  if (spcpar(ix)) then
311  ixlim=nser
312  else
313  ixlim=nsite
314  end if
315  return
316  end
317 
318 c
319 c------------------------------------------------------------------------
320 c =========================
321 c function GETPRM
322 c =========================
323 c>@brief Given a parameter index and a site/spectrum index, this function
324 c> returns the value of the parameter from the fparm array (for
325 c> site parameters) or from the spectral parameter arrays for spectrum
326 c> parameters.
327 c------------------------------------------------------------------------
328 c
329 
330  function getprm(ixparm,ixsite)
331 c
332  use nlsdim
333  use eprprm
334  use expdat
335  use tridag
336  use parcom
337  use symdef
338  use stdio
339 c
340  implicit none
341  integer ixparm,ixsite
342  double precision getprm
343 c
344  integer ixmax
345 c
346  logical spcpar,spcprm
347  external spcpar
348 c
349  spcprm=spcpar(ixparm)
350  if (spcprm) then
351  ixmax=nser
352  else
353  ixmax=nsite
354  endif
355 c
356 c --- Illegal index
357 c
358  if (ixsite.le.0.or.ixsite.gt.ixmax) then
359  write(luttyo,1000) ixsite
360  return
361  end if
362 c
363  if (spcprm) then
364  if (ixparm.eq.ib0) then
365  getprm=sb0(ixsite)
366  else if (ixparm.eq.iphase) then
367  getprm=sphs(ixsite)
368  else if (ixparm.eq.ipsi) then
369  getprm=spsi(ixsite)
370  else if (ixparm.eq.ilb) then
371  getprm=slb(ixsite)
372  else if (ixparm.eq.irange) then
373  getprm=srng(ixsite)
374  end if
375 c
376  else
377  getprm=fparm(abs(mod(ixparm,100)),ixsite)
378  end if
379 c
380  return
381 c
382  1000 format(' *** Illegal index: ',i2,' ***')
383  end
integer, save luout
Definition: stdio.f90:32
double precision, dimension(mxspc), save sbi
Definition: expdat.f90:40
integer, parameter irange
Definition: eprprm.f90:92
integer function ixlim(ix)
Definition: setprm.f90:300
double precision, dimension(mxspc), save sdb
Definition: expdat.f90:40
integer, parameter igib2
Definition: eprprm.f90:92
integer, parameter ilb
Definition: eprprm.f90:92
integer, dimension(mxspc), save nft
Definition: expdat.f90:45
subroutine setipr(ixparm, ixsite, ival)
Analogous routine to setprm for integer parameters There are only two user-settable integer spectrum ...
Definition: setprm.f90:142
double precision, dimension(mxspc), save srng
Definition: expdat.f90:40
integer, parameter iiwflg
Definition: eprprm.f90:101
integer, dimension(mxsite, mxspc), save modtd
Definition: tridag.f90:28
Definition: stdio.f90:26
character *30, dimension(mxspc), save dataid
Definition: expdat.f90:51
double precision, dimension(mxspc), save spsi
Definition: expdat.f90:40
double precision, dimension(nfprm, mxsite), target, save fparm
Definition: parcom.f90:54
Definition: basis.f90:19
logical function spcpar(ix)
Returns .true. if the index argument corresponds to a floating point or integer parameter that cannot...
Definition: setprm.f90:263
integer, parameter ipsi
Definition: eprprm.f90:92
integer, parameter igib0
Definition: eprprm.f90:92
integer, parameter mxspc
Definition: nlsdim.f90:39
integer, save nspc
Definition: expdat.f90:45
integer, save nsite
Definition: parcom.f90:62
double precision, dimension(mxspc), save sphs
Definition: expdat.f90:40
integer, dimension(niprm, mxsite), target, save iparm
Definition: parcom.f90:60
integer, parameter mxsite
Definition: nlsdim.f90:39
integer, save nser
Definition: parcom.f90:62
integer, parameter ib0
Definition: eprprm.f90:92
integer, parameter iwxx
Definition: eprprm.f90:92
integer, parameter iphase
Definition: eprprm.f90:92
integer, dimension(mxspc), save idrv
Definition: expdat.f90:45
integer, parameter iiderv
Definition: eprprm.f90:101
character *7 function ptype(ix)
Definition: setprm.f90:285
integer, parameter luttyo
Definition: stdio.f90:29
double precision, dimension(mxspc), save sb0
Definition: expdat.f90:40
integer, parameter idfld
Definition: eprprm.f90:92
integer, dimension(mxspc), save npts
Definition: expdat.f90:45
integer, parameter spherical
Definition: symdef.f90:14
double precision function getprm(ixparm, ixsite)
Given a parameter index and a site/spectrum index, this function returns the value of the parameter f...
Definition: setprm.f90:331
double precision, dimension(mxspc), save slb
Definition: expdat.f90:40
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
integer, parameter infld
Definition: eprprm.f90:101
integer, parameter ifldi
Definition: eprprm.f90:92