forked from afbarnard/go-lbfgsb
-
Notifications
You must be signed in to change notification settings - Fork 3
/
Copy pathlbfgsb.f
3952 lines (3536 loc) · 128 KB
/
lbfgsb.f
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
c
c L-BFGS-B is released under the “New BSD License” (aka “Modified BSD License”
c or “3-clause license”)
c Please read attached file License.txt
c
c=========== L-BFGS-B (version 3.0. April 25, 2011 ===================
c
c This is a modified version of L-BFGS-B. Minor changes in the updated
c code appear preceded by a line comment as follows
c
c c-jlm-jn
c
c Major changes are described in the accompanying paper:
c
c Jorge Nocedal and Jose Luis Morales, Remark on "Algorithm 778:
c L-BFGS-B: Fortran Subroutines for Large-Scale Bound Constrained
c Optimization" (2011). To appear in ACM Transactions on
c Mathematical Software,
c
c The paper describes an improvement and a correction to Algorithm 778.
c It is shown that the performance of the algorithm can be improved
c significantly by making a relatively simple modication to the subspace
c minimization phase. The correction concerns an error caused by the use
c of routine dpmeps to estimate machine precision.
c
c The total work space **wa** required by the new version is
c
c 2*m*n + 11m*m + 5*n + 8*m
c
c the old version required
c
c 2*m*n + 12m*m + 4*n + 12*m
c
c
c J. Nocedal Department of Electrical Engineering and
c Computer Science.
c Northwestern University. Evanston, IL. USA
c
c
c J.L Morales Departamento de Matematicas,
c Instituto Tecnologico Autonomo de Mexico
c Mexico D.F. Mexico.
c
c March 2011
c
c=============================================================================
subroutine setulb(n, m, x, l, u, nbd, f, g, factr, pgtol, wa, iwa,
+ task, iprint, csave, lsave, isave, dsave)
character*60 task, csave
logical lsave(4)
integer n, m, iprint,
+ nbd(n), iwa(3*n), isave(44)
double precision f, factr, pgtol, x(n), l(n), u(n), g(n),
c
c-jlm-jn
+ wa(2*m*n + 5*n + 11*m*m + 8*m), dsave(29)
c ************
c
c Subroutine setulb
c
c This subroutine partitions the working arrays wa and iwa, and
c then uses the limited memory BFGS method to solve the bound
c constrained optimization problem by calling mainlb.
c (The direct method will be used in the subspace minimization.)
c
c n is an integer variable.
c On entry n is the dimension of the problem.
c On exit n is unchanged.
c
c m is an integer variable.
c On entry m is the maximum number of variable metric corrections
c used to define the limited memory matrix.
c On exit m is unchanged.
c
c x is a double precision array of dimension n.
c On entry x is an approximation to the solution.
c On exit x is the current approximation.
c
c l is a double precision array of dimension n.
c On entry l is the lower bound on x.
c On exit l is unchanged.
c
c u is a double precision array of dimension n.
c On entry u is the upper bound on x.
c On exit u is unchanged.
c
c nbd is an integer array of dimension n.
c On entry nbd represents the type of bounds imposed on the
c variables, and must be specified as follows:
c nbd(i)=0 if x(i) is unbounded,
c 1 if x(i) has only a lower bound,
c 2 if x(i) has both lower and upper bounds, and
c 3 if x(i) has only an upper bound.
c On exit nbd is unchanged.
c
c f is a double precision variable.
c On first entry f is unspecified.
c On final exit f is the value of the function at x.
c
c g is a double precision array of dimension n.
c On first entry g is unspecified.
c On final exit g is the value of the gradient at x.
c
c factr is a double precision variable.
c On entry factr >= 0 is specified by the user. The iteration
c will stop when
c
c (f^k - f^{k+1})/max{|f^k|,|f^{k+1}|,1} <= factr*epsmch
c
c where epsmch is the machine precision, which is automatically
c generated by the code. Typical values for factr: 1.d+12 for
c low accuracy; 1.d+7 for moderate accuracy; 1.d+1 for extremely
c high accuracy.
c On exit factr is unchanged.
c
c pgtol is a double precision variable.
c On entry pgtol >= 0 is specified by the user. The iteration
c will stop when
c
c max{|proj g_i | i = 1, ..., n} <= pgtol
c
c where pg_i is the ith component of the projected gradient.
c On exit pgtol is unchanged.
c
c wa is a double precision working array of length
c (2mmax + 5)nmax + 12mmax^2 + 12mmax.
c
c iwa is an integer working array of length 3nmax.
c
c task is a working string of characters of length 60 indicating
c the current job when entering and quitting this subroutine.
c
c iprint is an integer variable that must be set by the user.
c It controls the frequency and type of output generated:
c iprint<0 no output is generated;
c iprint=0 print only one line at the last iteration;
c 0<iprint<99 print also f and |proj g| every iprint iterations;
c iprint=99 print details of every iteration except n-vectors;
c iprint=100 print also the changes of active set and final x;
c iprint>100 print details of every iteration including x and g;
c When iprint > 0, the file iterate.dat will be created to
c summarize the iteration.
c
c csave is a working string of characters of length 60.
c
c lsave is a logical working array of dimension 4.
c On exit with 'task' = NEW_X, the following information is
c available:
c If lsave(1) = .true. then the initial X has been replaced by
c its projection in the feasible set;
c If lsave(2) = .true. then the problem is constrained;
c If lsave(3) = .true. then each variable has upper and lower
c bounds;
c
c isave is an integer working array of dimension 44.
c On exit with 'task' = NEW_X, the following information is
c available:
c isave(22) = the total number of intervals explored in the
c search of Cauchy points;
c isave(26) = the total number of skipped BFGS updates before
c the current iteration;
c isave(30) = the number of current iteration;
c isave(31) = the total number of BFGS updates prior the current
c iteration;
c isave(33) = the number of intervals explored in the search of
c Cauchy point in the current iteration;
c isave(34) = the total number of function and gradient
c evaluations;
c isave(36) = the number of function value or gradient
c evaluations in the current iteration;
c if isave(37) = 0 then the subspace argmin is within the box;
c if isave(37) = 1 then the subspace argmin is beyond the box;
c isave(38) = the number of free variables in the current
c iteration;
c isave(39) = the number of active constraints in the current
c iteration;
c n + 1 - isave(40) = the number of variables leaving the set of
c active constraints in the current iteration;
c isave(41) = the number of variables entering the set of active
c constraints in the current iteration.
c
c dsave is a double precision working array of dimension 29.
c On exit with 'task' = NEW_X, the following information is
c available:
c dsave(1) = current 'theta' in the BFGS matrix;
c dsave(2) = f(x) in the previous iteration;
c dsave(3) = factr*epsmch;
c dsave(4) = 2-norm of the line search direction vector;
c dsave(5) = the machine precision epsmch generated by the code;
c dsave(7) = the accumulated time spent on searching for
c Cauchy points;
c dsave(8) = the accumulated time spent on
c subspace minimization;
c dsave(9) = the accumulated time spent on line search;
c dsave(11) = the slope of the line search function at
c the current point of line search;
c dsave(12) = the maximum relative step length imposed in
c line search;
c dsave(13) = the infinity norm of the projected gradient;
c dsave(14) = the relative step length in the line search;
c dsave(15) = the slope of the line search function at
c the starting point of the line search;
c dsave(16) = the square of the 2-norm of the line search
c direction vector.
c
c Subprograms called:
c
c L-BFGS-B Library ... mainlb.
c
c
c References:
c
c [1] R. H. Byrd, P. Lu, J. Nocedal and C. Zhu, ``A limited
c memory algorithm for bound constrained optimization'',
c SIAM J. Scientific Computing 16 (1995), no. 5, pp. 1190--1208.
c
c [2] C. Zhu, R.H. Byrd, P. Lu, J. Nocedal, ``L-BFGS-B: a
c limited memory FORTRAN code for solving bound constrained
c optimization problems'', Tech. Report, NAM-11, EECS Department,
c Northwestern University, 1994.
c
c (Postscript files of these papers are available via anonymous
c ftp to eecs.nwu.edu in the directory pub/lbfgs/lbfgs_bcm.)
c
c * * *
c
c NEOS, November 1994. (Latest revision June 1996.)
c Optimization Technology Center.
c Argonne National Laboratory and Northwestern University.
c Written by
c Ciyou Zhu
c in collaboration with R.H. Byrd, P. Lu-Chen and J. Nocedal.
c
c
c ************
c-jlm-jn
integer lws,lr,lz,lt,ld,lxp,lwa,
+ lwy,lsy,lss,lwt,lwn,lsnd
if (task .eq. 'START') then
isave(1) = m*n
isave(2) = m**2
isave(3) = 4*m**2
isave(4) = 1 ! ws m*n
isave(5) = isave(4) + isave(1) ! wy m*n
isave(6) = isave(5) + isave(1) ! wsy m**2
isave(7) = isave(6) + isave(2) ! wss m**2
isave(8) = isave(7) + isave(2) ! wt m**2
isave(9) = isave(8) + isave(2) ! wn 4*m**2
isave(10) = isave(9) + isave(3) ! wsnd 4*m**2
isave(11) = isave(10) + isave(3) ! wz n
isave(12) = isave(11) + n ! wr n
isave(13) = isave(12) + n ! wd n
isave(14) = isave(13) + n ! wt n
isave(15) = isave(14) + n ! wxp n
isave(16) = isave(15) + n ! wa 8*m
endif
lws = isave(4)
lwy = isave(5)
lsy = isave(6)
lss = isave(7)
lwt = isave(8)
lwn = isave(9)
lsnd = isave(10)
lz = isave(11)
lr = isave(12)
ld = isave(13)
lt = isave(14)
lxp = isave(15)
lwa = isave(16)
call mainlb(n,m,x,l,u,nbd,f,g,factr,pgtol,
+ wa(lws),wa(lwy),wa(lsy),wa(lss), wa(lwt),
+ wa(lwn),wa(lsnd),wa(lz),wa(lr),wa(ld),wa(lt),wa(lxp),
+ wa(lwa),
+ iwa(1),iwa(n+1),iwa(2*n+1),task,iprint,
+ csave,lsave,isave(22),dsave)
return
end
c======================= The end of setulb =============================
subroutine mainlb(n, m, x, l, u, nbd, f, g, factr, pgtol, ws, wy,
+ sy, ss, wt, wn, snd, z, r, d, t, xp, wa,
+ index, iwhere, indx2, task,
+ iprint, csave, lsave, isave, dsave)
implicit none
character*60 task, csave
logical lsave(4)
integer n, m, iprint, nbd(n), index(n),
+ iwhere(n), indx2(n), isave(23)
double precision f, factr, pgtol,
+ x(n), l(n), u(n), g(n), z(n), r(n), d(n), t(n),
c-jlm-jn
+ xp(n),
+ wa(8*m),
+ ws(n, m), wy(n, m), sy(m, m), ss(m, m),
+ wt(m, m), wn(2*m, 2*m), snd(2*m, 2*m), dsave(29)
c ************
c
c Subroutine mainlb
c
c This subroutine solves bound constrained optimization problems by
c using the compact formula of the limited memory BFGS updates.
c
c n is an integer variable.
c On entry n is the number of variables.
c On exit n is unchanged.
c
c m is an integer variable.
c On entry m is the maximum number of variable metric
c corrections allowed in the limited memory matrix.
c On exit m is unchanged.
c
c x is a double precision array of dimension n.
c On entry x is an approximation to the solution.
c On exit x is the current approximation.
c
c l is a double precision array of dimension n.
c On entry l is the lower bound of x.
c On exit l is unchanged.
c
c u is a double precision array of dimension n.
c On entry u is the upper bound of x.
c On exit u is unchanged.
c
c nbd is an integer array of dimension n.
c On entry nbd represents the type of bounds imposed on the
c variables, and must be specified as follows:
c nbd(i)=0 if x(i) is unbounded,
c 1 if x(i) has only a lower bound,
c 2 if x(i) has both lower and upper bounds,
c 3 if x(i) has only an upper bound.
c On exit nbd is unchanged.
c
c f is a double precision variable.
c On first entry f is unspecified.
c On final exit f is the value of the function at x.
c
c g is a double precision array of dimension n.
c On first entry g is unspecified.
c On final exit g is the value of the gradient at x.
c
c factr is a double precision variable.
c On entry factr >= 0 is specified by the user. The iteration
c will stop when
c
c (f^k - f^{k+1})/max{|f^k|,|f^{k+1}|,1} <= factr*epsmch
c
c where epsmch is the machine precision, which is automatically
c generated by the code.
c On exit factr is unchanged.
c
c pgtol is a double precision variable.
c On entry pgtol >= 0 is specified by the user. The iteration
c will stop when
c
c max{|proj g_i | i = 1, ..., n} <= pgtol
c
c where pg_i is the ith component of the projected gradient.
c On exit pgtol is unchanged.
c
c ws, wy, sy, and wt are double precision working arrays used to
c store the following information defining the limited memory
c BFGS matrix:
c ws, of dimension n x m, stores S, the matrix of s-vectors;
c wy, of dimension n x m, stores Y, the matrix of y-vectors;
c sy, of dimension m x m, stores S'Y;
c ss, of dimension m x m, stores S'S;
c yy, of dimension m x m, stores Y'Y;
c wt, of dimension m x m, stores the Cholesky factorization
c of (theta*S'S+LD^(-1)L'); see eq.
c (2.26) in [3].
c
c wn is a double precision working array of dimension 2m x 2m
c used to store the LEL^T factorization of the indefinite matrix
c K = [-D -Y'ZZ'Y/theta L_a'-R_z' ]
c [L_a -R_z theta*S'AA'S ]
c
c where E = [-I 0]
c [ 0 I]
c
c snd is a double precision working array of dimension 2m x 2m
c used to store the lower triangular part of
c N = [Y' ZZ'Y L_a'+R_z']
c [L_a +R_z S'AA'S ]
c
c z(n),r(n),d(n),t(n), xp(n),wa(8*m) are double precision working arrays.
c z is used at different times to store the Cauchy point and
c the Newton point.
c xp is used to safeguard the projected Newton direction
c
c sg(m),sgo(m),yg(m),ygo(m) are double precision working arrays.
c
c index is an integer working array of dimension n.
c In subroutine freev, index is used to store the free and fixed
c variables at the Generalized Cauchy Point (GCP).
c
c iwhere is an integer working array of dimension n used to record
c the status of the vector x for GCP computation.
c iwhere(i)=0 or -3 if x(i) is free and has bounds,
c 1 if x(i) is fixed at l(i), and l(i) .ne. u(i)
c 2 if x(i) is fixed at u(i), and u(i) .ne. l(i)
c 3 if x(i) is always fixed, i.e., u(i)=x(i)=l(i)
c -1 if x(i) is always free, i.e., no bounds on it.
c
c indx2 is an integer working array of dimension n.
c Within subroutine cauchy, indx2 corresponds to the array iorder.
c In subroutine freev, a list of variables entering and leaving
c the free set is stored in indx2, and it is passed on to
c subroutine formk with this information.
c
c task is a working string of characters of length 60 indicating
c the current job when entering and leaving this subroutine.
c
c iprint is an INTEGER variable that must be set by the user.
c It controls the frequency and type of output generated:
c iprint<0 no output is generated;
c iprint=0 print only one line at the last iteration;
c 0<iprint<99 print also f and |proj g| every iprint iterations;
c iprint=99 print details of every iteration except n-vectors;
c iprint=100 print also the changes of active set and final x;
c iprint>100 print details of every iteration including x and g;
c When iprint > 0, the file iterate.dat will be created to
c summarize the iteration.
c
c csave is a working string of characters of length 60.
c
c lsave is a logical working array of dimension 4.
c
c isave is an integer working array of dimension 23.
c
c dsave is a double precision working array of dimension 29.
c
c
c Subprograms called
c
c L-BFGS-B Library ... cauchy, subsm, lnsrlb, formk,
c
c errclb, prn1lb, prn2lb, prn3lb, active, projgr,
c
c freev, cmprlb, matupd, formt.
c
c Minpack2 Library ... timer
c
c Linpack Library ... dcopy, ddot.
c
c
c References:
c
c [1] R. H. Byrd, P. Lu, J. Nocedal and C. Zhu, ``A limited
c memory algorithm for bound constrained optimization'',
c SIAM J. Scientific Computing 16 (1995), no. 5, pp. 1190--1208.
c
c [2] C. Zhu, R.H. Byrd, P. Lu, J. Nocedal, ``L-BFGS-B: FORTRAN
c Subroutines for Large Scale Bound Constrained Optimization''
c Tech. Report, NAM-11, EECS Department, Northwestern University,
c 1994.
c
c [3] R. Byrd, J. Nocedal and R. Schnabel "Representations of
c Quasi-Newton Matrices and their use in Limited Memory Methods'',
c Mathematical Programming 63 (1994), no. 4, pp. 129-156.
c
c (Postscript files of these papers are available via anonymous
c ftp to eecs.nwu.edu in the directory pub/lbfgs/lbfgs_bcm.)
c
c * * *
c
c NEOS, November 1994. (Latest revision June 1996.)
c Optimization Technology Center.
c Argonne National Laboratory and Northwestern University.
c Written by
c Ciyou Zhu
c in collaboration with R.H. Byrd, P. Lu-Chen and J. Nocedal.
c
c
c ************
logical prjctd,cnstnd,boxed,updatd,wrk
character*3 word
integer i,k,nintol,itfile,iback,nskip,
+ head,col,iter,itail,iupdat,
+ nseg,nfgv,info,ifun,
+ iword,nfree,nact,ileave,nenter
double precision theta,fold,ddot,dr,rr,tol,
+ xstep,sbgnrm,ddum,dnorm,dtd,epsmch,
+ cpu1,cpu2,cachyt,sbtime,lnscht,time1,time2,
+ gd,gdold,stp,stpmx,time
double precision one,zero
parameter (one=1.0d0,zero=0.0d0)
if (task .eq. 'START') then
epsmch = epsilon(one)
call timer(time1)
c Initialize counters and scalars when task='START'.
c for the limited memory BFGS matrices:
col = 0
head = 1
theta = one
iupdat = 0
updatd = .false.
iback = 0
itail = 0
iword = 0
nact = 0
ileave = 0
nenter = 0
fold = zero
dnorm = zero
cpu1 = zero
gd = zero
stpmx = zero
sbgnrm = zero
stp = zero
gdold = zero
dtd = zero
c for operation counts:
iter = 0
nfgv = 0
nseg = 0
nintol = 0
nskip = 0
nfree = n
ifun = 0
c for stopping tolerance:
tol = factr*epsmch
c for measuring running time:
cachyt = 0
sbtime = 0
lnscht = 0
c 'word' records the status of subspace solutions.
word = '---'
c 'info' records the termination information.
info = 0
itfile = 8
if (iprint .ge. 1) then
c open a summary file 'iterate.dat'
open (8, file = 'iterate.dat', status = 'unknown')
endif
c Check the input arguments for errors.
call errclb(n,m,factr,l,u,nbd,task,info,k)
if (task(1:5) .eq. 'ERROR') then
call prn3lb(n,x,f,task,iprint,info,itfile,
+ iter,nfgv,nintol,nskip,nact,sbgnrm,
+ zero,nseg,word,iback,stp,xstep,k,
+ cachyt,sbtime,lnscht)
return
endif
call prn1lb(n,m,l,u,x,iprint,itfile,epsmch)
c Initialize iwhere & project x onto the feasible set.
call active(n,l,u,nbd,x,iwhere,iprint,prjctd,cnstnd,boxed)
c The end of the initialization.
else
c restore local variables.
prjctd = lsave(1)
cnstnd = lsave(2)
boxed = lsave(3)
updatd = lsave(4)
nintol = isave(1)
itfile = isave(3)
iback = isave(4)
nskip = isave(5)
head = isave(6)
col = isave(7)
itail = isave(8)
iter = isave(9)
iupdat = isave(10)
nseg = isave(12)
nfgv = isave(13)
info = isave(14)
ifun = isave(15)
iword = isave(16)
nfree = isave(17)
nact = isave(18)
ileave = isave(19)
nenter = isave(20)
theta = dsave(1)
fold = dsave(2)
tol = dsave(3)
dnorm = dsave(4)
epsmch = dsave(5)
cpu1 = dsave(6)
cachyt = dsave(7)
sbtime = dsave(8)
lnscht = dsave(9)
time1 = dsave(10)
gd = dsave(11)
stpmx = dsave(12)
sbgnrm = dsave(13)
stp = dsave(14)
gdold = dsave(15)
dtd = dsave(16)
c After returning from the driver go to the point where execution
c is to resume.
if (task(1:5) .eq. 'FG_LN') goto 666
if (task(1:5) .eq. 'NEW_X') goto 777
if (task(1:5) .eq. 'FG_ST') goto 111
if (task(1:4) .eq. 'STOP') then
if (task(7:9) .eq. 'CPU') then
c restore the previous iterate.
call dcopy(n,t,1,x,1)
call dcopy(n,r,1,g,1)
f = fold
endif
goto 999
endif
endif
c Compute f0 and g0.
task = 'FG_START'
c return to the driver to calculate f and g; reenter at 111.
goto 1000
111 continue
nfgv = 1
c Compute the infinity norm of the (-) projected gradient.
call projgr(n,l,u,nbd,x,g,sbgnrm)
if (iprint .ge. 1) then
write (6,1002) iter,f,sbgnrm
write (itfile,1003) iter,nfgv,sbgnrm,f
endif
if (sbgnrm .le. pgtol) then
c terminate the algorithm.
task = 'CONVERGENCE: NORM_OF_PROJECTED_GRADIENT_<=_PGTOL'
goto 999
endif
c ----------------- the beginning of the loop --------------------------
222 continue
if (iprint .ge. 99) write (6,1001) iter + 1
iword = -1
c
if (.not. cnstnd .and. col .gt. 0) then
c skip the search for GCP.
call dcopy(n,x,1,z,1)
wrk = updatd
nseg = 0
goto 333
endif
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
c
c Compute the Generalized Cauchy Point (GCP).
c
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
call timer(cpu1)
call cauchy(n,x,l,u,nbd,g,indx2,iwhere,t,d,z,
+ m,wy,ws,sy,wt,theta,col,head,
+ wa(1),wa(2*m+1),wa(4*m+1),wa(6*m+1),nseg,
+ iprint, sbgnrm, info, epsmch)
if (info .ne. 0) then
c singular triangular system detected; refresh the lbfgs memory.
if(iprint .ge. 1) write (6, 1005)
info = 0
col = 0
head = 1
theta = one
iupdat = 0
updatd = .false.
call timer(cpu2)
cachyt = cachyt + cpu2 - cpu1
goto 222
endif
call timer(cpu2)
cachyt = cachyt + cpu2 - cpu1
nintol = nintol + nseg
c Count the entering and leaving variables for iter > 0;
c find the index set of free and active variables at the GCP.
call freev(n,nfree,index,nenter,ileave,indx2,
+ iwhere,wrk,updatd,cnstnd,iprint,iter)
nact = n - nfree
333 continue
c If there are no free variables or B=theta*I, then
c skip the subspace minimization.
if (nfree .eq. 0 .or. col .eq. 0) goto 555
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
c
c Subspace minimization.
c
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
call timer(cpu1)
c Form the LEL^T factorization of the indefinite
c matrix K = [-D -Y'ZZ'Y/theta L_a'-R_z' ]
c [L_a -R_z theta*S'AA'S ]
c where E = [-I 0]
c [ 0 I]
if (wrk) call formk(n,nfree,index,nenter,ileave,indx2,iupdat,
+ updatd,wn,snd,m,ws,wy,sy,theta,col,head,info)
if (info .ne. 0) then
c nonpositive definiteness in Cholesky factorization;
c refresh the lbfgs memory and restart the iteration.
if(iprint .ge. 1) write (6, 1006)
info = 0
col = 0
head = 1
theta = one
iupdat = 0
updatd = .false.
call timer(cpu2)
sbtime = sbtime + cpu2 - cpu1
goto 222
endif
c compute r=-Z'B(xcp-xk)-Z'g (using wa(2m+1)=W'(xcp-x)
c from 'cauchy').
call cmprlb(n,m,x,g,ws,wy,sy,wt,z,r,wa,index,
+ theta,col,head,nfree,cnstnd,info)
if (info .ne. 0) goto 444
c-jlm-jn call the direct method.
call subsm( n, m, nfree, index, l, u, nbd, z, r, xp, ws, wy,
+ theta, x, g, col, head, iword, wa, wn, iprint, info)
444 continue
if (info .ne. 0) then
c singular triangular system detected;
c refresh the lbfgs memory and restart the iteration.
if(iprint .ge. 1) write (6, 1005)
info = 0
col = 0
head = 1
theta = one
iupdat = 0
updatd = .false.
call timer(cpu2)
sbtime = sbtime + cpu2 - cpu1
goto 222
endif
call timer(cpu2)
sbtime = sbtime + cpu2 - cpu1
555 continue
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
c
c Line search and optimality tests.
c
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
c Generate the search direction d:=z-x.
do 40 i = 1, n
d(i) = z(i) - x(i)
40 continue
call timer(cpu1)
666 continue
call lnsrlb(n,l,u,nbd,x,f,fold,gd,gdold,g,d,r,t,z,stp,dnorm,
+ dtd,xstep,stpmx,iter,ifun,iback,nfgv,info,task,
+ boxed,cnstnd,csave,isave(22),dsave(17))
if (info .ne. 0 .or. iback .ge. 20) then
c restore the previous iterate.
call dcopy(n,t,1,x,1)
call dcopy(n,r,1,g,1)
f = fold
if (col .eq. 0) then
c abnormal termination.
if (info .eq. 0) then
info = -9
c restore the actual number of f and g evaluations etc.
nfgv = nfgv - 1
ifun = ifun - 1
iback = iback - 1
endif
task = 'ABNORMAL_TERMINATION_IN_LNSRCH'
iter = iter + 1
goto 999
else
c refresh the lbfgs memory and restart the iteration.
if(iprint .ge. 1) write (6, 1008)
if (info .eq. 0) nfgv = nfgv - 1
info = 0
col = 0
head = 1
theta = one
iupdat = 0
updatd = .false.
task = 'RESTART_FROM_LNSRCH'
call timer(cpu2)
lnscht = lnscht + cpu2 - cpu1
goto 222
endif
else if (task(1:5) .eq. 'FG_LN') then
c return to the driver for calculating f and g; reenter at 666.
goto 1000
else
c calculate and print out the quantities related to the new X.
call timer(cpu2)
lnscht = lnscht + cpu2 - cpu1
iter = iter + 1
c Compute the infinity norm of the projected (-)gradient.
call projgr(n,l,u,nbd,x,g,sbgnrm)
c Print iteration information.
call prn2lb(n,x,f,g,iprint,itfile,iter,nfgv,nact,
+ sbgnrm,nseg,word,iword,iback,stp,xstep)
goto 1000
endif
777 continue
c Test for termination.
if (sbgnrm .le. pgtol) then
c terminate the algorithm.
task = 'CONVERGENCE: NORM_OF_PROJECTED_GRADIENT_<=_PGTOL'
goto 999
endif
ddum = max(abs(fold), abs(f), one)
if ((fold - f) .le. tol*ddum) then
c terminate the algorithm.
task = 'CONVERGENCE: REL_REDUCTION_OF_F_<=_FACTR*EPSMCH'
if (iback .ge. 10) info = -5
c i.e., to issue a warning if iback>10 in the line search.
goto 999
endif
c Compute d=newx-oldx, r=newg-oldg, rr=y'y and dr=y's.
do 42 i = 1, n
r(i) = g(i) - r(i)
42 continue
rr = ddot(n,r,1,r,1)
if (stp .eq. one) then
dr = gd - gdold
ddum = -gdold
else
dr = (gd - gdold)*stp
call dscal(n,stp,d,1)
ddum = -gdold*stp
endif
if (dr .le. epsmch*ddum) then
c skip the L-BFGS update.
nskip = nskip + 1
updatd = .false.
if (iprint .ge. 1) write (6,1004) dr, ddum
goto 888
endif
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
c
c Update the L-BFGS matrix.
c
cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
updatd = .true.
iupdat = iupdat + 1
c Update matrices WS and WY and form the middle matrix in B.
call matupd(n,m,ws,wy,sy,ss,d,r,itail,
+ iupdat,col,head,theta,rr,dr,stp,dtd)
c Form the upper half of the pds T = theta*SS + L*D^(-1)*L';
c Store T in the upper triangular of the array wt;
c Cholesky factorize T to J*J' with
c J' stored in the upper triangular of wt.
call formt(m,wt,sy,ss,col,theta,info)
if (info .ne. 0) then
c nonpositive definiteness in Cholesky factorization;
c refresh the lbfgs memory and restart the iteration.
if(iprint .ge. 1) write (6, 1007)
info = 0
col = 0
head = 1
theta = one
iupdat = 0
updatd = .false.
goto 222
endif
c Now the inverse of the middle matrix in B is
c [ D^(1/2) O ] [ -D^(1/2) D^(-1/2)*L' ]
c [ -L*D^(-1/2) J ] [ 0 J' ]
888 continue
c -------------------- the end of the loop -----------------------------
goto 222
999 continue
call timer(time2)
time = time2 - time1
call prn3lb(n,x,f,task,iprint,info,itfile,
+ iter,nfgv,nintol,nskip,nact,sbgnrm,
+ time,nseg,word,iback,stp,xstep,k,
+ cachyt,sbtime,lnscht)
1000 continue
c Save local variables.
lsave(1) = prjctd
lsave(2) = cnstnd
lsave(3) = boxed
lsave(4) = updatd
isave(1) = nintol
isave(3) = itfile
isave(4) = iback
isave(5) = nskip
isave(6) = head
isave(7) = col
isave(8) = itail
isave(9) = iter
isave(10) = iupdat
isave(12) = nseg
isave(13) = nfgv
isave(14) = info
isave(15) = ifun
isave(16) = iword
isave(17) = nfree
isave(18) = nact
isave(19) = ileave
isave(20) = nenter
dsave(1) = theta
dsave(2) = fold
dsave(3) = tol
dsave(4) = dnorm
dsave(5) = epsmch
dsave(6) = cpu1
dsave(7) = cachyt
dsave(8) = sbtime
dsave(9) = lnscht
dsave(10) = time1
dsave(11) = gd
dsave(12) = stpmx
dsave(13) = sbgnrm
dsave(14) = stp
dsave(15) = gdold
dsave(16) = dtd
1001 format (//,'ITERATION ',i5)
1002 format
+ (/,'At iterate',i5,4x,'f= ',1p,d12.5,4x,'|proj g|= ',1p,d12.5)
1003 format (2(1x,i4),5x,'-',5x,'-',3x,'-',5x,'-',5x,'-',8x,'-',3x,
+ 1p,2(1x,d10.3))
1004 format (' ys=',1p,e10.3,' -gs=',1p,e10.3,' BFGS update SKIPPED')
1005 format (/,
+' Singular triangular system detected;',/,
+' refresh the lbfgs memory and restart the iteration.')
1006 format (/,
+' Nonpositive definiteness in Cholesky factorization in formk;',/,
+' refresh the lbfgs memory and restart the iteration.')
1007 format (/,
+' Nonpositive definiteness in Cholesky factorization in formt;',/,
+' refresh the lbfgs memory and restart the iteration.')
1008 format (/,
+' Bad direction in the line search;',/,
+' refresh the lbfgs memory and restart the iteration.')
return
end