/[MITgcm]/MITgcm_contrib/mlosch/optim_m1qn3/testbed/m1qn3.F
ViewVC logotype

Contents of /MITgcm_contrib/mlosch/optim_m1qn3/testbed/m1qn3.F

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


Revision 1.1 - (show annotations) (download)
Wed May 2 12:27:42 2012 UTC (12 years ago) by mlosch
Branch: MAIN
CVS Tags: HEAD
simple testing environment for optim_m1qn3

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

  ViewVC Help
Powered by ViewVC 1.1.22