NLSL
addprm.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 ADDPRM
5 c =========================
6 c
7 c>@brief Add a parameter to list of parameters being varied for nonlinear
8 c> least-squares. Also maintain the lists of
9 c>
10 c> 1. ixst : secondary parameter index (for multiple sites/spectra)
11 c> 2. ibnd : boundary flag for variable parameter
12 c> 3. prmin : minimum for variable parameter
13 c> 4. prmax : maximum for variable parameter
14 c> 5. prscl : desired accuracy for given parameter
15 c> 6. xfdstp : Size of forward-differences step
16 c> 7. tag : Name of the parameter (character*9)
17 c> 8. ixx : index of each element in the fparm array into the x vector
18 c>
19 c> The values of these parameters for the given variable parameter
20 c> are passed to ADDPRM:
21 c> subroutine parameter |kept in list
22 c> ---------------------|------------
23 c> ix2 |ixst
24 c> ibd |ibnd
25 c> prmn |prmin
26 c> prmx |prmax
27 c> prsc |prscl
28 c> step |xfdstp
29 c> ident |tag
30 c>
31 c> Notes:
32 c>
33 c> - Each addition to the list is included in sort order with respect to
34 c> the index ix. This is an artifact of the original implementation
35 c> of the program, and no longer required.
36 c>
37 c> - If ix2 = -1, the parameter specified by ix is to be varied individually
38 c> for EACH existing site/spectrum
39 c> If ix2 = 0, the parameter is to be varied for ALL sites/spectra at once
40 c
41 c Includes:
42 c nlsdim.inc
43 c eprprm.inc
44 c expdat.inc
45 c parcom.inc
46 c lpnam.inc
47 c lmcom.inc
48 c stdio.inc
49 c rndoff.inc
50 c prmeqv.inc
51 c
52 c----------------------------------------------------------------------
53  subroutine addprm(ix,ix2,ibd,prmn,prmx,prsc,step,ident)
54 c
55  use nlsdim
56  use eprprm
57  use expdat
58  use parcom
59  use lpnam
60  use lmcom
61  use stdio
62  use rnddbl
63 c use prmeqv
64 c
65  implicit none
66  integer ix,ix2,ibd
67  double precision prmn,prmx,prsc,step
68  character*9 ident
69 c
70  integer i,isi,ixabs,j,jx,jx1,jx2,k,k1,l,lu,nadd,nmov
71 c
72  integer itrim
73  logical spcpar,tcheck
74  external itrim,spcpar,tcheck
75 c
76 c######################################################################
77 c
78 c
79 c --- Check index bounds
80 c
81  if ( (spcpar(ix) .and. ix2.gt.mxspc)
82  # .or. (.not.spcpar(ix) .and. ix2.gt.mxsite) ) then
83  write (luttyo,1002)
84  return
85  end if
86 c
87 c----------------------------------------------------------------------
88 c ix2= -1: vary the parameter individually for EACH exiting site/spectrum
89 c ix2= 0: vary the same parameter for ALL existing sites/spectra
90 c----------------------------------------------------------------------
91  if (ix2.lt.0) then
92  jx1=1
93  if (spcpar(ix)) then
94  jx2=nser
95  else
96  jx2=nsite
97  end if
98  else
99  jx1=ix2
100  jx2=ix2
101  end if
102 c
103  nadd=jx2-jx1+1
104 c
105 c----------------------------------------------------------------------
106 c See if there is enough room in the parameter list
107 c----------------------------------------------------------------------
108  if (nprm+nadd.gt.mxvar) then
109  write(luttyo,1001) mxvar
110  return
111  end if
112 c
113 c----------------------------------------------------------------------
114 c Loop over all parameters to be added
115 c----------------------------------------------------------------------
116  lu=luttyo
117  do 21 jx=jx1,jx2
118 c
119  if (.not.tcheck(ix,jx,ident,lu)) return
120  lu=0
121 c
122 c----------------------------------------------------------------------
123 c Search parameter list to see if parameter is already there
124 c----------------------------------------------------------------------
125  ixabs=abs(mod(ix,100))
126  do i=1,nprm
127 c *** parameter already in list
128  if (ixpr(i).eq.ixabs .and. ixst(i).eq.jx) then
129  write(luttyo,1000) tag(i)
130  return
131  end if
132 c *** parameter not in list
133  if ( ixpr(i).gt.ixabs .or.
134  # (ixpr(i).eq.ixabs.and.ixst(i).gt.jx) ) go to 14
135  end do
136 c
137 c----------------------------------------------------------------------
138 c Reached end of list: parameter will be added to the end
139 c----------------------------------------------------------------------
140  i=nprm+1
141  nprm=nprm+1
142  go to 18
143 c
144 c----------------------------------------------------------------------
145 c Found proper location for new parameter: move top of list up
146 c----------------------------------------------------------------------
147  14 nmov=nprm-i+1
148  do j=1,nmov
149  k=nprm-j+1
150  k1=k+1
151  x(k1)=x(k)
152  ixpr(k1)=ixpr(k)
153  ixst(k1)=ixst(k)
154  prmin(k1)=prmin(k)
155  prmax(k1)=prmax(k)
156  prscl(k1)=prscl(k)
157  xfdstp(k1)=xfdstp(k)
158  ibnd(k1)=ibnd(k)
159  tag(k1)=tag(k)
160  if (ixst(k1).le.0) then
161  do l=1,mxsite
162  ixx(ixpr(k1),l)=k1
163  end do
164  else
165  ixx(ixpr(k1),ixst(k1))=k1
166  end if
167 c
168  end do
169  nprm=nprm+1
170 c
171 c----------------------------------------------------------------------
172 c Add new parameter(s) to proper location in list
173 c----------------------------------------------------------------------
174  18 if (jx.eq.0) then
175  x(i)=fparm(ixabs,1)
176  do l=1,mxsite
177  ixx(ixabs,l)=i
178  end do
179  else
180  x(i)=fparm(ixabs,jx)
181  ixx(ixabs,jx)=i
182  end if
183  ixpr(i)=ixabs
184  ixst(i)=jx
185  prmin(i)=prmn
186  prmax(i)=prmx
187  prscl(i)=prsc
188  xfdstp(i)=step*x(i)
189  if(dabs(xfdstp(i)).lt.rndoff) xfdstp(i)=step
190  ibnd(i)=ibd
191  tag(i)=ident
192 c
193  if (jx.ne.0) then
194  k=itrim(tag(i))
195  write(tag(i)(k+1:),2000) jx
196  end if
197 c
198  21 continue
199 c
200  write(luttyo,1003) nprm
201  return
202 c
203 c #### format statements ###########################################
204 c
205  1000 format('*** Parameter ''',a,''' is already being varied ***')
206  1001 format('*** Limit of',i2,' variable parameters exceeded ***')
207  1002 format('*** Index out of range ***')
208  1003 format(i3,' parameters now being varied')
209  2000 format('(',i1,')')
210  end
integer, dimension(mxvar), save ixst
Definition: parcom.f90:62
integer, parameter mxvar
Definition: nlsdim.f90:39
Definition: stdio.f90:26
double precision, dimension(nfprm, mxsite), target, save fparm
Definition: parcom.f90:54
double precision, dimension(mxvar), save xfdstp
Definition: parcom.f90:56
integer, parameter mxspc
Definition: nlsdim.f90:39
integer, dimension(mxvar), save ixpr
Definition: parcom.f90:62
integer, save nprm
Definition: parcom.f90:62
double precision, dimension(mxvar), save prmin
Definition: parcom.f90:56
integer, save nsite
Definition: parcom.f90:62
integer, parameter mxsite
Definition: nlsdim.f90:39
integer, save nser
Definition: parcom.f90:62
double precision, dimension(mxvar), save prmax
Definition: parcom.f90:56
Definition: lpnam.f90:18
Definition: lmcom.f90:13
integer, parameter luttyo
Definition: stdio.f90:29
integer, dimension(mxvar), save ibnd
Definition: parcom.f90:62
integer, dimension(nfprm, mxsite), save ixx
Definition: parcom.f90:62
subroutine addprm(ix, ix2, ibd, prmn, prmx, prsc, step, ident)
Add a parameter to list of parameters being varied for nonlinear least-squares. Also maintain the lis...
Definition: addprm.f90:54
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67
double precision, parameter rndoff
Definition: rnddbl.f90:86
double precision, dimension(mxvar), save prscl
Definition: parcom.f90:56
double precision, dimension(mxjcol), save x
Definition: lmcom.f90:17