/[MITgcm]/MITgcm_contrib/ecco_utils/ecco_v4_release3_optimization/lsopt/lsline.F
ViewVC logotype

Annotation of /MITgcm_contrib/ecco_utils/ecco_v4_release3_optimization/lsopt/lsline.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.2 - (hide annotations) (download)
Wed Jan 17 18:01:05 2018 UTC (7 years, 6 months ago) by ou.wang
Branch: MAIN
CVS Tags: HEAD
Changes since 1.1: +48 -11 lines
Update README and add more printouts

1 ou.wang 1.1
2     subroutine lsline(
3     & simul
4     & , nn, ifail, lphprint
5     & , ifunc, nfunc
6     & , ff, dotdg
7     & , tmin, tmax, tact, epsx
8     & , dd, gg, xx, xdiff
9     & )
10    
11     c ==================================================================
12     c SUBROUTINE lsline
13     c ==================================================================
14     c
15     c o line search algorithm for determining control vector update;
16     c After computing updated control vector from given gradient,
17     c a forward and adjoint model run are performed (simul.F)
18     c using the updated control vector.
19     c Tests are then applied to see whether solution has improved.
20     c
21     c o Reference: J.C. Gilbert & C. Lemarechal
22     c Some numerical experiments with variable-storage
23     c quasi-Newton algorithms
24     c Mathematical Programming 45 (1989), pp. 407-435
25     c
26     c o started: ??? not reproducible
27     c
28     c o changed: Patrick Heimbach, MIT/EAPS
29     c
30     c o Version: 2.0, 24-Feb-2000: Patrick Heimbach, MIT/EAPS
31     c - severe changes in structure including various
32     c shifts of variables which are only used in this
33     c routine
34     c - number of 3 control flags for error handling
35     c (indic, moderl, ifail) reduced to 1 (ifail)
36     c and homogenized with routine lsoptv
37     c
38     c o Version: 2.1.0, 02-Mar-2000: Patrick Heimbach, MIT/EAPS
39     c - initial computation of tact and
40     c xdiff = xx + tact*dd
41     c moved to new routine lsupdxx
42     c tmin, tmax, tact needed as parameters
43     c
44     c ==================================================================
45     c SUBROUTINE lsline
46     c ==================================================================
47    
48     #include "blas1.h"
49    
50     implicit none
51    
52     c----------------------------------
53     c declare arguments
54     c----------------------------------
55     integer nn, ifail, ifunc, nfunc
56     double precision ff, dotdg, tmin, tmax, tact, epsx
57     double precision xx(nn), dd(nn), gg(nn), xdiff(nn)
58    
59     logical lphprint
60    
61     external simul
62    
63     c----------------------------------
64     c declare local variables
65     c----------------------------------
66    
67     double precision xpara1, xpara2
68     parameter( xpara1 = 0.000001, xpara2=0.99999 )
69    
70     double precision factor
71     parameter( factor = 0.2 )
72    
73     double precision barmax
74     parameter( barmax = 0.3 )
75     double precision barmul
76     parameter( barmul = 5.0 )
77     double precision barmin
78     parameter( barmin = 0.01 )
79    
80     integer i, indic
81    
82     double precision tg, fg, td, ta
83     double precision fa, fpa, fp
84     double precision fnew, fdiff
85     double precision z, test, barr
86     double precision left, right, told
87    
88     external DDOT
89     double precision DDOT
90    
91     c----------------------------------
92     c check parameters
93     c----------------------------------
94    
95     if ( (nn.le.0)
96     & .or. (dotdg.ge.0.0)
97     & .or. (xpara1.le.0.0) .or. (xpara1.ge.0.5)
98     & .or. (xpara2.le.xpara1) .or. (xpara2.ge.1.0) ) then
99     ifail = 9
100     go to 999
101     endif
102    
103     c----------------------------------
104     c initialization
105     c----------------------------------
106     indic = 0
107    
108     barr = barmin
109     fg = ff
110     fa = ff
111     fpa = dotdg
112    
113     td = 0.0
114     tg = 0.0
115     ta = 0.0
116    
117     c=======================================================================
118     c begin of simulation iter.
119     c=======================================================================
120    
121     do ifunc = 1, nfunc
122    
123     if (lphprint)
124     & print *, 'pathei-lsopt: ', ifunc, ' simul.'
125    
126     c------------------------------------
127 ou.wang 1.2 cow( for offline mode (loffline.eq..TRUE.), simul reads
128     c in ecco_ctrl and ecco_cost from iteration optimcycle.
129     cow) There are no actual forward/adjoint runs in simul.
130 ou.wang 1.1 c------------------------------------
131     call simul( indic, nn, xdiff, fnew, gg )
132    
133 ou.wang 1.2 cow( compute tact (could be different from 1)
134 ou.wang 1.1 do i = 1, nn
135     if(dd(i) .ne. 0.and.(xdiff(i) - xx(i)).ne.0) then
136     tact = (xdiff(i) - xx(i))/ dd(i)
137     goto 301
138     endif
139     end do
140     301 continue
141 ou.wang 1.2 c if tact is only slightly different from 1., then set it to 1.
142     c The difference is likely caused by trunction error.
143 ou.wang 1.1 if(abs(tact-1.).le.1.d-6) tact = 1.
144 ou.wang 1.2 write(*,'(a,e15.6)')
145     & ' pathei-lsopt: step size (tact) for the latest iteration
146     &: ',tact
147     cow) now we know tact (step size for the latest iteration).
148 ou.wang 1.1
149     C compute fp = direction(dd) . gradient(gg)
150     fp = DDOT( nn, dd, 1, gg, 1 )
151     fdiff = fnew - ff
152    
153     print*, 'pathei-lsopt: '
154 ou.wang 1.2 cow( print out info related to Wolfe tests
155     print*,'=================================================='
156     print*, 'Wolfe test 1 (sufficient cost decrease condition):'
157     write(*,'(a,3e15.6)') ' Cost (iter i, iter i-1, reduction): ',
158     & fnew, ff, fdiff
159     write(*,'(a,a,3e15.6)')' alpha1, step size for iter i, ',
160     & '<dd(i-1).gg(i-1)> ',
161     & xpara1, tact, dotdg
162 ou.wang 1.1 if(fdiff .gt. tact*xpara1*dotdg) then
163 ou.wang 1.2 write(*,'(a,a,e15.6,a,e15.6)')
164     & ' cost reduction GT alpha1*step size',
165     & '*<dd(i-1).gg(i-1)> ', fdiff, ' > ', tact*xpara1*dotdg
166 ou.wang 1.1 print*, 'Wolfe test 1: failed'
167     else
168 ou.wang 1.2 write(*,'(a,a,e15.6,a,e15.6)')
169     & ' cost reduction LE alpha1*step size',
170     & '*<dd(i-1).gg(i-1)> ', fdiff, ' <= ', tact*xpara1*dotdg
171 ou.wang 1.1 print*, 'Wolfe test 1: passed'
172     endif
173 ou.wang 1.2 print*,''
174 ou.wang 1.1
175 ou.wang 1.2 print*, 'Wolfe test 2 (curvature condition):'
176     write(*,'(a,3e15.6)')' <dd(i-1).gg(i)>, <dd(i-1).gg(i-1)>: ',
177     & fp, dotdg
178     write(*,'(a,a,3e15.6)')
179     & ' alpha2,', 'alpha2*<dd(i-1).gg(i-1)>: ',
180     & xpara2, xpara2*dotdg
181 ou.wang 1.1 if (fp .gt. xpara2*dotdg) then
182 ou.wang 1.2 write(*,'(a,e15.6,a,e15.6)')
183     & ' <dd(i-1).gg(i)> GT alpha2* <dd(i-1).gg(i-1)>: ',
184     & fp, ' > ', xpara2*dotdg
185 ou.wang 1.1 print*, 'Wolfe test 2: passed'
186     else
187 ou.wang 1.2 write(*,'(a,e15.6,a,e15.6)')
188     & ' <dd(i-1).gg(i)> LE alpha2* <dd(i-1).gg(i-1)>: ',
189     & fp, ' <= ', xpara2*dotdg
190 ou.wang 1.1 print*, 'Wolfe test 2: failed'
191     endif
192 ou.wang 1.2 print*, '=================================================='
193     cow)
194 ou.wang 1.1
195     c-----------------------------------------
196     c apply 1st Wolfe test
197     c-----------------------------------------
198    
199     if (fdiff .gt. tact*xpara1*dotdg) then
200 ou.wang 1.2 C if (lphprint)
201     C & print *, 'Wolfe test 1 (Armijo Rule) Failed'
202 ou.wang 1.1 if (lphprint)
203 ou.wang 1.2 & print *, 'pathei-lsopt: interpolate to get new step size'
204 ou.wang 1.1 td = tact
205     ifail = 0
206     go to 500
207     endif
208    
209     c-----------------------------------------
210     c 1st Wolfe test 1 ok, apply 2nd Wolf test
211     c-----------------------------------------
212     if (fp .gt. xpara2*dotdg) then
213 ou.wang 1.2 C if (lphprint)
214     C & print *, 'Pass Wolfe test 2 (Curvature condition)'
215 ou.wang 1.1 ifail = 0
216     go to 320
217     endif
218    
219     if (ifail.eq.0) go to 350
220    
221     c-----------------------------------------
222     c 2nd Wolfe test 2 ok, donc pas serieux, on sort
223     c-----------------------------------------
224    
225     320 continue
226    
227     ff = fnew
228     do i = 1, nn
229     xx(i) = xdiff(i)
230     end do
231     cph(
232     if (lphprint)
233     & print *, 'pathei-lsopt: no inter-/extrapolation in lsline'
234     cph)
235     go to 888
236    
237    
238     c-----------------------------------------
239     c extrapolation
240     c-----------------------------------------
241    
242     350 continue
243    
244 ou.wang 1.2 if (lphprint)
245     & print *, 'pathei-lsopt: extrapolate to get new step size'
246 ou.wang 1.1 tg = tact
247     fg = fnew
248     if (td .ne. 0.0) go to 500
249    
250     told = tact
251     left = (1.0+barmin)*tact
252     right = 10.0*tact
253     call cubic( tact, fnew, fp, ta, fa, fpa, left, right )
254     ta = told
255     if (tact.ge.tmax) then
256     ifail = 7
257     tact = tmax
258     endif
259    
260     if (lphprint)
261     & print *, 'pathei-lsopt: extrapolation: ',
262     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
263    
264     go to 900
265    
266     c-----------------------------------------
267     c interpolation
268     c-----------------------------------------
269    
270     500 continue
271    
272     C td: tact (step size)
273     C tg: was set to zero
274     C So the limits are between 0.01 * tg and tact-0.01*tg (0 and tact for isuml = 1)
275     test = barr*(td-tg)
276     left = tg+test
277     right = td-test
278     told = tact
279     C input:
280     C tact: was set to 1 in input, but in output store the new step size
281     C fnew: the new cost function in input
282     C fp: new dotdg = dd . gg ( direction DOT gradient)
283     C ta: set to zero (namely set x0 = 0 so tact is distance between x1 and x0)
284     C fa: old cost function
285     C fpa: old dotdg = dd . gg (diretion DOT gradient)
286     call cubic( tact, fnew, fp, ta, fa, fpa, left, right )
287     ta = told
288     if (tact.gt.left .and. tact.lt.right) then
289     C falling within left and right brackets. Keep the same barr
290     barr = dmax1( barmin, barr/barmul )
291     else
292     C outside of the left and right brackets. Reduce barr by half
293     barr = dmin1( barmul*barr, barmax )
294     endif
295    
296     if (lphprint)
297     & print *, 'pathei-lsopt: interpolation: ',
298     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
299    
300     c-----------------------------------------
301     c end of iteration loop
302     c - tact peut etre bloque sur tmax
303     c (venant de lextrapolation avec ifail=7)
304     c-----------------------------------------
305    
306     900 continue
307    
308     fa = fnew
309     fpa = fp
310    
311     c-----------------------------------------
312     c --- faut-il continuer ?
313     c-----------------------------------------
314     if (td .eq. 0.0) go to 950
315     if (td-tg .lt. tmin) go to 920
316    
317     c-----------------------------------------
318     c limit due to machine precision
319     c-----------------------------------------
320     do i = 1, nn
321     z = xx(i) + tact*dd(i)
322     if ((z.ne.xx(i)) .and. (z.ne.xdiff(i))) go to 950
323     end do
324    
325     c-----------------------------------------
326     c arret sur dxmin ou de secours
327     c-----------------------------------------
328     920 continue
329     ifail = 8
330    
331     c-----------------------------------------
332     c si tg=0, xx = xx_depart,
333     c sinon on prend xx=x_left qui fait decroitre f
334     c-----------------------------------------
335     if (tg .ne. 0.0) then
336     ff = fg
337     do i = 1, nn
338     xx(i) = xx(i) + tg*dd(i)
339     end do
340     endif
341    
342     go to 999
343    
344     c-----------------------------------------
345     c update vector for new simulation
346     c-----------------------------------------
347     950 continue
348    
349     do i = 1, nn
350     xdiff(i) = xx(i) + tact*dd(i)
351     end do
352    
353     c=======================================================================
354     c end of simulation iter.
355     c=======================================================================
356     if (lphprint)
357     & print *, 'pathei-lsopt: end of simulation iter: ',
358     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
359     end do
360    
361     c-----------------------------------------
362     c too many function evaluations
363     c-----------------------------------------
364     ifail = 6
365     ff = fg
366     do i = 1, nn
367     xx(i) = xx(i) + tg*dd(i)
368     end do
369    
370    
371     999 continue
372     C Set ifail=99 so we do not modify OPWARMD and OPWARMI
373 ou.wang 1.2 C We modify OPWARMD and OPWARMI Only if both Wolfe
374     C conditions are satisfied.
375 ou.wang 1.1
376     ifail = 99
377    
378     c-----------------------------------------
379     c end of routine
380     c-----------------------------------------
381     888 continue
382    
383     if (lphprint)
384     & print *, 'pathei-lsopt: end of lsline: ',
385     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
386     return
387    
388     end

  ViewVC Help
Powered by ViewVC 1.1.22