/[MITgcm]/MITgcm_contrib/mlosch/optim_m1qn3/m1qn3_offline.F
ViewVC logotype

Contents of /MITgcm_contrib/mlosch/optim_m1qn3/m1qn3_offline.F

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


Revision 1.1 - (show annotations) (download)
Thu Apr 26 11:10:06 2012 UTC (12 years ago) by mlosch
Branch: MAIN
First working version of a new optimization package that uses a slightly
modified version of m1qn3, v3.3
(https://who.rocq.inria.fr/Jean-Charles.Gilbert/modulopt/optimization-routines/m1qn3/m1qn3.html)
to work as an offline optimizer. The advantage of m1qn3_offline is, that
it is run in reverse communication control mode, so that it gives back
control to the call routine (here a script) to provide a new estimate of the
cost function and the gradient based on the control vector. This way we can
do complete line searches that are meaningful.

1 subroutine m1qn3_offline
2 & (simul,prosca,ctonb,ctcab,n,x,f,g,dxmin,df1,
3 & epsg,normtype,impres,io,imode,omode,niter,nsim,
4 & iz,dz,ndz,reverse,indic,izs,rzs,dzs)
5 c
6 c-----------------------------------------------------------------------
7 c
8 c M1QN3, Version 3.3, October 2009
9 c
10 c M1qn3 has two running modes: the SID (Scalar Initial Scaling) mode
11 c and the DIS (Diagonal Initial Scaling) mode. Both do not require
12 c the same amount of storage, the same subroutines, ...
13 c In the description below, items that differ in the DIS mode with
14 c respect to the SIS mode are given in brakets.
15 c
16 c Use the following subroutines:
17 c M1QN3A
18 c DDD, DDDS
19 c NLIS0 + DCUBE (Dec 88)
20 c MUPDTS, DYSTBL.
21 c
22 c The following routines are proposed to the user in case the
23 c Euclidean scalar product is used:
24 c DUCLID, DTONBE, DTCABE.
25 c
26 c La sous-routine M1QN3 est une interface entre le programme
27 c appelant et la sous-routine M1QN3A, le minimiseur proprement dit.
28 c
29 c Le module PROSCA est sense realiser le produit scalaire de deux
30 c vecteurs de Rn; le module DTONB est sense realiser le changement
31 c de coordonnees correspondant au changement de bases: base
32 c euclidienne -> base orthonormale (pour le produit scalaire
33 c PROSCA); le module CTBAB fait la transformation inverse: base
34 c orthonormale -> base euclidienne.
35 c
36 c Iz is an integer working zone for M1QN3A, its dimension is 5.
37 c It is formed of 5 scalars that are set by the optimizer:
38 c - the dimension of the problem,
39 c - an identifier of the scaling mode,
40 c - the number of updates,
41 c - two pointers.
42 c
43 c Dz est la zone de travail pour M1QN3A, de dimension ndz.
44 c Elle est subdivisee en
45 c 3 [ou 4] vecteurs de dimension n: d,gg,[diag,]aux
46 c m scalaires: alpha
47 c m vecteurs de dimension n: ybar
48 c m vecteurs de dimension n: sbar
49 c
50 c m est alors le plus grand entier tel que
51 c m*(2*n+1)+3*n .le. ndz [m*(2*n+1)+4*n .le. ndz)]
52 c soit m := (ndz-3*n) / (2*n+1) [m := (ndz-4*n) / (2*n+1)].
53 c Il faut avoir m >= 1, donc ndz >= 5n+1 [ndz >= 6n+1].
54 c
55 c A chaque iteration la metrique est formee a partir d un multiple
56 c de l'identite [d'une matrice diagonale] D qui est mise a jour m
57 c fois par la formule de BFGS en utilisant les m couples {y,s} les
58 c plus recents.
59 c
60 c-----------------------------------------------------------------------
61 c
62 c Authors: Jean Charles Gilbert, Claude Lemarechal, INRIA.
63 c
64 c Copyright 2008, 2009, INRIA.
65 c
66 c M1QN3 is distributed under the terms of the GNU General Public
67 c License.
68 c
69 c-----------------------------------------------------------------------
70 c
71 c This program is free software: you can redistribute it and/or
72 c modify it under the terms of the GNU General Public License as
73 c published by the Free Software Foundation, either version 3 of
74 c the License, or (at your option) any later version.
75 c
76 c This program is distributed in the hope that it will be useful,
77 c but WITHOUT ANY WARRANTY; without even the implied warranty of
78 c MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
79 c General Public License for more details.
80 c
81 c You should have received a copy of the GNU General Public License
82 c along with this program. If not, see
83 c <http://www.gnu.org/licenses/>.
84 c
85 c-----------------------------------------------------------------------
86 c
87 implicit none
88 c
89 c arguments
90 c
91 character*3 normtype
92 integer n,impres,io,imode(3),omode,niter,nsim,iz(5),ndz,indic,
93 & reverse,izs(*)
94 real rzs(*)
95 double precision x(n),f,g(n),dxmin,df1,epsg,dz(ndz),dzs(*)
96 external simul,prosca,ctonb,ctcab
97 c
98 c --- local variables
99 c
100 c The variable 'reentry' is used to know where to jump when entering
101 c the subroutine in reverse commnunication (it is better than using
102 c 'reverse', which is known outside m1qn3 and could be changed by
103 c the user):
104 c = 0: no jump, start at the begining of the subroutine
105 c = 1: jump inside m1qn3a, where the simulator was called with
106 c indic=1
107 c = 2: jump inside m1qn3a/mlis3, where the simulator was called with
108 c indic=4
109 c The variable 'reentry' is set to 0 by the compiler and when m1qn3
110 c has completed its job (see below); since the value is saved, it
111 c has always the value 0 when entering the solver for the first
112 c time when solving a new problem, even if it is in the same
113 c program.
114 c
115 logical sscale
116 integer ntravu,m,mmemo
117 double precision gnorm,ps
118 #include "m1qn3_common.h"
119 CML logical inmemo,sscale
120 CML integer ntravu,id,igg,idiag,iaux,ialpha,iybar,isbar,m,mmemo,
121 CML & reentry
122 CML double precision gnorm,ps
123 CML save inmemo,sscale,ntravu,id,igg,idiag,iaux,ialpha,iybar,isbar,m,
124 CML & mmemo,reentry,gnorm,ps
125 c
126 c data
127 c
128 CML this needs to be done outside of m1qn3_offline
129 CML data reentry /0/
130 c
131 c --- function
132 c
133 double precision ddot,dnrmi
134 c
135 c --- stop if reverse < 0 (m1qn3 should not be called with reverse < 0)
136 c
137 if (reverse.lt.0) then
138 write (io,'(/a,a,i0,a/)')
139 & " >>> m1qn3 should not be called with a negative reverse",
140 & " (=", reverse, ")"
141 stop
142 endif
143 c
144 c --- possible jumps
145 c 9999: go directly in m1qn3a
146 c
147 if (reentry.gt.0) goto 9999
148 c
149 c---- license notice
150 c
151 if (impres.ge.5) then
152 write (io,934)
153 endif
154 934 format (1x,79("-")
155 & /1x,"M1QN3 Copyright (C) 2008, J. Ch. Gilbert, Cl. ",
156 & "Lemarechal."
157 & /1x,79("-")
158 & /1x,"This program comes with ABSOLUTELY NO WARRANTY. This is",
159 & " free software, and you"
160 & /1x,"are welcome to redistribute it under certain ",
161 & "conditions. See the file COPYING "
162 & /1x,"in the root directory of the M1QN3 distribution for ",
163 & "details."
164 & /1x,79("-"))
165 c
166 c---- impressions initiales et controle des arguments
167 c
168 if (impres.ge.1) then
169 write (io,900) n,dxmin,df1,epsg,normtype,niter,nsim,impres
170 if (reverse.gt.0) then
171 write (io,'(5x,a)') "reverse communication"
172 else
173 write (io,'(5x,a)') "direct communication"
174 endif
175 endif
176 900 format (/" M1QN3 (Version 3.3, October 2009): entry point"/
177 & 5x,"dimension of the problem (n):",i14/
178 & 5x,"absolute precision on x (dxmin):",9x,1pd9.2/
179 & 5x,"expected decrease for f (df1):",11x,1pd9.2/
180 & 5x,"relative precision on g (epsg):",10x,1pd9.2,
181 & " (",a3,"-norm)"/
182 & 5x,"maximal number of iterations (niter):",i6/
183 & 5x,"maximal number of simulations (nsim):",i6/
184 & 5x,"printing level (impres):",15x,i4)
185 c
186 c---- check the arguments
187 c
188 if (n.le.0) then
189 omode=2
190 if (reverse.gt.0) reverse = -1
191 if (impres.ge.1) write (io,'(/a)')
192 & " >>> m1qn3: n should be > 0"
193 return
194 endif
195 if (niter.le.0) then
196 omode=2
197 if (reverse.gt.0) reverse = -1
198 if (impres.ge.1) write (io,'(/a)')
199 & " >>> m1qn3: niter should be > 0"
200 return
201 endif
202 if (nsim.le.0) then
203 omode=2
204 if (reverse.gt.0) reverse = -1
205 if (impres.ge.1) write (io,'(/a)')
206 & " >>> m1qn3: nsim should be > 0"
207 return
208 endif
209 if (dxmin.le.0.d0) then
210 omode=2
211 if (reverse.gt.0) reverse = -1
212 if (impres.ge.1) write (io,'(/a)')
213 & " >>> m1qn3: dxmin should be > 0.d0"
214 return
215 endif
216 c if (epsg.le.0.d0 .or. epsg.gt.1.d0) then
217 if (epsg.le.0.d0) then
218 omode=2
219 if (reverse.gt.0) reverse = -1
220 if (impres.ge.1) write (io,'(/a)')
221 & " >>> m1qn3: epsg should be > 0.d0"
222 return
223 endif
224 if (epsg.ge.1.d0) then
225 omode=1
226 niter=0
227 nsim=0
228 epsg=1.d0
229 if (reverse.gt.0) reverse = -1
230 if (impres.ge.1) write (io,'(/a)')
231 & " >>> m1qn3: epsg is >= 1.d0, no need to make progress"
232 goto 1000
233 endif
234 if ((normtype.ne.'two') .and.
235 & (normtype.ne.'sup') .and.
236 & (normtype.ne.'dfn')) then
237 omode=2
238 if (reverse.gt.0) reverse = -1
239 write (io,'(/a,a,a/)') " >>> m1qn3: unknown norm type '",
240 & normtype, "'"
241 return
242 endif
243 if (impres.lt.0) then
244 omode=2
245 if (reverse.gt.0) reverse = -1
246 write (io,'(/a,i0/)')
247 & " >>> m1qn3: impres should be >= 0 and has the value ",
248 & impres
249 return
250 endif
251 c
252 c---- what method
253 c
254 if (imode(1).eq.0) then
255 if (impres.ge.1) write (io,920)
256 920 format (/" m1qn3: Diagonal Initial Scaling mode")
257 sscale=.false.
258 else
259 if (impres.ge.1) write (io,921)
260 921 format (/" m1qn3: Scalar Initial Scaling mode")
261 sscale=.true.
262 endif
263 c
264 if ((ndz.lt.5*n+1).or.((.not.sscale).and.(ndz.lt.6*n+1))) then
265 omode=2
266 if (reverse.gt.0) reverse = -1
267 if (impres.ge.1) write (io,922)
268 922 format (/" >>> m1qn3: not enough memory allocated")
269 return
270 endif
271 c
272 c---- Compute m
273 c
274 call mupdts (sscale,inmemo,n,m,ndz)
275 c
276 c --- Check the value of m (if (y,s) pairs in core, m will be >= 1)
277 c
278 if (m.lt.1) then
279 omode=2
280 if (reverse.gt.0) reverse = -1
281 if (impres.ge.1) write (io,930)
282 930 format (/" >>> m1qn3: m is set too small in mupdts")
283 return
284 endif
285 c
286 c --- mmemo = number of (y,s) pairs in core memory
287 c
288 mmemo=1
289 if (inmemo) mmemo=m
290 c
291 ntravu=2*(2+mmemo)*n+m
292 if (sscale) ntravu=ntravu-n
293 if (impres.ge.1) write (io,931) ndz,ntravu,m
294 931 format (/5x,"allocated memory (ndz) :",i9/
295 & 5x,"used memory : ",i9/
296 & 5x,"number of updates : ",i9)
297 if (ndz.lt.ntravu) then
298 omode=2
299 if (reverse.gt.0) reverse = -1
300 if (impres.ge.1) write (io,922)
301 return
302 endif
303 c
304 if (impres.ge.1) then
305 if (inmemo) then
306 write (io,932)
307 else
308 write (io,933)
309 endif
310 endif
311 932 format (5x,"(y,s) pairs are stored in core memory")
312 933 format (5x,"(y,s) pairs are stored by the user")
313 c
314 c---- cold start or warm restart ?
315 c check iz: iz(1)=n, iz(2)=(0 if DIS, 1 if SIS),
316 c iz(3)=m, iz(4)=jmin, iz(5)=jmax
317 c
318 if (imode(2).eq.0) then
319 if (impres.ge.1) write (io,940)
320 else
321 if (iz(1).ne.n .or. iz(2).ne.imode(1) .or. iz(3).ne.m .or.
322 & iz(4).lt.1 .or. iz(5).lt.0 .or. iz(4).gt.iz(3) .or.
323 & iz(5).gt.iz(3))
324 & then
325 omode=2
326 if (reverse.gt.0) reverse = -1
327 if (impres.ge.1) then
328 write (io,941)
329 if (iz(1).ne.n) write (io,942)
330 if (iz(2).ne.imode(1)) write (io,943)
331 if (iz(3).ne.m) write (io,944)
332 if (iz(4).lt.1 .or. iz(5).lt.0 .or. iz(4).gt.iz(3)
333 & .or. iz(5).gt.iz(3)) write (io,945)
334 endif
335 return
336 endif
337 if (impres.ge.1) write (io,946)
338 endif
339 940 format (/" m1qn3: cold start"/1x)
340 941 format (/" >>> m1qn3: inconsistent warm restart ")
341 942 format (" >>> m1qn3: (the number of variables has changed)")
342 943 format (" >>> m1qn3: (the scaling mode has changed)")
343 944 format (" >>> m1qn3: (the number of updates has changed)")
344 945 format (" >>> m1qn3: (wrong pointers)")
345 946 format (/" m1qn3: warm restart"/1x)
346 iz(1)=n
347 iz(2)=0
348 if (sscale) iz(2)=1
349 iz(3)=m
350 c
351 c---- split the working zone dz
352 c
353 idiag=1
354 iybar=idiag+n
355 if (sscale) iybar=1
356 isbar=iybar+n*mmemo
357 id=isbar+n*mmemo
358 igg=id+n
359 iaux=igg+n
360 ialpha=iaux+n
361 c
362 c---- call the optimization code
363 c
364 9999 continue
365 call m1qn3a (simul,prosca,ctonb,ctcab,n,x,f,g,dxmin,df1,epsg,
366 & normtype,impres,io,imode,omode,niter,nsim,inmemo,
367 & iz(3),iz(4),iz(5),dz(id),dz(igg),dz(idiag),dz(iaux),
368 & dz(ialpha),dz(iybar),dz(isbar),reverse,reentry,indic,
369 & izs,rzs,dzs)
370 if (reentry.gt.0) return
371 c
372 c---- impressions finales
373 c
374 1000 continue
375 if (impres.ge.1) write (io,960) omode,niter,nsim,epsg
376 960 format (/1x,79("-")/
377 & /" m1qn3: output mode is ",i2
378 & /5x,"number of iterations: ",i14
379 & /5x,"number of simulations: ",i13
380 & /5x,"realized relative precision on g: ",1pd9.2)
381 if (normtype.eq.'two') then
382 gnorm = sqrt(ddot(n,g,1,g,1))
383 elseif (normtype.eq.'sup') then
384 gnorm = dnrmi(n,g)
385 elseif (normtype.eq.'dfn') then
386 call prosca (n,g,g,ps,izs,rzs,dzs)
387 gnorm=dsqrt(ps)
388 endif
389
390 if (impres.ge.1) write (io,961) f,normtype,gnorm
391 961 format (5x,"f = ",1pd15.8
392 & /5x,a3,"-norm of g = ",1pd15.8)
393
394 return
395 end
396 c
397 c--------0---------0---------0---------0---------0---------0---------0--
398 c
399 subroutine m1qn3a (simul,prosca,ctonb,ctcab,n,x,f,g,dxmin,df1,
400 & epsg,normtype,impres,io,imode,omode,niter,nsim,
401 & inmemo,m,jmin,jmax,d,gg,diag,aux,alpha,ybar,
402 & sbar,reverse,reentry,indic,izs,rzs,dzs)
403 c----
404 c
405 c Code d optimisation proprement dit.
406 c
407 c----
408 c
409 implicit none
410 c
411 c arguments
412 c
413 logical inmemo
414 character*3 normtype
415 integer n,impres,io,imode(3),omode,niter,nsim,m,jmin,jmax,indic,
416 & reverse,reentry,izs(*)
417 real rzs(*)
418 double precision x(n),f,g(n),dxmin,df1,epsg,d(n),gg(n),diag(n),
419 & aux(n),alpha(m),ybar(n,1),sbar(n,1),dzs(*)
420 external simul,prosca,ctonb,ctcab
421 c
422 c variables locales
423 c
424 logical skip_update
425 integer i,moderl
426 double precision preco,precos,ys,den,dk1,ps,ps2,hp0
427 #include "m1qn3a_common.h"
428 CML logical sscale,cold,warm,skip_update
429 CML integer i,itmax,moderl,isim,jcour
430 CML double precision d1,t,tmin,tmax,gnorm,gnorms,eps1,ff,preco,precos,
431 CML & ys,den,dk,dk1,ps,ps2,hp0
432 CML save sscale,cold,warm,i,itmax,moderl,isim,jcour,d1,t,tmin,tmax,
433 CML & gnorm,gnorms,eps1,ff,preco,precos,ys,den,dk,dk1,ps,ps2,hp0
434 c
435 c parametres
436 c
437 double precision rm1,rm2
438 parameter (rm1=0.0001d+0,rm2=0.99d+0)
439 double precision pi
440 parameter (pi=3.1415927d+0)
441 double precision rmin
442 c
443 c function
444 c
445 double precision ddot,dnrmi
446 c
447 c --- possible jumps
448 c 9998: call of the simulator in m1qn3a with indic = 1
449 c 9999: call of the simulator in mlis3 with indic = 4
450 c
451 if (reentry.eq.1) goto 9998
452 if (reentry.eq.2) goto 9999
453 c
454 c---- initialisation
455 c
456 rmin=1.d-20
457 c
458 sscale=.true.
459 if (imode(1).eq.0) sscale=.false.
460 c
461 warm=.false.
462 if (imode(2).eq.1) warm=.true.
463 cold=.not.warm
464 c
465 skip_update = .false.
466 c
467 itmax=niter
468 niter=0
469 isim=1
470 eps1=1.d+0
471 c
472 call prosca (n,g,g,ps,izs,rzs,dzs)
473 gnorm = dsqrt(ps)
474 if (normtype.eq.'two') then
475 gnorms = sqrt(ddot(n,g,1,g,1))
476 elseif (normtype.eq.'sup') then
477 gnorms = dnrmi(n,g)
478 elseif (normtype.eq.'dfn') then
479 gnorms = gnorm
480 endif
481 if (impres.ge.1) write (io,900) f,normtype,gnorms
482 900 format (5x,"f = ",1pd15.8
483 & /5x,a3,"-norm of g = ",1pd15.8)
484 if (gnorms.lt.rmin) then
485 omode=2
486 if (impres.ge.1) write (io,901)
487 goto 1000
488 endif
489 901 format (/" >>> m1qn3a: initial gradient is too small")
490 c
491 c --- initialisation pour dd
492 c
493 if (cold) then
494 jmin=1
495 jmax=0
496 endif
497 jcour=1
498 if (inmemo) jcour=jmax
499 c
500 c --- mise a l'echelle de la premiere direction de descente
501 c
502 if (cold) then
503 c
504 c --- use Fletcher's scaling and initialize diag to 1.
505 c
506 precos=2.d+0*df1/gnorm**2
507 do 10 i=1,n
508 d(i)=-g(i)*precos
509 diag(i)=1.d+0
510 10 continue
511 if (impres.ge.5) write(io,902) precos
512 902 format (/" m1qn3a: descent direction -g: precon = ",d10.3)
513 else
514 c
515 c --- use the matrix stored in [diag and] the (y,s) pairs
516 c
517 if (sscale) then
518 call prosca (n,ybar(1,jcour),ybar(1,jcour),ps,izs,rzs,dzs)
519 precos=1.d+0/ps
520 endif
521 do 11 i=1,n
522 d(i)=-g(i)
523 11 continue
524 if (inmemo) then
525 call dd (prosca,ctonb,ctcab,n,sscale,m,d,aux,jmin,jmax,
526 & precos,diag,alpha,ybar,sbar,izs,rzs,dzs)
527 else
528 call dds (prosca,ctonb,ctcab,n,sscale,m,d,aux,jmin,jmax,
529 & precos,diag,alpha,ybar,sbar,izs,rzs,dzs)
530 endif
531 endif
532 c
533 if (impres.eq.3) write(io,903)
534 if (impres.eq.4) write(io,903)
535 903 format (/1x,79("-"))
536 904 format (1x)
537 c
538 c --- initialisation pour mlis3
539 c
540 tmax=1.d+20
541 call prosca (n,d,g,hp0,izs,rzs,dzs)
542 if (hp0.ge.0.d+0) then
543 omode=7
544 if (impres.ge.1) write (io,905) niter,hp0
545 goto 1000
546 endif
547 905 format (/" >>> m1qn3 (iteration ",i2,"): "
548 & /5x," the search direction d is not a ",
549 & "descent direction: (g,d) = ",d12.5)
550 c
551 c --- compute the angle (-g,d)
552 c
553 if (warm.and.impres.ge.5) then
554 call prosca (n,g,g,ps,izs,rzs,dzs)
555 ps=dsqrt(ps)
556 call prosca (n,d,d,ps2,izs,rzs,dzs)
557 ps2=dsqrt(ps2)
558 ps=hp0/ps/ps2
559 ps=dmin1(-ps,1.d+0)
560 ps=dacos(ps)
561 d1=ps*180.d+0/pi
562 write (io,906) sngl(d1)
563 endif
564 906 format (/" m1qn3: descent direction d: ",
565 & "angle(-g,d) = ",f5.1," degrees")
566 c
567 c---- Debut de l iteration. on cherche x(k+1) de la forme x(k) + t*d,
568 c avec t > 0. On connait d.
569 c
570 c Debut de la boucle: etiquette 100,
571 c Sortie de la boucle: goto 1000.
572 c
573 100 niter=niter+1
574 if (impres.ge.5) write(io,903)
575 if (impres.ge.4) write(io,904)
576 if (impres.ge.4) write (io,910) niter,isim,f,hp0
577 910 format (" m1qn3: iter ",i0,", simul ",i0,
578 & ", f=",1pd15.8,", h'(0)=",d12.5)
579 c
580 c --- free simulation if desired
581 c
582 if (imode(3).gt.0) then
583 if (mod(niter-1,imode(3)).eq.0) then
584 indic=1
585 if (reverse.gt.0) then
586 reentry = 1
587 return
588 else
589 call simul(indic,n,x,f,g,izs,rzs,dzs)
590 endif
591 endif
592 endif
593 9998 continue
594 c
595 c --- recherche lineaire et nouveau point x(k+1)
596 c
597 do 101 i=1,n
598 gg(i)=g(i)
599 101 continue
600 ff=f
601 if (impres.ge.5) write (io,911)
602 911 format (/" m1qn3: line search")
603 c
604 c --- calcul de tmin
605 c
606 tmin=0.d+0
607 do 200 i=1,n
608 tmin=dmax1(tmin,dabs(d(i)))
609 200 continue
610 tmin=dxmin/tmin
611 t=1.d+0
612 d1=hp0
613 c
614 9999 continue
615 call mlis3 (n,simul,prosca,x,f,d1,t,tmin,tmax,d,g,rm2,rm1,impres,
616 & io,moderl,isim,nsim,aux,reverse,reentry,indic,izs,rzs,
617 & dzs)
618 if (reentry.gt.0) return
619 c
620 c --- mlis3 renvoie les nouvelles valeurs de x, f et g
621 c
622 if (moderl.ne.0) then
623 if (moderl.lt.0) then
624 c
625 c --- calcul impossible
626 c t, g: ou les calculs sont impossibles
627 c x, f: ceux du t_gauche (donc f <= ff)
628 c
629 omode=moderl
630 goto 1000
631 elseif (moderl.eq.1) then
632 c
633 c --- descente bloquee sur tmax
634 c
635 skip_update = .true.
636 c omode=3
637 c if (impres.ge.1) write(io,912) niter
638 c 912 format (/" >>> m1qn3 (iteration ",i0,
639 c & "): line search blocked on tmax: "/
640 c & " >>> possible reasons: bad scaling,",
641 c & " unbounded problem")
642 elseif (moderl.eq.4) then
643 c
644 c --- nsim atteint
645 c x, f: ceux du t_gauche (donc f <= ff)
646 c
647 omode=5
648 goto 1000
649 elseif (moderl.eq.5) then
650 c
651 c --- arret demande par l utilisateur (indic = 0)
652 c x, f: ceux en sortie du simulateur
653 c
654 omode=0
655 goto 1000
656 elseif (moderl.eq.6) then
657 c
658 c --- arret sur dxmin ou appel incoherent
659 c x, f: ceux du t_gauche (donc f <= ff)
660 c
661 omode=6
662 goto 1000
663 endif
664 else
665 skip_update = .false.
666 endif
667 c
668 c NOTE: stopping tests are now done after having updated the matrix, so
669 c that update information can be stored in case of a later warm restart
670 c
671 c --- mise a jour de la matrice
672 c
673 if (skip_update) then
674 if (impres.ge.5) write(io,'(/a)')
675 & " m1qn3: matrix update is skipped"
676 elseif (m.gt.0) then
677 c
678 c --- mise a jour des pointeurs
679 c
680 jmax=jmax+1
681 if (jmax.gt.m) jmax=jmax-m
682 if ((cold.and.niter.gt.m).or.(warm.and.jmin.eq.jmax)) then
683 jmin=jmin+1
684 if (jmin.gt.m) jmin=jmin-m
685 endif
686 if (inmemo) jcour=jmax
687 c
688 c --- y, s et (y,s)
689 c
690 do 400 i=1,n
691 sbar(i,jcour)=t*d(i)
692 ybar(i,jcour)=g(i)-gg(i)
693 400 continue
694 if (impres.ge.5) then
695 call prosca (n,sbar(1,jcour),sbar(1,jcour),ps,izs,rzs,dzs)
696 dk1=dsqrt(ps)
697 if (niter.gt.1) write (io,930) dk1/dk
698 930 format (/" m1qn3: convergence rate, s(k)/s(k-1) = ",
699 & 1pd12.5)
700 dk=dk1
701 endif
702 call prosca (n,ybar(1,jcour),sbar(1,jcour),ys,izs,rzs,dzs)
703 if (ys.le.0.d+0) then
704 omode=7
705 if (impres.ge.1) write (io,931) niter,ys
706 931 format (/" >>> m1qn3 (iteration ",i2,
707 & "): the scalar product (y,s) = ",d12.5
708 & /27x,"is not positive")
709 goto 1000
710 endif
711 c
712 c --- ybar et sbar
713 c
714 d1=dsqrt(1.d+0/ys)
715 do 410 i=1,n
716 sbar(i,jcour)=d1*sbar(i,jcour)
717 ybar(i,jcour)=d1*ybar(i,jcour)
718 410 continue
719 if (.not.inmemo) call ystbl (.true.,ybar,sbar,n,jmax)
720 c
721 c --- compute the scalar or diagonal preconditioner
722 c
723 if (impres.ge.5) write(io,932)
724 932 format (/" m1qn3: matrix update:")
725 c
726 c --- Here is the Oren-Spedicato factor, for scalar scaling
727 c
728 if (sscale) then
729 call prosca (n,ybar(1,jcour),ybar(1,jcour),ps,izs,rzs,dzs)
730 precos=1.d+0/ps
731 c
732 if (impres.ge.5) write (io,933) precos
733 933 format (5x,"Oren-Spedicato factor = ",d10.3)
734 c
735 c --- Scale the diagonal to Rayleigh s ellipsoid.
736 c Initially (niter.eq.1) and for a cold start, this is
737 c equivalent to an Oren-Spedicato scaling of the
738 c identity matrix.
739 c
740 else
741 call ctonb (n,ybar(1,jcour),aux,izs,rzs,dzs)
742 ps=0.d0
743 do 420 i=1,n
744 ps=ps+diag(i)*aux(i)*aux(i)
745 420 continue
746 d1=1.d0/ps
747 if (impres.ge.5) then
748 write (io,934) d1
749 934 format(5x,"fitting the ellipsoid: factor = ",1pd10.3)
750 endif
751 do 421 i=1,n
752 diag(i)=diag(i)*d1
753 421 continue
754 c
755 c --- update the diagonal
756 c (gg is used as an auxiliary vector)
757 c
758 call ctonb (n,sbar(1,jcour),gg,izs,rzs,dzs)
759 ps=0.d0
760 do 430 i=1,n
761 ps=ps+gg(i)*gg(i)/diag(i)
762 430 continue
763 den=ps
764 do 431 i=1,n
765 diag(i)=1.d0/
766 & (1.d0/diag(i)+aux(i)**2-(gg(i)/diag(i))**2/den)
767 if (diag(i).le.0.d0) then
768 if (impres.ge.5) write (io,935) i,diag(i),rmin
769 diag(i)=rmin
770 endif
771 431 continue
772 935 format (/" >>> m1qn3-WARNING: diagonal element ",i8,
773 & " is negative (",d10.3,"), reset to ",d10.3)
774 c
775 if (impres.ge.5) then
776 ps=0.d0
777 do 440 i=1,n
778 ps=ps+diag(i)
779 440 continue
780 ps=ps/n
781 preco=ps
782 c
783 ps2=0.d0
784 do 441 i=1,n
785 ps2=ps2+(diag(i)-ps)**2
786 441 continue
787 ps2=dsqrt(ps2/n)
788 write (io,936) preco,ps2
789 936 format (5x,"updated diagonal: average value = ",
790 & 1pd10.3,", sqrt(variance) = ",d10.3)
791 endif
792 endif
793 endif
794 c
795 c --- printings
796 c
797 c
798 c --- tests d arret
799 c
800 call prosca(n,g,g,ps,izs,rzs,dzs)
801 if (normtype.eq.'two') then
802 gnorm = sqrt(ddot(n,g,1,g,1))
803 elseif (normtype.eq.'sup') then
804 gnorm = dnrmi(n,g)
805 elseif (normtype.eq.'dfn') then
806 gnorm = dsqrt(ps)
807 endif
808 eps1 = gnorm/gnorms
809 c
810 if (impres.eq.3) then
811 if (mod(niter-1,50).eq.0) write(io,'(/a,a)')
812 & " iter simul stepsize f |g|",
813 & " |g|/|g0|"
814 write(io,
815 & '(1x,i5,2x,i5,2x,1pd8.2,2x,d21.14,2x,d11.5,2x,d10.4)')
816 & niter, isim, t, f, gnorm, eps1
817 endif
818 if (impres.ge.5) write (io,940) eps1
819 940 format (/" m1qn3: stopping criterion on g: ",1pd12.5)
820 if (eps1.lt.epsg) then
821 omode=1
822 goto 1000
823 endif
824 if (niter.eq.itmax) then
825 omode=4
826 if (impres.ge.1) write (io,941) niter
827 941 format (/" >>> m1qn3 (iteration ",i0,
828 & "): maximal number of iterations")
829 goto 1000
830 endif
831 if (isim.gt.nsim) then
832 omode=5
833 if (impres.ge.1) write (io,942) niter,isim
834 942 format (/" >>> m1qn3 (iteration ",i3,"): ",i6,
835 & " simulations (maximal number reached)")
836 goto 1000
837 endif
838 c
839 c --- calcul de la nouvelle direction de descente d = - H.g
840 c
841 if (m.eq.0) then
842 preco=2.d0*(ff-f)/ps
843 do 500 i=1,n
844 d(i)=-g(i)*preco
845 500 continue
846 else
847 do 510 i=1,n
848 d(i)=-g(i)
849 510 continue
850 if (inmemo) then
851 call dd (prosca,ctonb,ctcab,n,sscale,m,d,aux,jmin,jmax,
852 & precos,diag,alpha,ybar,sbar,izs,rzs,dzs)
853 else
854 call dds (prosca,ctonb,ctcab,n,sscale,m,d,aux,jmin,jmax,
855 & precos,diag,alpha,ybar,sbar,izs,rzs,dzs)
856 endif
857 endif
858 c
859 c --- test: la direction d est-elle de descente ?
860 c hp0 sera utilise par mlis3
861 c
862 call prosca (n,d,g,hp0,izs,rzs,dzs)
863 if (hp0.ge.0.d+0) then
864 omode=7
865 if (impres.ge.1) write (io,905) niter,hp0
866 goto 1000
867 endif
868 if (impres.ge.5) then
869 call prosca (n,g,g,ps,izs,rzs,dzs)
870 ps=dsqrt(ps)
871 call prosca (n,d,d,ps2,izs,rzs,dzs)
872 ps2=dsqrt(ps2)
873 ps=hp0/ps/ps2
874 ps=dmin1(-ps,1.d+0)
875 ps=dacos(ps)
876 d1=ps
877 d1=d1*180.d0/pi
878 write (io,906) sngl(d1)
879 endif
880 c
881 c---- on poursuit les iterations
882 c
883 goto 100
884 c
885 c --- n1qn3 has finished for ever
886 c
887 1000 continue
888 if (reverse.ne.0) reverse = -1
889 reentry = 0
890 nsim=isim
891 epsg=eps1
892 return
893 end
894 c
895 c--------0---------0---------0---------0---------0---------0---------0--
896 c
897 subroutine dd (prosca,ctonb,ctcab,n,sscale,nm,depl,aux,jmin,jmax,
898 & precos,diag,alpha,ybar,sbar,izs,rzs,dzs)
899 c----
900 c
901 c calcule le produit H.g ou
902 c . H est une matrice construite par la formule de bfgs inverse
903 c a nm memoires a partir de la matrice diagonale diag
904 c dans un espace hilbertien dont le produit scalaire
905 c est donne par prosca
906 c (cf. J. Nocedal, Math. of Comp. 35/151 (1980) 773-782)
907 c . g est un vecteur de dimension n (en general le gradient)
908 c
909 c la matrice diag apparait donc comme un preconditionneur diagonal
910 c
911 c depl = g (en entree), = H g (en sortie)
912 c
913 c la matrice H est memorisee par les vecteurs des tableaux
914 c ybar, sbar et les pointeurs jmin, jmax
915 c
916 c alpha(nm) est une zone de travail
917 c
918 c izs(*),rzs(*),dzs(*) sont des zones de travail pour prosca
919 c
920 c----
921 c
922 implicit none
923 c
924 c arguments
925 c
926 logical sscale
927 integer n,nm,jmin,jmax,izs(*)
928 real rzs(*)
929 double precision depl(n),precos,diag(n),alpha(nm),ybar(n,1),
930 & sbar(n,1),aux(n),dzs(*)
931 external prosca,ctonb,ctcab
932 c
933 c variables locales
934 c
935 integer jfin,i,j,jp
936 double precision r,ps
937 c
938 jfin=jmax
939 if (jfin.lt.jmin) jfin=jmax+nm
940 c
941 c phase de descente
942 c
943 do 100 j=jfin,jmin,-1
944 jp=j
945 if (jp.gt.nm) jp=jp-nm
946 call prosca (n,depl,sbar(1,jp),ps,izs,rzs,dzs)
947 r=ps
948 alpha(jp)=r
949 do 20 i=1,n
950 depl(i)=depl(i)-r*ybar(i,jp)
951 20 continue
952 100 continue
953 c
954 c preconditionnement
955 c
956 if (sscale) then
957 do 150 i=1,n
958 depl(i)=depl(i)*precos
959 150 continue
960 else
961 call ctonb (n,depl,aux,izs,rzs,dzs)
962 do 151 i=1,n
963 aux(i)=aux(i)*diag(i)
964 151 continue
965 call ctcab (n,aux,depl,izs,rzs,dzs)
966 endif
967 c
968 c remontee
969 c
970 do 200 j=jmin,jfin
971 jp=j
972 if (jp.gt.nm) jp=jp-nm
973 call prosca (n,depl,ybar(1,jp),ps,izs,rzs,dzs)
974 r=alpha(jp)-ps
975 do 120 i=1,n
976 depl(i)=depl(i)+r*sbar(i,jp)
977 120 continue
978 200 continue
979 return
980 end
981 c
982 c--------0---------0---------0---------0---------0---------0---------0--
983 c
984 subroutine dds (prosca,ctonb,ctcab,n,sscale,nm,depl,aux,jmin,
985 & jmax,precos,diag,alpha,ybar,sbar,izs,rzs,dzs)
986 c----
987 c
988 c This subroutine has the same role as dd (computation of the
989 c product H.g). It supposes however that the (y,s) pairs are not
990 c stored in core memory, but on a devise chosen by the user.
991 c The access to this devise is performed via the subroutine ystbl.
992 c
993 c----
994 c
995 implicit none
996 c
997 c arguments
998 c
999 logical sscale
1000 integer n,nm,jmin,jmax,izs(*)
1001 real rzs(*)
1002 double precision depl(n),precos,diag(n),alpha(nm),ybar(n),sbar(n),
1003 & aux(n),dzs(*)
1004 external prosca,ctonb,ctcab
1005 c
1006 c variables locales
1007 c
1008 integer jfin,i,j,jp
1009 double precision r,ps
1010 c
1011 jfin=jmax
1012 if (jfin.lt.jmin) jfin=jmax+nm
1013 c
1014 c phase de descente
1015 c
1016 do 100 j=jfin,jmin,-1
1017 jp=j
1018 if (jp.gt.nm) jp=jp-nm
1019 call ystbl (.false.,ybar,sbar,n,jp)
1020 call prosca (n,depl,sbar,ps,izs,rzs,dzs)
1021 r=ps
1022 alpha(jp)=r
1023 do 20 i=1,n
1024 depl(i)=depl(i)-r*ybar(i)
1025 20 continue
1026 100 continue
1027 c
1028 c preconditionnement
1029 c
1030 if (sscale) then
1031 do 150 i=1,n
1032 depl(i)=depl(i)*precos
1033 150 continue
1034 else
1035 call ctonb (n,depl,aux,izs,rzs,dzs)
1036 do 151 i=1,n
1037 aux(i)=aux(i)*diag(i)
1038 151 continue
1039 call ctcab (n,aux,depl,izs,rzs,dzs)
1040 endif
1041 c
1042 c remontee
1043 c
1044 do 200 j=jmin,jfin
1045 jp=j
1046 if (jp.gt.nm) jp=jp-nm
1047 call ystbl (.false.,ybar,sbar,n,jp)
1048 call prosca (n,depl,ybar(1),ps,izs,rzs,dzs)
1049 r=alpha(jp)-ps
1050 do 120 i=1,n
1051 depl(i)=depl(i)+r*sbar(i)
1052 120 continue
1053 200 continue
1054 return
1055 end
1056 c
1057 c--------0---------0---------0---------0---------0---------0---------0--
1058 c
1059 subroutine mlis3 (n,simul,prosca,x,f,fpn,t,tmin,tmax,d,g,
1060 1 amd,amf,imp,io,logic,nap,napmax,xn,
1061 1 reverse,reentry,indic,izs,rzs,dzs)
1062 c
1063 c ----
1064 c
1065 c mlis3 + minuscules + commentaires
1066 c + version amelioree (XII 88): interpolation cubique systematique
1067 c et anti-overflows
1068 c + declaration variables (II/89, JCG).
1069 c + barr is also progressively decreased (12/93, CL & JChG).
1070 c barmul is set to 5.
1071 c
1072 c ----------------------------------------------------------------
1073 c
1074 c en sortie logic =
1075 c
1076 c 0 descente serieuse
1077 c 1 descente bloquee sur tmax
1078 c 4 nap > napmax
1079 c 5 arret demande par le simulateur
1080 c 6 fonction et gradient pas d accord
1081 c < 0 contrainte implicite active
1082 c
1083 c indic on entry (in reverse communication)
1084 c
1085 c < 0: the simulator cannot compute f and g
1086 c = 0: the simulator wants to stop
1087 c > 0: the simulator has done its work
1088 c
1089 c indic on return (in reverse communication)
1090 c
1091 c = 4: the simulator is asked to compute f and g
1092 c
1093 c reverse
1094 c
1095 c = 0: direct communication
1096 c = 1: reverse communication
1097 c
1098 c reentry on entry
1099 c
1100 c = 2: reverse communication, return from a simulation, skip to
1101 c the place where the simulator was called (9999)
1102 c
1103 c reentry return
1104 c
1105 c = 0: reverse communication, the linesearch has finished its job
1106 c = 2: reverse communication, a simulation is required
1107 c
1108 c ----
1109 c
1110 implicit none
1111 c
1112 c --- arguments
1113 c
1114 integer n,imp,io,logic,nap,napmax,indic,reverse,reentry,izs(*)
1115 real rzs(*)
1116 double precision x(n),f,fpn,t,tmin,tmax,d(n),g(n),amd,amf,xn(n),
1117 & dzs(*)
1118 external simul,prosca
1119 c
1120 c --- variables locales
1121 c
1122 #include "mlis3_common.h"
1123 logical t_increased
1124 integer i,indica,indicd
1125 double precision tg,fg,fpg,td,ta,d2,fp,ffn,fd,
1126 & fpd,z,test,gauche,droite,taa,ps
1127 CML save t_increased,i,indica,indicd,tesf,tesd,tg,fg,fpg,td,ta,fa,
1128 CML & fpa,d2,fn,fp,ffn,fd,fpd,z,test,barmin,barmul,barmax,barr,
1129 CML & gauche,droite,taa,ps
1130 c
1131 1000 format (/4x," mlis3 ",4x,"fpn=",1pd10.3," d2=",d9.2,
1132 1 " tmin=",d9.2," tmax=",d9.2)
1133 1001 format (/4x," mlis3",3x,"stop on tmin",5x,
1134 1 "stepsizes",11x,"functions",8x,"derivatives")
1135 1002 format (4x," mlis3",37x,1pd10.3,2d11.3)
1136 1003 format (4x," mlis3",1pd14.3,2d11.3)
1137 1004 format (4x," mlis3",37x,1pd10.3," indic=",i3)
1138 1005 format (4x," mlis3",14x,1pd18.8,2x,d21.14,2x,d11.4)
1139 1006 format (4x," mlis3",14x,1pd18.8," indic=",i3)
1140 1007 format (/4x," mlis3",10x,"tmin forced to tmax")
1141 1008 format (/4x," mlis3",10x,"inconsistent call")
1142 c
1143 c --- possible jump
1144 c
1145 if (reentry.eq.2) goto 9999
1146 c
1147 if (n.gt.0 .and. fpn.lt.0.d0 .and. t.gt.0.d0
1148 1 .and. tmax.gt.0.d0 .and. amf.gt.0.d0
1149 1 .and. amd.gt.amf .and. amd.lt.1.d0) go to 5
1150 logic=6
1151 go to 999
1152 5 tesf=amf*fpn
1153 tesd=amd*fpn
1154 barmin=0.01d0
1155 barmul=5.d0
1156 barmax=0.3d0
1157 barr=barmin
1158 td=0.d0
1159 tg=0.d0
1160 fn=f
1161 fg=fn
1162 fpg=fpn
1163 ta=0.d0
1164 fa=fn
1165 fpa=fpn
1166 call prosca (n,d,d,ps,izs,rzs,dzs)
1167 d2=ps
1168 c
1169 c elimination d un t initial ridiculement petit
1170 c
1171 c<
1172 c if (t.gt.tmin) go to 20
1173 c t=tmin
1174 c if (t.le.tmax) go to 20
1175 c if (imp.gt.0) write (io,1007)
1176 c tmin=tmax
1177 c 20 if (fn+t*fpn.lt.fn+0.9d0*t*fpn) go to 30
1178 c t=2.d0*t
1179 c go to 20
1180 c changed into
1181 if (t.lt.tmin) then
1182 t=tmin
1183 if (imp.ge.4) write (io,'(a)')
1184 & ' mlis3: initial step-size increased to tmin'
1185 if (t.gt.tmax) then
1186 if (imp.gt.0) write (io,1007)
1187 tmin=tmax
1188 endif
1189 endif
1190 t_increased = .false.
1191 do while (fn+t*fpn.ge.fn+0.9d0*t*fpn)
1192 t_increased = .true.
1193 t = 2.d0*t
1194 enddo
1195 if (t_increased .and. (imp.ge.4)) write (io,'(a,1pd10.3)')
1196 & ' mlis3: initial step-size increased to ',t
1197 c>
1198 30 indica=1
1199 logic=0
1200 if (t.gt.tmax) then
1201 t=tmax
1202 logic=1
1203 endif
1204 if (imp.ge.4) write (io,1000) fpn,d2,tmin,tmax
1205 c
1206 c --- nouveau x
1207 c initialize xn to the current iterate
1208 c use x as the trial iterate
1209 c
1210 do i=1,n
1211 xn(i)=x(i)
1212 x(i)=xn(i)+t*d(i)
1213 enddo
1214 c
1215 c --- boucle
1216 c
1217 100 nap=nap+1
1218 if(nap.gt.napmax) then
1219 logic=4
1220 fn=fg
1221 do i=1,n
1222 xn(i)=xn(i)+tg*d(i)
1223 enddo
1224 go to 999
1225 endif
1226 c
1227 c --- appel simulateur
1228 c
1229 indic=4
1230 if (reverse.gt.0) then
1231 reentry = 2
1232 return
1233 else
1234 call simul(indic,n,x,f,g,izs,rzs,dzs)
1235 endif
1236 9999 continue
1237 if (indic.eq.0) then
1238 c
1239 c --- arret demande par l utilisateur
1240 c
1241 logic=5
1242 fn=f
1243 do 170 i=1,n
1244 xn(i)=x(i)
1245 170 continue
1246 go to 999
1247 endif
1248 if(indic.lt.0) then
1249 c
1250 c --- les calculs n ont pas pu etre effectues par le simulateur
1251 c
1252 td=t
1253 indicd=indic
1254 logic=0
1255 if (imp.ge.4) write (io,1004) t,indic
1256 t=tg+0.1d0*(td-tg)
1257 go to 905
1258 endif
1259 c
1260 c --- les tests elementaires sont faits, on y va
1261 c
1262 call prosca (n,d,g,ps,izs,rzs,dzs)
1263 fp=ps
1264 c
1265 c --- premier test de Wolfe
1266 c
1267 ffn=f-fn
1268 if(ffn.gt.t*tesf) then
1269 td=t
1270 fd=f
1271 fpd=fp
1272 indicd=indic
1273 logic=0
1274 if(imp.ge.4) write (io,1002) t,ffn,fp
1275 go to 500
1276 endif
1277 c
1278 c --- test 1 ok, donc deuxieme test de Wolfe
1279 c
1280 if(imp.ge.4) write (io,1003) t,ffn,fp
1281 if(fp.gt.tesd) then
1282 logic=0
1283 go to 320
1284 endif
1285 if (logic.eq.0) go to 350
1286 c
1287 c --- test 2 ok, donc pas serieux, on sort
1288 c
1289 320 fn=f
1290 do 330 i=1,n
1291 xn(i)=x(i)
1292 330 continue
1293 go to 999
1294 c
1295 c
1296 c
1297 350 tg=t
1298 fg=f
1299 fpg=fp
1300 if(td.ne.0.d0) go to 500
1301 c
1302 c extrapolation
1303 c
1304 taa=t
1305 gauche=(1.d0+barmin)*t
1306 droite=10.d0*t
1307 call ecube (t,f,fp,ta,fa,fpa,gauche,droite)
1308 ta=taa
1309 if(t.lt.tmax) go to 900
1310 logic=1
1311 t=tmax
1312 go to 900
1313 c
1314 c interpolation
1315 c
1316 500 if(indica.le.0) then
1317 ta=t
1318 t=0.9d0*tg+0.1d0*td
1319 go to 900
1320 endif
1321 test=barr*(td-tg)
1322 gauche=tg+test
1323 droite=td-test
1324 taa=t
1325 call ecube (t,f,fp,ta,fa,fpa,gauche,droite)
1326 ta=taa
1327 if (t.gt.gauche .and. t.lt.droite) then
1328 barr=dmax1(barmin,barr/barmul)
1329 c barr=barmin
1330 else
1331 barr=dmin1(barmul*barr,barmax)
1332 endif
1333 c
1334 c --- fin de boucle
1335 c - t peut etre bloque sur tmax
1336 c (venant de l extrapolation avec logic=1)
1337 c
1338 900 fa=f
1339 fpa=fp
1340 905 indica=indic
1341 c
1342 c --- faut-il continuer ?
1343 c
1344 if (td.eq.0.d0) go to 950
1345 if (td-tg.lt.tmin) go to 920
1346 c
1347 c --- limite de precision machine (arret de secours) ?
1348 c
1349 do 910 i=1,n
1350 z=xn(i)+t*d(i)
1351 if (z.ne.xn(i).and.z.ne.x(i)) go to 950
1352 910 continue
1353 if (imp.gt.3) write (io,'(5x,a)') "mlis3 no change in x"
1354 c
1355 c --- arret sur dxmin ou de secours
1356 c
1357 920 logic=6
1358 c
1359 c si indicd<0, derniers calculs non faits par simul
1360 c
1361 if (indicd.lt.0) logic=indicd
1362 c
1363 c si tg=0, xn = xn_depart,
1364 c sinon on prend xn=x_gauche qui fait decroitre f
1365 c
1366 if (tg.eq.0.d0) go to 940
1367 fn=fg
1368 do 930 i=1,n
1369 930 xn(i)=xn(i)+tg*d(i)
1370 940 if (imp.le.3) go to 999
1371 write (io,1001)
1372 write (io,1005) tg,fg,fpg
1373 if (logic.eq.6) write (io,1005) td,fd,fpd
1374 if (logic.eq.7) write (io,1006) td,indicd
1375 go to 999
1376 c
1377 c recopiage de x et boucle
1378 c
1379 950 do 960 i=1,n
1380 960 x(i)=xn(i)+t*d(i)
1381 go to 100
1382 c --- linesearch finished, no skip at next entry in mlis3
1383 999 if (reverse.ne.0) reentry = 0
1384 return
1385 end
1386 c
1387 c
1388 c--------0---------0---------0---------0---------0---------0---------0--
1389 c
1390 subroutine ecube (t,f,fp,ta,fa,fpa,tlower,tupper)
1391 c
1392 implicit none
1393 c
1394 c --- arguments
1395 c
1396 double precision sign,den,anum,t,f,fp,ta,fa,fpa,tlower,tupper
1397 c
1398 c --- variables locales
1399 c
1400 double precision z1,b,discri
1401 c
1402 c Using f and fp at t and ta, computes new t by cubic formula
1403 c safeguarded inside [tlower,tupper].
1404 c
1405 z1=fp+fpa-3.d0*(fa-f)/(ta-t)
1406 b=z1+fp
1407 c
1408 c first compute the discriminant (without overflow)
1409 c
1410 if (dabs(z1).le.1.d0) then
1411 discri=z1*z1-fp*fpa
1412 else
1413 discri=fp/z1
1414 discri=discri*fpa
1415 discri=z1-discri
1416 if (z1.ge.0.d0 .and. discri.ge.0.d0) then
1417 discri=dsqrt(z1)*dsqrt(discri)
1418 go to 120
1419 endif
1420 if (z1.le.0.d0 .and. discri.le.0.d0) then
1421 discri=dsqrt(-z1)*dsqrt(-discri)
1422 go to 120
1423 endif
1424 discri=-1.d0
1425 endif
1426 if (discri.lt.0.d0) then
1427 if (fp.lt.0.d0) t=tupper
1428 if (fp.ge.0.d0) t=tlower
1429 go to 900
1430 endif
1431 c
1432 c discriminant nonnegative, compute solution (without overflow)
1433 c
1434 discri=dsqrt(discri)
1435 120 if (t-ta.lt.0.d0) discri=-discri
1436 sign=(t-ta)/dabs(t-ta)
1437 if (b*sign.gt.0.d+0) then
1438 t=t+fp*(ta-t)/(b+discri)
1439 else
1440 den=z1+b+fpa
1441 anum=b-discri
1442 if (dabs((t-ta)*anum).lt.(tupper-tlower)*dabs(den)) then
1443 t=t+anum*(ta-t)/den
1444 else
1445 t=tupper
1446 endif
1447 endif
1448 900 t=dmax1(t,tlower)
1449 t=dmin1(t,tupper)
1450 return
1451 end
1452 c
1453 c--------0---------0---------0---------0---------0---------0---------0--
1454 c
1455 subroutine mupdts (sscale,inmemo,n,m,nrz)
1456 c
1457 implicit none
1458 c
1459 c arguments
1460 c
1461 logical sscale,inmemo
1462 integer n,m,nrz
1463 c----
1464 c
1465 c On entry:
1466 c sscale: .true. if scalar initial scaling,
1467 c .false. if diagonal initial scaling
1468 c n: number of variables
1469 c
1470 c This routine has to return:
1471 c m: the number of updates to form the approximate Hessien H,
1472 c inmemo: .true., if the vectors y and s used to form H are stored
1473 c in core memory,
1474 c .false. otherwise (storage of y and s on disk, for
1475 c instance).
1476 c When inmemo=.false., the routine "ystbl", which stores and
1477 c restores (y,s) pairs, has to be rewritten.
1478 c
1479 c----
1480 c
1481 if (sscale) then
1482 m=(nrz-3*n)/(2*n+1)
1483 else
1484 m=(nrz-4*n)/(2*n+1)
1485 endif
1486 inmemo=.true.
1487 return
1488 end
1489 c
1490 c--------0---------0---------0---------0---------0---------0---------0--
1491 c
1492 subroutine ystbl (store,ybar,sbar,n,j)
1493 c----
1494 c
1495 c This subroutine should store (if store = .true.) or restore
1496 c (if store = .false.) a pair (ybar,sbar) at or from position
1497 c j in memory. Be sure to have 1 <= j <= m, where m in the number
1498 c of updates specified by subroutine mupdts.
1499 c
1500 c The subroutine is used only when the (y,s) pairs are not
1501 c stored in core memory in the arrays ybar(.,.) and sbar(.,.).
1502 c In this case, the subroutine has to be written by the user.
1503 c
1504 c----
1505 c
1506 implicit none
1507 c
1508 c arguments
1509 c
1510 logical store
1511 integer n,j
1512 double precision ybar(n),sbar(n)
1513 c
1514 return
1515 end
1516 c
1517 c--------0---------0---------0---------0---------0---------0---------0--
1518 c
1519 subroutine ctonbe (n,u,v,izs,rzs,dzs)
1520 c
1521 implicit none
1522 integer n,izs(*)
1523 real rzs(*)
1524 double precision u(1),v(1),dzs(*)
1525 c
1526 integer i
1527 c
1528 do i=1,n
1529 v(i)=u(i)
1530 enddo
1531 return
1532 end
1533 c
1534 c--------0---------0---------0---------0---------0---------0---------0--
1535 c
1536 subroutine ctcabe (n,u,v,izs,rzs,dzs)
1537 c
1538 implicit none
1539 integer n,izs(*)
1540 real rzs(*)
1541 double precision u(1),v(1),dzs(*)
1542 c
1543 integer i
1544 c
1545 do i=1,n
1546 v(i)=u(i)
1547 enddo
1548 return
1549 end
1550 c
1551 c--------0---------0---------0---------0---------0---------0---------0--
1552 c
1553 subroutine euclid (n,x,y,ps,izs,rzs,dzs)
1554 c
1555 implicit none
1556 integer n,izs(*)
1557 real rzs(*)
1558 double precision x(n),y(n),ps,dzs(*)
1559 c
1560 integer i
1561 c
1562 ps=0.d0
1563 do i=1,n
1564 ps=ps+x(i)*y(i)
1565 enddo
1566 return
1567 end
1568 c
1569 c--------0---------0---------0---------0---------0---------0---------0--
1570 c
1571 subroutine simul_rc (indic,n,x,f,g,izs,rzs,dzs)
1572 c
1573 implicit none
1574 integer indic,n,izs(*)
1575 real rzs(*)
1576 double precision x(n),f,g(n),dzs(*)
1577 c
1578 return
1579 end
1580 c
1581 c--------0---------0---------0---------0---------0---------0---------0--
1582 c
1583 double precision function dnrmi (n,v)
1584 c
1585 integer n
1586 double precision v(n)
1587 c
1588 c----
1589 c
1590 c Commputes the infinty-norm of the vector v(n)
1591 c
1592 c----
1593 c
1594 c --- local variables
1595 c
1596 integer i
1597 double precision norm
1598 c
1599 c --- compute
1600 c
1601 norm = 0.d0
1602 if (n.gt.0) then
1603 do i=1,n
1604 norm = max(norm,abs(v(i)))
1605 enddo
1606 endif
1607 dnrmi = norm
1608 c
1609 return
1610 end

  ViewVC Help
Powered by ViewVC 1.1.22