NLSL
srchc.f90
Go to the documentation of this file.
1 c Version 1.5.1b 11/6/95
2 c----------------------------------------------------------------------
3 c =========================
4 c subroutine SRCHC
5 c =========================
6 c
7 c search <parameter> { xtol <xtol> step <step> bound <bound> }
8 c
9 c xtol: Convergence tolerance for fitting parameter
10 c step: Size of first step to take away from initial parameter value
11 c bound: Boundary for search (may not exceed initial value
12 c +/- bound
13 c maxfun: Maximum number of function evaluations
14 c
15 c Uses
16 c mnbrak
17 c brent
18 c l1pfun
19 c setprm
20 c gettkn
21 c ftoken
22 c ipfind
23 c spcpar
24 c tcheck
25 c indtkn
26 c itrim
27 c----------------------------------------------------------------------
28  subroutine srchc( line )
29 c
30  use nlsdim
31  use eprprm
32  use expdat
33  use lmcom
34  use parcom
35  use iterat
36  use errmsg
37  use lpnam
38  use stdio
39 c
40  implicit none
41  character line*80
42 c
43  integer i,itmp,bflag,jx,jx1,jx2,lth,lu
44  double precision ax,bx,cx,fa,fb,fc,fmin,fval,pmin
45  character token*30,prmID*9,tagtmp*9
46 c
47  integer NKEYWD
48  parameter(nkeywd=6)
49 c
50  integer ipfind,indtkn,itrim
51  double precision l1pfun,brent,getprm,gammq,residx,corrl,wtdres
52  logical ftoken,tcheck,spcpar
53  external l1pfun,brent,ftoken,spcpar,tcheck,
54  # ipfind,indtkn,itrim,getprm,residx,corrl,wtdres
55 c
56  character*8 keywrd(nkeywd)
57  data keywrd / 'XTOL','FTOL','STEP','BOUND','MAXFUN','SRANGE'/
58 c
59 c######################################################################
60 c
61 c----------------------------------------------------------------------
62 c Get the name of the parameter
63 c----------------------------------------------------------------------
64  call gettkn(line,token,lth)
65  lth=min(lth,6)
66  call touppr(token,lth)
67 c
68  ixp1p=ipfind(token,lth)
69 c
70 c----------------------------------------------------------------------
71 c Check whether parameter may be varied
72 c----------------------------------------------------------------------
73  if (ixp1p.eq.0 .or. ixp1p.gt.nfprm) then
74  write(luttyo,1002) token(:lth)
75  return
76  end if
77 c
78  if (ixp1p.eq.iser) then
79  write(luttyo,1011) token(:lth)
80  return
81  end if
82 c
83  if (ixp1p.lt.-100) then
84  prmid=alias2( -99-(iwxx+ixp1p) )
85  else if (ixp1p.lt.0) then
86  prmid=alias1( 1-(iwxx+ixp1p) )
87  else
88  prmid=parnam(ixp1p)
89  end if
90 c
91 c --- Get secondary index
92 c
93  ixs1p=indtkn( line )
94 c
95 c --- Set bounds for site index of search parameter
96 c
97  if (ixs1p.le.0) then
98 c
99 c --- if parameter is to be searched globally for all sites,
100 c check that the initial values of the parameter are the
101 c same for all sites
102 c
103  do i=2,nsite
104  if (getprm(ixp1p,i).ne.getprm(ixp1p,1)) then
105  write (luout,1004) prmid(:itrim(prmid))
106  if (luout.ne.luttyo)
107  # write (luttyo,1004) prmid(:itrim(prmid))
108  return
109  end if
110  end do
111 c
112  jx1=1
113  jx2=mxsite
114  if (spcpar(ixp1p)) jx2=mxspc
115  else
116  jx1=ixs1p
117  jx2=ixs1p
118  write (prmid(itrim(prmid)+1:),1005) ixs1p
119  end if
120 c
121 c --- Check whether proper symmetry has been specified for
122 c tensor components
123 c
124  lu=luttyo
125  do jx=jx1,jx2
126  if (.not.tcheck(ixp1p,ixs1p,prmid,lu)) return
127  lu=0
128  end do
129 c
130  ixp1p=iabs( mod(ixp1p,100) )
131 c
132 c----------------------------------------------------------------------
133 c Look for a keyword
134 c----------------------------------------------------------------------
135 c
136  13 call gettkn(line,token,lth)
137  lth=min(lth,8)
138 c
139 c------------------------------------------------
140 c************************************************
141 c No more keywords: call the line search routine
142 c************************************************
143 c------------------------------------------------
144  if (lth.eq.0) then
145 c
146 c ---------------------------------------------------
147 c Find a starting place for the bracketing procedure
148 c ---------------------------------------------------
149 c
150  if (ixs1p.le.0) then
151  ax=fparm(ixp1p,1)
152  else
153  ax=fparm(ixp1p,ixs1p)
154  end if
155  bx=ax+pstep
156 c
157 c -------------------------------------------------------------------
158 c Swap indices for the search parameter into the first elements of
159 c the parameter and site index arrays for the NLS x-vector
160 c (This is so x can be used by the lfun routine as in a NLS procedure)
161 c -------------------------------------------------------------------
162 c NOTE: Nonsense is saved if no paramaters have been varied!
163 c The swap-out should occur only if nprm .gt. 1
164 c
165  itmp=ixpr(1)
166  ixpr(1)=ixp1p
167  ixp1p=itmp
168  itmp=ixst(1)
169  ixst(1)=ixs1p
170  ixs1p=itmp
171  tagtmp=tag(1)
172  tag(1)=prmid
173 c
174 c --------------------
175 c Bracket the minimum
176 c --------------------
177 c
178  call catchc( hltfit )
179  warn=.true.
180 c
181  write (luout,1009) prmid(:itrim(prmid))
182  if (luout.ne.luttyo) write (luttyo,1009) prmid(:itrim(prmid))
183 c
184  bflag=1
185  iter=0
186  fnmin=0.0d0
187  call mnbrak(ax,bx,cx,fa,fb,fc,l1pfun,bflag,pbound)
188 c
189 c ---------------------------------
190 c User-terminated during bracketing
191 c ---------------------------------
192  if (bflag.lt.0) then
193  write (luttyo,1007)
194  go to 15
195  else if (bflag.eq.nomin) then
196  write (luttyo,1012) prmid(:itrim(prmid)),cx
197  call setprm( ixpr(1), ixst(1), cx )
198  fnorm=fc
199  go to 14
200  end if
201 c
202 c -----------------------------------
203 c Use Brent's method for minimization
204 c -----------------------------------
205 c
206  write (luout,1010) ax,cx
207  if (luout.ne.luttyo) write (luttyo,1010) ax,cx
208  bflag=1
209  iter=0
210  fnmin=fb
211  fmin=brent(ax,bx,cx,fb,l1pfun,ptol,pftol,pmin,bflag)
212 c
213  if (bflag.lt.0) then
214 c
215 c ---------------------------------
216 c User-terminated:report best value
217 c ---------------------------------
218  write (luttyo,1008) prmid,pmin
219  if (luout.ne.luttyo) write (luout,1008) prmid,pmin
220  else
221 c
222 c ------------------------
223 c Minimum found: report it
224 c ------------------------
225  write (luttyo,1006) prmid,pmin
226  if (luout.ne.luttyo) write (luout,1006) prmid,pmin
227  end if
228 c
229  fnorm=fmin
230  call setprm( ixpr(1), ixst(1), pmin )
231 c
232  14 if (iwflag.ne.0) then
234  rdchsq=chisqr/float(ndatot-nprm)
235  qfit=gammq( 0.5d0*(ndatot-nprm), 0.5*chisqr )
236  else
237  chisqr=wtdres( fvec,ndatot,nspc,ixsp,npts,rmsn )**2
238  rdchsq=chisqr/float(ndatot)
239  end if
240 c
241  write(luout,2036) fnorm,chisqr,rdchsq,corrl(),residx()
242  if (iwflag.ne.0) write (luout,2038) qfit
243 c
244  if (luout.ne.luttyo)
245  # write(luttyo,2036) fnorm,chisqr,rdchsq,corrl(),
246  # residx()
247 c
248 c --------------------------------------------------------
249 c Restore the parameter/site index arrays and I.D. string
250 c --------------------------------------------------------
251  15 itmp=ixpr(1)
252  ixpr(1)=ixp1p
253  ixp1p=itmp
254 c
255  itmp=ixst(1)
256  ixst(1)=ixs1p
257  ixs1p=itmp
258  tag(1)=tagtmp
259 c
260  call uncatchc( hltfit )
261 c
262  return
263  end if
264 c
265 c -------------------
266 c Check keyword list
267 c -------------------
268  call touppr(token,lth)
269  do i=1,nkeywd
270  if (token(:lth).eq.keywrd(i)(:lth)) goto 16
271  end do
272 c *** Unrecognized keyword
273  write (luttyo,1000) token(:lth)
274  return
275 c
276 c ----------------------------------------------------------
277 c Keyword found: for keywords requiring an argument, convert
278 c next token and assign appropriate value
279 c ------------------------------------------------------------
280  16 call gettkn(line,token,lth)
281 c *** No value given
282  if (lth.eq.0) then
283  write(luttyo,1003) keywrd(i)
284  return
285  end if
286 c
287  if (ftoken(token,lth,fval)) then
288 c *** XTOL keyword
289  if (i.eq.1) then
290  ptol=fval
291 c *** FTOL keyword
292  else if (i.eq.2) then
293  pftol=fval
294 c *** STEP keyword
295  else if (i.eq.3) then
296  pstep=fval
297 c *** BOUND keyword
298  else if (i.eq.4) then
299  pbound=fval
300 c *** MAXFUN keyword
301  else if (i.eq.5) then
302  mxpitr=int(fval)
303 c *** MAXFUN keyword
304  else if (i.eq.6) then
305  srange=fval/100.0
306  end if
307 c *** Illegal numeric value
308  else
309  write(luttyo,1001) token(:lth)
310  end if
311  go to 13
312 c
313 c ##### Format statements ########################################
314 c
315  1000 format('*** Unrecognized SEARCH keyword: ''',a,''' ***')
316  1001 format('*** Numeric value expected: ''',a,''' ***')
317  1002 format('*** ''',a,''' is not a variable parameter ***')
318  1003 format('*** No value given for ''',a,''' ***')
319  1004 format('*** ',a,' is not the same for all currently defined',
320  # ' sites ***')
321  1005 format('(',i1,')')
322  1006 format(/2x,'Minimum found at ',a,'= ',g12.6/)
323  1007 format('*** Terminated by user during bracketing of minimum ***')
324  1008 format('*** Terminated by user ***'/'Best point is ',a,'= ',f9.5/)
325  1009 format(/2x,'Bracketing the minimum in ',a)
326  1010 format(/2x,'Minimum is between ',f9.5,' and ',f9.5/)
327  1011 format('*** ',a,' is series parameter: cannot be searched ***')
328  1012 format('*** Minimum is at step bound ***'/2x,a,'=',f9.5/)
329  2036 format(10x,'Residual norm= ',g13.6/
330  # 10x,'Chi-squared=',g13.6,5x,'Reduced Chi-sq=',g13.6/
331  # 10x,'Correlation = ',f8.5,5x,'Residual index =',f8.5/)
332  2038 format(12x,'Goodness of fit from chi-squared (Q)=',g13.6)
333  end
integer, save luout
Definition: stdio.f90:32
integer, dimension(mxvar), save ixst
Definition: parcom.f90:62
subroutine srchc(line)
Definition: srchc.f90:29
double precision, save srange
Definition: parcom.f90:56
integer, save ndatot
Definition: expdat.f90:45
integer, dimension(mxspc), save ixsp
Definition: expdat.f90:45
void FORTRAN() uncatchc(int *flag)
Definition: catch.c:49
character *6, dimension(nalias), save alias1
Definition: lpnam.f90:45
double precision, save fnmin
Definition: iterat.f90:15
integer, save iter
Definition: iterat.f90:14
double precision, save qfit
Definition: iterat.f90:15
integer, parameter nomin
Definition: errmsg.f90:51
integer, save hltfit
Definition: stdio.f90:32
double precision, save pstep
Definition: parcom.f90:56
double precision, save rdchsq
Definition: iterat.f90:15
Definition: stdio.f90:26
double precision, dimension(nfprm, mxsite), target, save fparm
Definition: parcom.f90:54
double precision, save pftol
Definition: parcom.f90:56
subroutine touppr(string, lth)
Definition: strutl2.f90:22
double precision, save fnorm
Definition: iterat.f90:15
integer, save ixs1p
Definition: parcom.f90:62
integer, save iwflag
Definition: iterat.f90:14
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 nspc
Definition: expdat.f90:45
integer, save nsite
Definition: parcom.f90:62
integer, save ixp1p
Definition: parcom.f90:62
integer, parameter mxsite
Definition: nlsdim.f90:39
double precision, dimension(mxpt), save fvec
Definition: lmcom.f90:17
integer, parameter iwxx
Definition: eprprm.f90:92
void FORTRAN() catchc(int *flag)
Definition: catch.c:33
Definition: lpnam.f90:18
Definition: lmcom.f90:13
integer, save iser
Definition: parcom.f90:62
integer, parameter luttyo
Definition: stdio.f90:29
double precision, save pbound
Definition: parcom.f90:56
double precision, save ptol
Definition: parcom.f90:56
double precision, dimension(mxspc), save rmsn
Definition: expdat.f90:40
integer, save mxpitr
Definition: parcom.f90:62
logical, save warn
Definition: stdio.f90:33
subroutine mnbrak(ax, bx, cx, fa, fb, fc, func, iflag, bound)
Definition: mnbrak.f90:2
integer, dimension(mxspc), save npts
Definition: expdat.f90:45
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
character *6, dimension(nfprm), save parnam
Definition: lpnam.f90:27
character *9, dimension(mxjcol), save tag
Definition: parcom.f90:67
integer, parameter nfprm
Definition: nlsdim.f90:57
double precision, save chisqr
Definition: iterat.f90:15
double precision, dimension(mxjcol), save x
Definition: lmcom.f90:17
character *6, dimension(nalias), save alias2
Definition: lpnam.f90:51