/[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.1 - (hide annotations) (download)
Wed Jan 3 17:13:47 2018 UTC (7 years, 6 months ago) by ou.wang
Branch: MAIN
Check in the optimization used in ECCO v4r3

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     c compute cost function and gradient
128     c------------------------------------
129     call simul( indic, nn, xdiff, fnew, gg )
130    
131     C compute tact (could be different from 1)
132     do i = 1, nn
133     if(dd(i) .ne. 0.and.(xdiff(i) - xx(i)).ne.0) then
134     tact = (xdiff(i) - xx(i))/ dd(i)
135     goto 301
136     endif
137     end do
138     301 continue
139     C if tact is only slightly different from 1., then it
140     C is set to 1. The difference is likely caused by trunction error.
141     if(abs(tact-1.).le.1.d-6) tact = 1.
142    
143    
144     C compute fp = direction(dd) . gradient(gg)
145     fp = DDOT( nn, dd, 1, gg, 1 )
146     fdiff = fnew - ff
147    
148     print*, 'pathei-lsopt: '
149     if(fdiff .gt. tact*xpara1*dotdg) then
150     print*, 'Wolfe test 1: failed'
151     else
152     print*, 'Wolfe test 1: passed'
153     endif
154    
155     if (fp .gt. xpara2*dotdg) then
156     print*, 'Wolfe test 2: passed'
157     else
158     print*, 'Wolfe test 2: failed'
159     endif
160    
161    
162     c-----------------------------------------
163     c apply 1st Wolfe test
164     c-----------------------------------------
165    
166     if (fdiff .gt. tact*xpara1*dotdg) then
167     if (lphprint)
168     & print *, 'Wolfe test 1 (Armijo Rule) Failed'
169     td = tact
170     ifail = 0
171     go to 500
172     endif
173    
174     c-----------------------------------------
175     c 1st Wolfe test 1 ok, apply 2nd Wolf test
176     c-----------------------------------------
177     if (fp .gt. xpara2*dotdg) then
178     if (lphprint)
179     & print *, 'Pass Wolfe test 2 (Curvature condition)'
180     ifail = 0
181     go to 320
182     endif
183    
184     if (ifail.eq.0) go to 350
185    
186     c-----------------------------------------
187     c 2nd Wolfe test 2 ok, donc pas serieux, on sort
188     c-----------------------------------------
189    
190     320 continue
191    
192     ff = fnew
193     do i = 1, nn
194     xx(i) = xdiff(i)
195     end do
196     cph(
197     if (lphprint)
198     & print *, 'pathei-lsopt: no inter-/extrapolation in lsline'
199     cph)
200     go to 888
201    
202    
203     c-----------------------------------------
204     c extrapolation
205     c-----------------------------------------
206    
207     350 continue
208    
209     tg = tact
210     fg = fnew
211     if (td .ne. 0.0) go to 500
212    
213     told = tact
214     left = (1.0+barmin)*tact
215     right = 10.0*tact
216     call cubic( tact, fnew, fp, ta, fa, fpa, left, right )
217     ta = told
218     if (tact.ge.tmax) then
219     ifail = 7
220     tact = tmax
221     endif
222    
223     if (lphprint)
224     & print *, 'pathei-lsopt: extrapolation: ',
225     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
226    
227     go to 900
228    
229     c-----------------------------------------
230     c interpolation
231     c-----------------------------------------
232    
233     500 continue
234    
235     C td: tact (step size)
236     C tg: was set to zero
237     C So the limits are between 0.01 * tg and tact-0.01*tg (0 and tact for isuml = 1)
238     test = barr*(td-tg)
239     left = tg+test
240     right = td-test
241     told = tact
242     C input:
243     C tact: was set to 1 in input, but in output store the new step size
244     C fnew: the new cost function in input
245     C fp: new dotdg = dd . gg ( direction DOT gradient)
246     C ta: set to zero (namely set x0 = 0 so tact is distance between x1 and x0)
247     C fa: old cost function
248     C fpa: old dotdg = dd . gg (diretion DOT gradient)
249     call cubic( tact, fnew, fp, ta, fa, fpa, left, right )
250     ta = told
251     if (tact.gt.left .and. tact.lt.right) then
252     C falling within left and right brackets. Keep the same barr
253     barr = dmax1( barmin, barr/barmul )
254     else
255     C outside of the left and right brackets. Reduce barr by half
256     barr = dmin1( barmul*barr, barmax )
257     endif
258    
259     if (lphprint)
260     & print *, 'pathei-lsopt: interpolation: ',
261     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
262    
263     c-----------------------------------------
264     c end of iteration loop
265     c - tact peut etre bloque sur tmax
266     c (venant de lextrapolation avec ifail=7)
267     c-----------------------------------------
268    
269     900 continue
270    
271     fa = fnew
272     fpa = fp
273    
274     c-----------------------------------------
275     c --- faut-il continuer ?
276     c-----------------------------------------
277     if (td .eq. 0.0) go to 950
278     if (td-tg .lt. tmin) go to 920
279    
280     c-----------------------------------------
281     c limit due to machine precision
282     c-----------------------------------------
283     do i = 1, nn
284     z = xx(i) + tact*dd(i)
285     if ((z.ne.xx(i)) .and. (z.ne.xdiff(i))) go to 950
286     end do
287    
288     c-----------------------------------------
289     c arret sur dxmin ou de secours
290     c-----------------------------------------
291     920 continue
292     ifail = 8
293    
294     c-----------------------------------------
295     c si tg=0, xx = xx_depart,
296     c sinon on prend xx=x_left qui fait decroitre f
297     c-----------------------------------------
298     if (tg .ne. 0.0) then
299     ff = fg
300     do i = 1, nn
301     xx(i) = xx(i) + tg*dd(i)
302     end do
303     endif
304    
305     go to 999
306    
307     c-----------------------------------------
308     c update vector for new simulation
309     c-----------------------------------------
310     950 continue
311    
312     do i = 1, nn
313     xdiff(i) = xx(i) + tact*dd(i)
314     end do
315    
316     c=======================================================================
317     c end of simulation iter.
318     c=======================================================================
319     if (lphprint)
320     & print *, 'pathei-lsopt: end of simulation iter: ',
321     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
322     end do
323    
324     c-----------------------------------------
325     c too many function evaluations
326     c-----------------------------------------
327     ifail = 6
328     ff = fg
329     do i = 1, nn
330     xx(i) = xx(i) + tg*dd(i)
331     end do
332    
333    
334     999 continue
335     C Set ifail=99 so we do not modify OPWARMD and OPWARMI
336     C We modify OPWARMD and OPWARMI Only if the result
337     C passes both Wolfe condtions
338    
339     ifail = 99
340    
341     c-----------------------------------------
342     c end of routine
343     c-----------------------------------------
344     888 continue
345    
346     if (lphprint)
347     & print *, 'pathei-lsopt: end of lsline: ',
348     & 'td, tg, tact, ifail = ', td, tg, tact, ifail
349     return
350    
351     end

  ViewVC Help
Powered by ViewVC 1.1.22