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

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

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


Revision 1.4 - (hide annotations) (download)
Wed May 2 10:10:57 2012 UTC (11 years, 11 months ago) by mlosch
Branch: MAIN
Changes since 1.3: +1 -2 lines
store more variables of m1qn3a in m1qn3a_common.h

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

  ViewVC Help
Powered by ViewVC 1.1.22