-
Notifications
You must be signed in to change notification settings - Fork 1
/
tn.f
1794 lines (1781 loc) · 53.2 KB
/
tn.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
subroutine lmqn (ifail, n, x, f, g, w, lw, sfun,
* msglvl, maxit, maxfun, eta, stepmx, accrcy, xtol)
implicit double precision (a-h,o-z)
integer msglvl, n, maxfun, ifail, lw, ipivot(1)
double precision x(n), g(n), w(lw), eta, xtol, stepmx, f, accrcy
c
c this routine solves the optimization problem
c
c minimize f(x)
c x
c
c where x is a vector of n real variables. the method used is
c a truncated-newton algorithm (see "newton-type minimization via
c the lanczos method" by s.g. nash (siam j. numer. anal. 21 (1984),
c pp. 770-778). this algorithm finds a local minimum of f(x). it does
c not assume that the function f is convex (and so cannot guarantee a
c global solution), but does assume that the function is bounded below.
c it can solve problems having any number of variables, but it is
c especially useful when the number of variables (n) is large.
c
c subroutine parameters:
c
c ierror - (integer) error code
c ( 0 => normal return)
c ( 2 => more than maxfun evaluations)
c ( 3 => line search failed to find
c ( lower point (may not be serious)
c (-1 => error in input parameters)
c n - (integer) number of variables
c x - (real*8) vector of length at least n; on input, an initial
c estimate of the solution; on output, the computed solution.
c g - (real*8) vector of length at least n; on output, the final
c value of the gradient
c f - (real*8) on input, a rough estimate of the value of the
c objective function at the solution; on output, the value
c of the objective function at the solution
c w - (real*8) work vector of length at least 14*n
c lw - (integer) the declared dimension of w
c sfun - a user-specified subroutine that computes the function
c and gradient of the objective function. it must have
c the calling sequence
c subroutine sfun (n, x, f, g)
c integer n
c double precision x(n), g(n), f
c
c eta - severity of the linesearch
c maxfun - maximum allowable number of function evaluations
c xtol - desired accuracy for the solution x*
c stepmx - maximum allowable step in the linesearch
c accrcy - accuracy of computed function values
c msglvl - determines quantity of printed output
c 0 = none, 1 = one line per major iteration.
c maxit - maximum number of inner iterations per step
c
c this routine is a truncated-newton method.
c the truncated-newton method is preconditioned by a limited-memory
c quasi-newton method (this preconditioning strategy is developed
c in this routine) with a further diagonal scaling (see routine ndia3).
c for further details on the parameters, see routine tn.
c
integer i, icycle, ioldg, ipk, iyk, loldg, lpk, lsr,
* lwtest, lyk, lyr, nftotl, niter, nm1, numf, nwhy
double precision abstol, alpha, difnew, difold, epsmch,
* epsred, fkeep, fm, fnew, fold, fstop, ftest, gnorm, gsk,
* gtg, gtpnew, oldf, oldgtp, one, pe, peps, pnorm, reltol,
* rteps, rtleps, rtol, rtolsq, small, spe, tiny,
* tnytol, toleps, xnorm, yksk, yrsr, zero
logical lreset, upd1
c
c the following imsl and standard functions are used
c
double precision dabs, ddot, dsqrt, step1, dnrm2
external sfun
common /subscr/ lgv,lz1,lzk,lv,lsk,lyk,ldiagb,lsr,lyr,
* loldg,lhg,lhyk,lpk,lemat,lwtest
c
c initialize parameters and constants
c
if (msglvl .ge. 1) write(*,800)
call setpar(n)
upd1 = .true.
ireset = 0
nfeval = 0
nmodif = 0
nlincg = 0
fstop = f
zero = 0.d0
one = 1.d0
nm1 = n - 1
c
c within this routine the array w(loldg) is shared by w(lhyr)
c
lhyr = loldg
c
c check parameters and set constants
c
call chkucp(lwtest,maxfun,nwhy,n,alpha,epsmch,
* eta,peps,rteps,rtol,rtolsq,stepmx,ftest,
* xtol,xnorm,x,lw,small,tiny,accrcy)
if (nwhy .lt. 0) go to 120
call setucr(small,nftotl,niter,n,f,fnew,
* fm,gtg,oldf,sfun,g,x)
fold = fnew
if (msglvl .ge. 1) write(*,810) niter,nftotl,nlincg,fnew,gtg
c
c check for small gradient at the starting point.
c
ftest = one + dabs(fnew)
if (gtg .lt. 1.d-4*epsmch*ftest*ftest) go to 90
c
c set initial values to other parameters
c
icycle = nm1
toleps = rtol + rteps
rtleps = rtolsq + epsmch
gnorm = dsqrt(gtg)
difnew = zero
epsred = 5.0d-2
fkeep = fnew
c
c set the diagonal of the approximate hessian to unity.
c
idiagb = ldiagb
do 10 i = 1,n
w(idiagb) = one
idiagb = idiagb + 1
10 continue
c
c ..................start of main iterative loop..........
c
c compute the new search direction
c
modet = msglvl - 3
call modlnp(modet,w(lpk),w(lgv),w(lz1),w(lv),
* w(ldiagb),w(lemat),x,g,w(lzk),
* n,w,lw,niter,maxit,nfeval,nmodif,
* nlincg,upd1,yksk,gsk,yrsr,lreset,sfun,.false.,ipivot,
* accrcy,gtpnew,gnorm,xnorm)
20 continue
call dcopy(n,g,1,w(loldg),1)
pnorm = dnrm2(n,w(lpk),1)
oldf = fnew
oldgtp = gtpnew
c
c prepare to compute the step length
c
pe = pnorm + epsmch
c
c compute the absolute and relative tolerances for the linear search
c
reltol = rteps*(xnorm + one)/pe
abstol = - epsmch*ftest/(oldgtp - epsmch)
c
c compute the smallest allowable spacing between points in
c the linear search
c
tnytol = epsmch*(xnorm + one)/pe
spe = stepmx/pe
c
c set the initial step length.
c
alpha = step1(fnew,fm,oldgtp,spe)
c
c perform the linear search
c
call linder(n,sfun,small,epsmch,reltol,abstol,tnytol,
* eta,zero,spe,w(lpk),oldgtp,x,fnew,alpha,g,numf,
* nwhy,w,lw)
c
fold = fnew
niter = niter + 1
nftotl = nftotl + numf
gtg = ddot(n,g,1,g,1)
if (msglvl .ge. 1) write(*,810) niter,nftotl,nlincg,fnew,gtg
if (nwhy .lt. 0) go to 120
if (nwhy .eq. 0 .or. nwhy .eq. 2) go to 30
c
c the linear search has failed to find a lower point
c
nwhy = 3
go to 100
30 if (nwhy .le. 1) go to 40
call sfun(n,x,fnew,g)
nftotl = nftotl + 1
c
c terminate if more than maxfun evalutations have been made
c
40 nwhy = 2
if (nftotl .gt. maxfun) go to 110
nwhy = 0
c
c set up parameters used in convergence and resetting tests
c
difold = difnew
difnew = oldf - fnew
c
c if this is the first iteration of a new cycle, compute the
c percentage reduction factor for the resetting test.
c
if (icycle .ne. 1) go to 50
if (difnew .gt. 2.0d0 *difold) epsred = epsred + epsred
if (difnew .lt. 5.0d-1*difold) epsred = 5.0d-1*epsred
50 continue
gnorm = dsqrt(gtg)
ftest = one + dabs(fnew)
xnorm = dnrm2(n,x,1)
c
c test for convergence
c
if ((alpha*pnorm .lt. toleps*(one + xnorm)
* .and. dabs(difnew) .lt. rtleps*ftest
* .and. gtg .lt. peps*ftest*ftest)
* .or. gtg .lt. 1.d-4*accrcy*ftest*ftest) go to 90
c
c compute the change in the iterates and the corresponding change
c in the gradients
c
isk = lsk
ipk = lpk
iyk = lyk
ioldg = loldg
do 60 i = 1,n
w(iyk) = g(i) - w(ioldg)
w(isk) = alpha*w(ipk)
ipk = ipk + 1
isk = isk + 1
iyk = iyk + 1
ioldg = ioldg + 1
60 continue
c
c set up parameters used in updating the direction of search.
c
yksk = ddot(n,w(lyk),1,w(lsk),1)
lreset = .false.
if (icycle .eq. nm1 .or. difnew .lt.
* epsred*(fkeep-fnew)) lreset = .true.
if (lreset) go to 70
yrsr = ddot(n,w(lyr),1,w(lsr),1)
if (yrsr .le. zero) lreset = .true.
70 continue
upd1 = .false.
c
c compute the new search direction
c
modet = msglvl - 3
call modlnp(modet,w(lpk),w(lgv),w(lz1),w(lv),
* w(ldiagb),w(lemat),x,g,w(lzk),
* n,w,lw,niter,maxit,nfeval,nmodif,
* nlincg,upd1,yksk,gsk,yrsr,lreset,sfun,.false.,ipivot,
* accrcy,gtpnew,gnorm,xnorm)
if (lreset) go to 80
c
c store the accumulated change in the point and gradient as an
c "average" direction for preconditioning.
c
call dxpy(n,w(lsk),1,w(lsr),1)
call dxpy(n,w(lyk),1,w(lyr),1)
icycle = icycle + 1
goto 20
c
c reset
c
80 ireset = ireset + 1
c
c initialize the sum of all the changes in x.
c
call dcopy(n,w(lsk),1,w(lsr),1)
call dcopy(n,w(lyk),1,w(lyr),1)
fkeep = fnew
icycle = 1
go to 20
c
c ...............end of main iteration.......................
c
90 ifail = 0
f = fnew
go to 130
100 oldf = fnew
c
c local search here could be installed here
c
110 f = oldf
c
c set ifail
c
120 ifail = nwhy
130 if (msglvl .ge. 1 .and. ifail .lt. 0) then
if (ifail .ne. 0) write(*,820) ifail
return
endif
if (msglvl .ge. 1) then
if (ifail .ne. 0) write(*,820) ifail
write(*,830) f
write(*,840)
nmax = 10
if (n .lt. nmax) nmax = n
write(*,850) (i,x(i),i=1,nmax)
endif
return
800 format(//' nit nf cg', 9x, 'f', 21x, 'gtg',//)
810 format(' ',i3,1x,i4,1x,i4,1x,1pd22.15,2x,1pd15.8)
820 format(//,' error code =', i3)
830 format(//,' optimal function value = ', 1pd22.15)
840 format(10x, 'current solution is (at most 10 components)', /,
* 14x, 'i', 11x, 'x(i)')
850 format(10x, i5, 2x, 1pd22.15)
end
c
c
subroutine lmqnbc (ifail, n, x, f, g, w, lw, sfun, low, up,
* ipivot, msglvl, maxit, maxfun, eta, stepmx, accrcy, xtol)
implicit double precision (a-h,o-z)
integer msglvl,n,maxfun,ifail,lw
integer ipivot(n)
double precision eta,xtol,stepmx,f,accrcy
double precision x(n),g(n),w(lw),low(n),up(n)
c this routine solves the optimization problem
c
c minimize f(x)
c x
c subject to low <= x <= up
c
c where x is a vector of n real variables. the method used is
c a truncated-newton algorithm (see "newton-type minimization via
c the lanczos algorithm" by s.g. nash (technical report 378, math.
c the lanczos method" by s.g. nash (siam j. numer. anal. 21 (1984),
c pp. 770-778). this algorithm finds a local minimum of f(x). it does
c not assume that the function f is convex (and so cannot guarantee a
c global solution), but does assume that the function is bounded below.
c it can solve problems having any number of variables, but it is
c especially useful when the number of variables (n) is large.
c
c subroutine parameters:
c
c ierror - (integer) error code
c ( 0 => normal return
c ( 2 => more than maxfun evaluations
c ( 3 => line search failed to find lower
c ( point (may not be serious)
c (-1 => error in input parameters
c n - (integer) number of variables
c x - (real*8) vector of length at least n; on input, an initial
c estimate of the solution; on output, the computed solution.
c g - (real*8) vector of length at least n; on output, the final
c value of the gradient
c f - (real*8) on input, a rough estimate of the value of the
c objective function at the solution; on output, the value
c of the objective function at the solution
c w - (real*8) work vector of length at least 14*n
c lw - (integer) the declared dimension of w
c sfun - a user-specified subroutine that computes the function
c and gradient of the objective function. it must have
c the calling sequence
c subroutine sfun (n, x, f, g)
c integer n
c double precision x(n), g(n), f
c low, up - (real*8) vectors of length at least n containing
c the lower and upper bounds on the variables. if
c there are no bounds on a particular variable, set
c the bounds to -1.d38 and 1.d38, respectively.
c ipivot - (integer) work vector of length at least n, used
c to record which variables are at their bounds.
c
c eta - severity of the linesearch
c maxfun - maximum allowable number of function evaluations
c xtol - desired accuracy for the solution x*
c stepmx - maximum allowable step in the linesearch
c accrcy - accuracy of computed function values
c msglvl - controls quantity of printed output
c 0 = none, 1 = one line per major iteration.
c maxit - maximum number of inner iterations per step
c
c this routine is a bounds-constrained truncated-newton method.
c the truncated-newton method is preconditioned by a limited-memory
c quasi-newton method (this preconditioning strategy is developed
c in this routine) with a further diagonal scaling (see routine ndia3).
c for further details on the parameters, see routine tnbc.
c
integer i, icycle, ioldg, ipk, iyk, loldg, lpk, lsr,
* lwtest, lyk, lyr, nftotl, niter, nm1, numf, nwhy
double precision abstol, alpha, difnew, difold, epsmch, epsred,
* fkeep, flast, fm, fnew, fold, fstop, ftest, gnorm, gsk,
* gtg, gtpnew, oldf, oldgtp, one, pe, peps, pnorm, reltol,
* rteps, rtleps, rtol, rtolsq, small, spe, tiny,
* tnytol, toleps, xnorm, yksk, yrsr, zero
logical conv, lreset, upd1, newcon
c
c the following standard functions and system functions are used
c
double precision dabs, ddot, dnrm2, dsqrt, step1
external sfun
common/subscr/ lgv, lz1, lzk, lv, lsk, lyk, ldiagb, lsr, lyr,
* loldg, lhg, lhyk, lpk, lemat, lwtest
c
c check that initial x is feasible and that the bounds are consistent
c
call crash(n,x,ipivot,low,up,ier)
if (ier .ne. 0) write(*,800)
if (ier .ne. 0) return
if (msglvl .ge. 1) write(*,810)
c
c initialize variables
c
call setpar(n)
upd1 = .true.
ireset = 0
nfeval = 0
nmodif = 0
nlincg = 0
fstop = f
conv = .false.
zero = 0.d0
one = 1.d0
nm1 = n - 1
c
c within this routine the array w(loldg) is shared by w(lhyr)
c
lhyr = loldg
c
c check parameters and set constants
c
call chkucp(lwtest,maxfun,nwhy,n,alpha,epsmch,
* eta,peps,rteps,rtol,rtolsq,stepmx,ftest,
* xtol,xnorm,x,lw,small,tiny,accrcy)
if (nwhy .lt. 0) go to 160
call setucr(small,nftotl,niter,n,f,fnew,
* fm,gtg,oldf,sfun,g,x)
fold = fnew
flast = fnew
c
c test the lagrange multipliers to see if they are non-negative.
c because the constraints are only lower bounds, the components
c of the gradient corresponding to the active constraints are the
c lagrange multipliers. afterwords, the projected gradient is formed.
c
do 10 i = 1,n
if (ipivot(i) .eq. 2) go to 10
if (-ipivot(i)*g(i) .ge. 0.d0) go to 10
ipivot(i) = 0
10 continue
call ztime(n,g,ipivot)
gtg = ddot(n,g,1,g,1)
if (msglvl .ge. 1)
* call monit(n,x,fnew,g,niter,nftotl,nfeval,lreset,ipivot)
c
c check if the initial point is a local minimum.
c
ftest = one + dabs(fnew)
if (gtg .lt. 1.d-4*epsmch*ftest*ftest) go to 130
c
c set initial values to other parameters
c
icycle = nm1
toleps = rtol + rteps
rtleps = rtolsq + epsmch
gnorm = dsqrt(gtg)
difnew = zero
epsred = 5.0d-2
fkeep = fnew
c
c set the diagonal of the approximate hessian to unity.
c
idiagb = ldiagb
do 15 i = 1,n
w(idiagb) = one
idiagb = idiagb + 1
15 continue
c
c ..................start of main iterative loop..........
c
c compute the new search direction
c
modet = msglvl - 3
call modlnp(modet,w(lpk),w(lgv),w(lz1),w(lv),
* w(ldiagb),w(lemat),x,g,w(lzk),
* n,w,lw,niter,maxit,nfeval,nmodif,
* nlincg,upd1,yksk,gsk,yrsr,lreset,sfun,.true.,ipivot,
* accrcy,gtpnew,gnorm,xnorm)
20 continue
call dcopy(n,g,1,w(loldg),1)
pnorm = dnrm2(n,w(lpk),1)
oldf = fnew
oldgtp = gtpnew
c
c prepare to compute the step length
c
pe = pnorm + epsmch
c
c compute the absolute and relative tolerances for the linear search
c
reltol = rteps*(xnorm + one)/pe
abstol = - epsmch*ftest/(oldgtp - epsmch)
c
c compute the smallest allowable spacing between points in
c the linear search
c
tnytol = epsmch*(xnorm + one)/pe
call stpmax(stepmx,pe,spe,n,x,w(lpk),ipivot,low,up)
c
c set the initial step length.
c
alpha = step1(fnew,fm,oldgtp,spe)
c
c perform the linear search
c
call linder(n,sfun,small,epsmch,reltol,abstol,tnytol,
* eta,zero,spe,w(lpk),oldgtp,x,fnew,alpha,g,numf,
* nwhy,w,lw)
newcon = .false.
if (dabs(alpha-spe) .gt. 1.d1*epsmch) go to 30
newcon = .true.
nwhy = 0
call modz(n,x,w(lpk),ipivot,epsmch,low,up,flast,fnew)
flast = fnew
c
30 if (msglvl .ge. 2) write(*,820) alpha,pnorm
fold = fnew
niter = niter + 1
nftotl = nftotl + numf
c
c if required, print the details of this iteration
c
if (msglvl .ge. 1)
* call monit(n,x,fnew,g,niter,nftotl,nfeval,lreset,ipivot)
if (nwhy .lt. 0) go to 160
if (nwhy .eq. 0 .or. nwhy .eq. 2) go to 40
c
c the linear search has failed to find a lower point
c
nwhy = 3
go to 140
40 if (nwhy .le. 1) go to 50
call sfun(n,x,fnew,g)
nftotl = nftotl + 1
c
c terminate if more than maxfun evaluations have been made
c
50 nwhy = 2
if (nftotl .gt. maxfun) go to 150
nwhy = 0
c
c set up parameters used in convergence and resetting tests
c
difold = difnew
difnew = oldf - fnew
c
c if this is the first iteration of a new cycle, compute the
c percentage reduction factor for the resetting test.
c
if (icycle .ne. 1) go to 60
if (difnew .gt. 2.d0*difold) epsred = epsred + epsred
if (difnew .lt. 5.0d-1*difold) epsred = 5.0d-1*epsred
60 call dcopy(n,g,1,w(lgv),1)
call ztime(n,w(lgv),ipivot)
gtg = ddot(n,w(lgv),1,w(lgv),1)
gnorm = dsqrt(gtg)
ftest = one + dabs(fnew)
xnorm = dnrm2(n,x,1)
c
c test for convergence
c
call cnvtst(conv,alpha,pnorm,toleps,xnorm,difnew,rtleps,
* ftest,gtg,peps,epsmch,gtpnew,fnew,flast,g,ipivot,n,accrcy)
if (conv) go to 130
call ztime(n,g,ipivot)
c
c compute the change in the iterates and the corresponding change
c in the gradients
c
if (newcon) go to 90
isk = lsk
ipk = lpk
iyk = lyk
ioldg = loldg
do 70 i = 1,n
w(iyk) = g(i) - w(ioldg)
w(isk) = alpha*w(ipk)
ipk = ipk + 1
isk = isk + 1
iyk = iyk + 1
ioldg = ioldg + 1
70 continue
c
c set up parameters used in updating the preconditioning strategy.
c
yksk = ddot(n,w(lyk),1,w(lsk),1)
lreset = .false.
if (icycle .eq. nm1 .or. difnew .lt.
* epsred*(fkeep-fnew)) lreset = .true.
if (lreset) go to 80
yrsr = ddot(n,w(lyr),1,w(lsr),1)
if (yrsr .le. zero) lreset = .true.
80 continue
upd1 = .false.
c
c compute the new search direction
c
90 if (upd1 .and. msglvl .ge. 2) write(*,830)
if (newcon .and. msglvl .ge. 2) write(*,840)
modet = msglvl - 3
call modlnp(modet,w(lpk),w(lgv),w(lz1),w(lv),
* w(ldiagb),w(lemat),x,g,w(lzk),
* n,w,lw,niter,maxit,nfeval,nmodif,
* nlincg,upd1,yksk,gsk,yrsr,lreset,sfun,.true.,ipivot,
* accrcy,gtpnew,gnorm,xnorm)
if (newcon) go to 20
if (lreset) go to 110
c
c compute the accumulated step and its corresponding
c gradient difference.
c
call dxpy(n,w(lsk),1,w(lsr),1)
call dxpy(n,w(lyk),1,w(lyr),1)
icycle = icycle + 1
goto 20
c
c reset
c
110 ireset = ireset + 1
c
c initialize the sum of all the changes in x.
c
call dcopy(n,w(lsk),1,w(lsr),1)
call dcopy(n,w(lyk),1,w(lyr),1)
fkeep = fnew
icycle = 1
go to 20
c
c ...............end of main iteration.......................
c
130 ifail = 0
f = fnew
goto 170
140 oldf = fnew
c
c local search could be installed here
c
150 f = oldf
if (msglvl .ge. 1) call monit(n,x,
* f,g,niter,nftotl,nfeval,.true.,ipivot)
c
c set ifail
c
160 ifail = nwhy
170 if (msglvl .ge. 1 .and. ifail .lt. 0) then
if (ifail .ne. 0) write(*,850) ifail
return
endif
if (msglvl .ge. 1) then
if (ifail .ne. 0) write(*,850) ifail
write(*,860) f
write(*,870)
nmax = 10
if (n .lt. nmax) nmax = n
write(*,880) (i,x(i),i=1,nmax)
endif
return
800 format(' there is no feasible point; terminating algorithm')
810 format(//' nit nf cg', 9x, 'f', 21x, 'gtg',//)
820 format(' linesearch results: alpha,pnorm',2(1pd12.4))
830 format(' upd1 is true - trivial preconditioning')
840 format(' newcon is true - constraint added in linesearch')
850 format(//,' error code =', i3)
860 format(//,' optimal function value = ', 1pd22.15)
870 format(10x, 'current solution is (at most 10 components)', /,
* 14x, 'i', 11x, 'x(i)')
880 format(10x, i5, 2x, 1pd22.15)
end
c
c
subroutine monit(n,x,f,g,niter,nftotl,nfeval,lreset,ipivot)
c
c print results of current iteration
c
implicit double precision (a-h,o-z)
double precision x(n),f,g(n),gtg
integer ipivot(n)
logical lreset
c
gtg = 0.d0
do 10 i = 1,n
if (ipivot(i) .ne. 0) go to 10
gtg = gtg + g(i)*g(i)
10 continue
write(*,800) niter,nftotl,nfeval,f,gtg
return
800 format(' ',i4,1x,i4,1x,i4,1x,1pd22.15,2x,1pd15.8)
end
c
c
subroutine ztime(n,x,ipivot)
implicit double precision (a-h,o-z)
double precision x(n)
integer ipivot(n)
c
c this routine multiplies the vector x by the constraint matrix z
c
do 10 i = 1,n
if (ipivot(i) .ne. 0) x(i) = 0.d0
10 continue
return
end
c
c
subroutine stpmax(stepmx,pe,spe,n,x,p,ipivot,low,up)
implicit double precision (a-h,o-z)
double precision low(n),up(n),x(n),p(n),stepmx,pe,spe,t
integer ipivot(n)
c
c compute the maximum allowable step length
c
spe = stepmx / pe
c spe is the standard (unconstrained) max step
do 10 i = 1,n
if (ipivot(i) .ne. 0) go to 10
if (p(i) .eq. 0.d0) go to 10
if (p(i) .gt. 0.d0) go to 5
t = low(i) - x(i)
if (t .gt. spe*p(i)) spe = t / p(i)
go to 10
5 t = up(i) - x(i)
if (t .lt. spe*p(i)) spe = t / p(i)
10 continue
return
end
c
c
subroutine modz(n,x,p,ipivot,epsmch,low,up,flast,fnew)
implicit double precision (a-h,o-z)
double precision x(n), p(n), epsmch, dabs, tol, low(n), up(n),
* flast, fnew
integer ipivot(n)
c
c update the constraint matrix if a new constraint is encountered
c
do 10 i = 1,n
if (ipivot(i) .ne. 0) go to 10
if (p(i) .eq. 0.d0) go to 10
if (p(i) .gt. 0.d0) go to 5
tol = 1.d1 * epsmch * (dabs(low(i)) + 1.d0)
if (x(i)-low(i) .gt. tol) go to 10
flast = fnew
ipivot(i) = -1
x(i) = low(i)
go to 10
5 tol = 1.d1 * epsmch * (dabs(up(i)) + 1.d0)
if (up(i)-x(i) .gt. tol) go to 10
flast = fnew
ipivot(i) = 1
x(i) = up(i)
10 continue
return
end
c
c
subroutine cnvtst(conv,alpha,pnorm,toleps,xnorm,difnew,rtleps,
* ftest,gtg,peps,epsmch,gtpnew,fnew,flast,g,ipivot,n,accrcy)
implicit double precision (a-h,o-z)
logical conv,ltest
integer ipivot(n)
double precision g(n), alpha, pnorm, toleps, xnorm, difnew,
* rtleps, ftest, gtg, peps, epsmch, gtpnew, fnew, flast, one,
* cmax, t, accrcy
c
c test for convergence
c
imax = 0
cmax = 0.d0
ltest = flast - fnew .le. -5.d-1*gtpnew
do 10 i = 1,n
if (ipivot(i) .eq. 0 .or. ipivot(i) .eq. 2) go to 10
t = -ipivot(i)*g(i)
if (t .ge. 0.d0) go to 10
conv = .false.
if (ltest) go to 10
if (cmax .le. t) go to 10
cmax = t
imax = i
10 continue
if (imax .eq. 0) go to 15
ipivot(imax) = 0
flast = fnew
return
15 continue
conv = .false.
one = 1.d0
if ((alpha*pnorm .ge. toleps*(one + xnorm)
* .or. dabs(difnew) .ge. rtleps*ftest
* .or. gtg .ge. peps*ftest*ftest)
* .and. gtg .ge. 1.d-4*accrcy*ftest*ftest) return
conv = .true.
c
c for details, see gill, murray, and wright (1981, p. 308) and
c fletcher (1981, p. 116). the multiplier tests (here, testing
c the sign of the components of the gradient) may still need to
c modified to incorporate tolerances for zero.
c
return
end
c
c
subroutine crash(n,x,ipivot,low,up,ier)
implicit double precision (a-h,o-z)
double precision x(n),low(n),up(n)
integer ipivot(n)
c
c this initializes the constraint information, and ensures that the
c initial point satisfies low <= x <= up.
c the constraints are checked for consistency.
c
ier = 0
do 30 i = 1,n
if (x(i) .lt. low(i)) x(i) = low(i)
if (x(i) .gt. up(i)) x(i) = up(i)
ipivot(i) = 0
if (x(i) .eq. low(i)) ipivot(i) = -1
if (x(i) .eq. up(i)) ipivot(i) = 1
if (up(i) .eq. low(i)) ipivot(i) = 2
if (low(i) .gt. up(i)) ier = -i
30 continue
return
end
c
c the vectors sk and yk, although not in the call,
c are used (via their position in w) by the routine msolve.
c
subroutine modlnp(modet,zsol,gv,r,v,diagb,emat,
* x,g,zk,n,w,lw,niter,maxit,nfeval,nmodif,nlincg,
* upd1,yksk,gsk,yrsr,lreset,sfun,bounds,ipivot,accrcy,
* gtp,gnorm,xnorm)
implicit double precision (a-h,o-z)
integer modet,n,niter,ipivot(1)
double precision zsol(n),g(n),gv(n),r(n),v(n),diagb(n),w(lw)
double precision emat(n),zk(n),x(n),accrcy
double precision alpha,beta,delta,gsk,gtp,pr,
* qold,qnew,qtest,rhsnrm,rnorm,rz,rzold,tol,vgv,yksk,yrsr
double precision gnorm,xnorm
double precision ddot,dnrm2
logical first,upd1,lreset,bounds
external sfun
c
c this routine performs a preconditioned conjugate-gradient
c iteration in order to solve the newton equations for a search
c direction for a truncated-newton algorithm. when the value of the
c quadratic model is sufficiently reduced,
c the iteration is terminated.
c
c parameters
c
c modet - integer which controls amount of output
c zsol - computed search direction
c g - current gradient
c gv,gz1,v - scratch vectors
c r - residual
c diagb,emat - diagonal preconditoning matrix
c niter - nonlinear iteration #
c feval - value of quadratic function
c
c *************************************************************
c initialization
c *************************************************************
c
c general initialization
c
if (modet .gt. 0) write(*,800)
if (maxit .eq. 0) return
first = .true.
rhsnrm = gnorm
tol = 1.d-12
qold = 0.d0
c
c initialization for preconditioned conjugate-gradient algorithm
c
call initpc(diagb,emat,n,w,lw,modet,
* upd1,yksk,gsk,yrsr,lreset)
do 10 i = 1,n
r(i) = -g(i)
v(i) = 0.d0
zsol(i) = 0.d0
10 continue
c
c ************************************************************
c main iteration
c ************************************************************
c
do 30 k = 1,maxit
nlincg = nlincg + 1
if (modet .gt. 1) write(*,810) k
c
c cg iteration to solve system of equations
c
if (bounds) call ztime(n,r,ipivot)
call msolve(r,zk,n,w,lw,upd1,yksk,gsk,
* yrsr,lreset,first)
if (bounds) call ztime(n,zk,ipivot)
rz = ddot(n,r,1,zk,1)
if (rz/rhsnrm .lt. tol) go to 80
if (k .eq. 1) beta = 0.d0
if (k .gt. 1) beta = rz/rzold
do 20 i = 1,n
v(i) = zk(i) + beta*v(i)
20 continue
if (bounds) call ztime(n,v,ipivot)
call gtims(v,gv,n,x,g,w,lw,sfun,first,delta,accrcy,xnorm)
if (bounds) call ztime(n,gv,ipivot)
nfeval = nfeval + 1
vgv = ddot(n,v,1,gv,1)
if (vgv/rhsnrm .lt. tol) go to 50
call ndia3(n,emat,v,gv,r,vgv,modet)
c
c compute linear step length
c
alpha = rz / vgv
if (modet .ge. 1) write(*,820) alpha
c
c compute current solution and related vectors
c
call daxpy(n,alpha,v,1,zsol,1)
call daxpy(n,-alpha,gv,1,r,1)
c
c test for convergence
c
gtp = ddot(n,zsol,1,g,1)
pr = ddot(n,r,1,zsol,1)
qnew = 5.d-1 * (gtp + pr)
qtest = k * (1.d0 - qold/qnew)
if (qtest .lt. 0.d0) go to 70
qold = qnew
if (qtest .le. 5.d-1) go to 70
c
c perform cautionary test
c
if (gtp .gt. 0) go to 40
rzold = rz
30 continue
c
c terminate algorithm
c
k = k-1
go to 70
c
c truncate algorithm in case of an emergency
c
40 if (modet .ge. -1) write(*,830) k
call daxpy(n,-alpha,v,1,zsol,1)
gtp = ddot(n,zsol,1,g,1)
go to 90
50 continue
if (modet .gt. -2) write(*,840)
60 if (k .gt. 1) go to 70
call msolve(g,zsol,n,w,lw,upd1,yksk,gsk,yrsr,lreset,first)
call negvec(n,zsol)
if (bounds) call ztime(n,zsol,ipivot)
gtp = ddot(n,zsol,1,g,1)
70 continue
if (modet .ge. -1) write(*,850) k,rnorm
go to 90
80 continue
if (modet .ge. -1) write(*,860)
if (k .gt. 1) go to 70
call dcopy(n,g,1,zsol,1)
call negvec(n,zsol)
if (bounds) call ztime(n,zsol,ipivot)
gtp = ddot(n,zsol,1,g,1)
go to 70
c
c store (or restore) diagonal preconditioning
c
90 continue
call dcopy(n,emat,1,diagb,1)
return
800 format(' ',//,' entering modlnp')
810 format(' ',//,' ### iteration ',i2,' ###')
820 format(' alpha',1pd16.8)
830 format(' g(t)z positive at iteration ',i2,
* ' - truncating method',/)
840 format(' ',10x,'hessian not positive-definite')
850 format(' ',/,8x,'modlan truncated after ',i3,' iterations',
* ' rnorm = ',1pd14.6)
860 format(' preconditioning not positive-definite')
end
c
c
subroutine ndia3(n,e,v,gv,r,vgv,modet)
implicit double precision (a-h,o-z)
double precision e(n),v(n),gv(n),r(n),vgv,vr,ddot
c