NLSL
sitec.f90
Go to the documentation of this file.
1 c NLSL Version 1.5 beta 11/24/95
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine SITEC
5 c =========================
6 c
7 c sites <n>
8 c
9 c n : number of different "sites" or components in each spectrum
10 c
11 c----------------------------------------------------------------------
12  subroutine sitec( line )
13 c
14  use nlsdim
15  use expdat
16  use parcom
17  use tridag
18  use basis
19  use stdio
20 c
21  implicit none
22  character*80 line
23 c
24  character token*30
25 c
26  integer i,ixpar,ixsite,j,lth,ival
27 c
28  integer itrim
29  logical itoken
30  external itrim,itoken
31 c
32 c----------------------------------------------------------------------
33 c Get the number of sites required
34 c----------------------------------------------------------------------
35 c
36  call gettkn(line,token,lth)
37 c
38  if (lth.gt.0 .and. itoken(token,lth,ival)) then
39 
40  if (ival.le.mxsite.and.ival.ge.1) then
41 c
42 c ------------------------------------------------------------
43 c If the number of sites is being reduced, mark all basis
44 c sets and tridiagonal matrices belonging to the sites that
45 c are being removed
46 c ------------------------------------------------------------
47  if (ival.lt.nsite) then
48  do i=ival+1,nsite
49  do j=1,mxspc
50  bsused( basno(i,j) )=0
51  modtd(i,j)=1
52  end do
53  end do
54  end if
55 c
56 c --------------------------------------------------
57 c Remove any variable parameter belonging to a site
58 c that is no longer defined
59 c --------------------------------------------------
60  do i=1,nprm
61  if (ixst(i).gt.ival) then
62  ixpar=ixpr(i)
63  ixsite=ixst(i)
64  token=tag(i)
65  call rmvprm(ixpar,ixsite,token)
66  end if
67  end do
68 c
69  nsite=ival
70 c
71 c
72 c ------------------------------------------------------------
73 c Issue warning message if some of the data in a multiple-site
74 c multiple spectrum fit are not normalized
75 c ------------------------------------------------------------
76  if (nsite.gt.1 .and. nser.gt.1) then
77  do i=1,nser
78  if (nrmlz(i).eq.0)
79  # write (luttyo,1002) i,dataid(i)(:itrim(dataid(i)))
80  end do
81  end if
82 c
83 c *** Bad # of sites
84  else
85  write(luttyo,1000) mxsite
86  end if
87 c
88 c *** Expecting integer value
89  else
90  write (luttyo,1001)
91  end if
92  return
93 c
94 c ###### format statements ########################################
95 c
96  1000 format('*** Number of sites must be ',i2,' or less ***')
97  1001 format('*** Expecting integer value ***')
98  1002 format('*** WARNING: spectrum ',i2,' (',a,
99  # ') not normalized ***')
100  end
integer, dimension(mxvar), save ixst
Definition: parcom.f90:62
integer, dimension(mxspc), save nrmlz
Definition: expdat.f90:45
integer, dimension(mxtdm), save bsused
Definition: basis.f90:23
integer, dimension(mxsite, mxspc), save basno
Definition: basis.f90:23
integer, dimension(mxsite, mxspc), save modtd
Definition: tridag.f90:28
Definition: stdio.f90:26
character *30, dimension(mxspc), save dataid
Definition: expdat.f90:51
Definition: basis.f90:19
subroutine gettkn(line, token, lth)
Written for free-form input of parameters for slow-motional calculations. Returns a token consisting ...
Definition: strutl1.f90:75
integer, parameter mxspc
Definition: nlsdim.f90:39
integer, dimension(mxvar), save ixpr
Definition: parcom.f90:62
integer, save nprm
Definition: parcom.f90:62
integer, save nsite
Definition: parcom.f90:62
integer, parameter mxsite
Definition: nlsdim.f90:39
integer, save nser
Definition: parcom.f90:62
subroutine sitec(line)
Definition: sitec.f90:13
integer, parameter luttyo
Definition: stdio.f90:29
subroutine rmvprm(ix, ix2, ident)
Definition: rmvprm.f90:25
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67