-
Notifications
You must be signed in to change notification settings - Fork 370
/
powell.F
2056 lines (1984 loc) · 68.9 KB
/
powell.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
!--------------------------------------------------------------------------------------------------!
! CP2K: A general program to perform molecular dynamics simulations !
! Copyright 2000-2024 CP2K developers group <https://cp2k.org> !
! !
! SPDX-License-Identifier: GPL-2.0-or-later !
!--------------------------------------------------------------------------------------------------!
!******************************************************************************
MODULE powell
USE kinds, ONLY: dp
USE mathconstants, ONLY: twopi
#include "../base/base_uses.f90"
IMPLICIT NONE
CHARACTER(len=*), PARAMETER, PRIVATE :: moduleN = 'powell'
TYPE opt_state_type
INTEGER :: state = -1
INTEGER :: nvar = -1
INTEGER :: iprint = -1
INTEGER :: unit = -1
INTEGER :: maxfun = -1
REAL(dp) :: rhobeg = 0.0_dp, rhoend = 0.0_dp
REAL(dp), DIMENSION(:), POINTER :: w => NULL()
REAL(dp), DIMENSION(:), POINTER :: xopt => NULL()
! local variables
INTEGER :: np = -1, nh = -1, nptm = -1, nftest = -1, idz = -1, itest = -1, nf = -1, nfm = -1, nfmm = -1, &
nfsav = -1, knew = -1, kopt = -1, ksave = -1, ktemp = -1
REAL(dp) :: rhosq = 0.0_dp, recip = 0.0_dp, reciq = 0.0_dp, fbeg = 0.0_dp, &
fopt = 0.0_dp, diffa = 0.0_dp, xoptsq = 0.0_dp, &
rho = 0.0_dp, delta = 0.0_dp, dsq = 0.0_dp, dnorm = 0.0_dp, &
ratio = 0.0_dp, temp = 0.0_dp, tempq = 0.0_dp, beta = 0.0_dp, &
dx = 0.0_dp, vquad = 0.0_dp, diff = 0.0_dp, diffc = 0.0_dp, &
diffb = 0.0_dp, fsave = 0.0_dp, detrat = 0.0_dp, hdiag = 0.0_dp, &
distsq = 0.0_dp, gisq = 0.0_dp, gqsq = 0.0_dp, f = 0.0_dp, &
bstep = 0.0_dp, alpha = 0.0_dp, dstep = 0.0_dp
END TYPE opt_state_type
PRIVATE
PUBLIC :: powell_optimize, opt_state_type
!******************************************************************************
CONTAINS
!******************************************************************************
! **************************************************************************************************
!> \brief ...
!> \param n ...
!> \param x ...
!> \param optstate ...
! **************************************************************************************************
SUBROUTINE powell_optimize(n, x, optstate)
INTEGER, INTENT(IN) :: n
REAL(dp), DIMENSION(*), INTENT(INOUT) :: x
TYPE(opt_state_type), INTENT(INOUT) :: optstate
CHARACTER(len=*), PARAMETER :: routineN = 'powell_optimize'
INTEGER :: handle, npt
CALL timeset(routineN, handle)
SELECT CASE (optstate%state)
CASE (0)
npt = 2*n + 1
ALLOCATE (optstate%w((npt + 13)*(npt + n) + 3*n*(n + 3)/2))
ALLOCATE (optstate%xopt(n))
! Initialize w
optstate%w = 0.0_dp
optstate%state = 1
CALL newuoa(n, x, optstate)
CASE (1, 2)
CALL newuoa(n, x, optstate)
CASE (3)
IF (optstate%unit > 0) THEN
WRITE (optstate%unit, *) "POWELL| Exceeding maximum number of steps"
END IF
optstate%state = -1
CASE (4)
IF (optstate%unit > 0) THEN
WRITE (optstate%unit, *) "POWELL| Error in trust region"
END IF
optstate%state = -1
CASE (5)
IF (optstate%unit > 0) THEN
WRITE (optstate%unit, *) "POWELL| N out of range"
END IF
optstate%state = -1
CASE (6, 7)
optstate%state = -1
CASE (8)
x(1:n) = optstate%xopt(1:n)
DEALLOCATE (optstate%w)
DEALLOCATE (optstate%xopt)
optstate%state = -1
CASE DEFAULT
CPABORT("")
END SELECT
CALL timestop(handle)
END SUBROUTINE powell_optimize
!******************************************************************************
! **************************************************************************************************
!> \brief ...
!> \param n ...
!> \param x ...
!> \param optstate ...
! **************************************************************************************************
SUBROUTINE newuoa(n, x, optstate)
INTEGER, INTENT(IN) :: n
REAL(dp), DIMENSION(*), INTENT(INOUT) :: x
TYPE(opt_state_type), INTENT(INOUT) :: optstate
INTEGER :: ibmat, id, ifv, igq, ihq, ipq, ivl, iw, &
ixb, ixn, ixo, ixp, izmat, maxfun, &
ndim, np, npt, nptm
REAL(dp) :: rhobeg, rhoend
maxfun = optstate%maxfun
rhobeg = optstate%rhobeg
rhoend = optstate%rhoend
!
! This subroutine seeks the least value of a function of many variab
! by a trust region method that forms quadratic models by interpolat
! There can be some freedom in the interpolation conditions, which i
! taken up by minimizing the Frobenius norm of the change to the sec
! derivative of the quadratic model, beginning with a zero matrix. T
! arguments of the subroutine are as follows.
!
! N must be set to the number of variables and must be at least two.
! NPT is the number of interpolation conditions. Its value must be i
! interval [N+2,(N+1)(N+2)/2].
! Initial values of the variables must be set in X(1),X(2),...,X(N).
! will be changed to the values that give the least calculated F.
! RHOBEG and RHOEND must be set to the initial and final values of a
! region radius, so both must be positive with RHOEND<=RHOBEG. Typ
! RHOBEG should be about one tenth of the greatest expected change
! variable, and RHOEND should indicate the accuracy that is requir
! the final values of the variables.
! The value of IPRINT should be set to 0, 1, 2 or 3, which controls
! amount of printing. Specifically, there is no output if IPRINT=0
! there is output only at the return if IPRINT=1. Otherwise, each
! value of RHO is printed, with the best vector of variables so fa
! the corresponding value of the objective function. Further, each
! value of F with its variables are output if IPRINT=3.
! MAXFUN must be set to an upper bound on the number of calls of CAL
! The array W will be used for working space. Its length must be at
! (NPT+13)*(NPT+N)+3*N*(N+3)/2.
!
! SUBROUTINE CALFUN (N,X,F) must be provided by the user. It must se
! the value of the objective function for the variables X(1),X(2),..
!
! Partition the working space array, so that different parts of it c
! treated separately by the subroutine that performs the main calcul
!
np = n + 1
npt = 2*n + 1
nptm = npt - np
IF (npt < n + 2 .OR. npt > ((n + 2)*np)/2) THEN
optstate%state = 5
RETURN
END IF
ndim = npt + n
ixb = 1
ixo = ixb + n
ixn = ixo + n
ixp = ixn + n
ifv = ixp + n*npt
igq = ifv + npt
ihq = igq + n
ipq = ihq + (n*np)/2
ibmat = ipq + npt
izmat = ibmat + ndim*n
id = izmat + npt*nptm
ivl = id + n
iw = ivl + ndim
!
! The above settings provide a partition of W for subroutine NEWUOB.
! The partition requires the first NPT*(NPT+N)+5*N*(N+3)/2 elements
! W plus the space that is needed by the last array of NEWUOB.
!
CALL newuob(n, npt, x, rhobeg, rhoend, maxfun, optstate%w(ixb:), optstate%w(ixo:), &
optstate%w(ixn:), optstate%w(ixp:), optstate%w(ifv:), optstate%w(igq:), optstate%w(ihq:), &
optstate%w(ipq:), optstate%w(ibmat:), optstate%w(izmat:), ndim, optstate%w(id:), &
optstate%w(ivl:), optstate%w(iw:), optstate)
optstate%xopt(1:n) = optstate%w(ixb:ixb + n - 1) + optstate%w(ixo:ixo + n - 1)
END SUBROUTINE newuoa
!******************************************************************************
! **************************************************************************************************
!> \brief ...
!> \param n ...
!> \param npt ...
!> \param x ...
!> \param rhobeg ...
!> \param rhoend ...
!> \param maxfun ...
!> \param xbase ...
!> \param xopt ...
!> \param xnew ...
!> \param xpt ...
!> \param fval ...
!> \param gq ...
!> \param hq ...
!> \param pq ...
!> \param bmat ...
!> \param zmat ...
!> \param ndim ...
!> \param d ...
!> \param vlag ...
!> \param w ...
!> \param opt ...
! **************************************************************************************************
SUBROUTINE newuob(n, npt, x, rhobeg, rhoend, maxfun, xbase, &
xopt, xnew, xpt, fval, gq, hq, pq, bmat, zmat, ndim, d, vlag, w, opt)
INTEGER, INTENT(in) :: n, npt
REAL(dp), DIMENSION(1:n), INTENT(inout) :: x
REAL(dp), INTENT(in) :: rhobeg, rhoend
INTEGER, INTENT(in) :: maxfun
REAL(dp), DIMENSION(*), INTENT(inout) :: xbase, xopt, xnew
REAL(dp), DIMENSION(npt, *), &
INTENT(inout) :: xpt
REAL(dp), DIMENSION(*), INTENT(inout) :: fval, gq, hq, pq
INTEGER, INTENT(in) :: ndim
REAL(dp), DIMENSION(npt, *), &
INTENT(inout) :: zmat
REAL(dp), DIMENSION(ndim, *), &
INTENT(inout) :: bmat
REAL(dp), DIMENSION(*), INTENT(inout) :: d, vlag, w
TYPE(opt_state_type) :: opt
INTEGER :: i, idz, ih, ip, ipt, itemp, &
itest, j, jp, jpt, k, knew, &
kopt, ksave, ktemp, nf, nfm, &
nfmm, nfsav, nftest, nh, np, &
nptm
REAL(dp) :: alpha, beta, bstep, bsum, crvmin, delta, detrat, diff, diffa, &
diffb, diffc, distsq, dnorm, dsq, dstep, dx, f, fbeg, fopt, fsave, &
gisq, gqsq, half, hdiag, one, ratio, recip, reciq, rho, rhosq, sum, &
suma, sumb, sumz, temp, tempq, tenth, vquad, xipt, xjpt, xoptsq, zero
!
! The arguments N, NPT, X, RHOBEG, RHOEND, IPRINT and MAXFUN are ide
! to the corresponding arguments in SUBROUTINE NEWUOA.
! XBASE will hold a shift of origin that should reduce the contribut
! from rounding errors to values of the model and Lagrange functio
! XOPT will be set to the displacement from XBASE of the vector of
! variables that provides the least calculated F so far.
! XNEW will be set to the displacement from XBASE of the vector of
! variables for the current calculation of F.
! XPT will contain the interpolation point coordinates relative to X
! FVAL will hold the values of F at the interpolation points.
! GQ will hold the gradient of the quadratic model at XBASE.
! HQ will hold the explicit second derivatives of the quadratic mode
! PQ will contain the parameters of the implicit second derivatives
! the quadratic model.
! BMAT will hold the last N columns of H.
! ZMAT will hold the factorization of the leading NPT by NPT submatr
! H, this factorization being ZMAT times Diag(DZ) times ZMAT^T, wh
! the elements of DZ are plus or minus one, as specified by IDZ.
! NDIM is the first dimension of BMAT and has the value NPT+N.
! D is reserved for trial steps from XOPT.
! VLAG will contain the values of the Lagrange functions at a new po
! They are part of a product that requires VLAG to be of length ND
! The array W will be used for working space. Its length must be at
! 10*NDIM = 10*(NPT+N).
IF (opt%state == 1) THEN
! initialize all variable that will be stored
np = 0
nh = 0
nptm = 0
nftest = 0
idz = 0
itest = 0
nf = 0
nfm = 0
nfmm = 0
nfsav = 0
knew = 0
kopt = 0
ksave = 0
ktemp = 0
rhosq = 0._dp
recip = 0._dp
reciq = 0._dp
fbeg = 0._dp
fopt = 0._dp
diffa = 0._dp
xoptsq = 0._dp
rho = 0._dp
delta = 0._dp
dsq = 0._dp
dnorm = 0._dp
ratio = 0._dp
temp = 0._dp
tempq = 0._dp
beta = 0._dp
dx = 0._dp
vquad = 0._dp
diff = 0._dp
diffc = 0._dp
diffb = 0._dp
fsave = 0._dp
detrat = 0._dp
hdiag = 0._dp
distsq = 0._dp
gisq = 0._dp
gqsq = 0._dp
f = 0._dp
bstep = 0._dp
alpha = 0._dp
dstep = 0._dp
!
END IF
ipt = 0
jpt = 0
xipt = 0._dp
xjpt = 0._dp
half = 0.5_dp
one = 1.0_dp
tenth = 0.1_dp
zero = 0.0_dp
np = n + 1
nh = (n*np)/2
nptm = npt - np
nftest = MAX(maxfun, 1)
IF (opt%state == 2) GOTO 1000
!
! Set the initial elements of XPT, BMAT, HQ, PQ and ZMAT to zero.
!
DO j = 1, n
xbase(j) = x(j)
DO k = 1, npt
xpt(k, j) = zero
END DO
DO i = 1, ndim
bmat(i, j) = zero
END DO
END DO
DO ih = 1, nh
hq(ih) = zero
END DO
DO k = 1, npt
pq(k) = zero
DO j = 1, nptm
zmat(k, j) = zero
END DO
END DO
!
! Begin the initialization procedure. NF becomes one more than the n
! of function values so far. The coordinates of the displacement of
! next initial interpolation point from XBASE are set in XPT(NF,.).
!
rhosq = rhobeg*rhobeg
recip = one/rhosq
reciq = SQRT(half)/rhosq
nf = 0
50 nfm = nf
nfmm = nf - n
nf = nf + 1
IF (nfm <= 2*n) THEN
IF (nfm >= 1 .AND. nfm <= N) THEN
xpt(nf, nfm) = rhobeg
ELSE IF (nfm > n) THEN
xpt(nf, nfmm) = -rhobeg
END IF
ELSE
itemp = (nfmm - 1)/n
jpt = nfm - itemp*n - n
ipt = jpt + itemp
IF (ipt > n) THEN
itemp = jpt
jpt = ipt - n
ipt = itemp
END IF
xipt = rhobeg
IF (fval(ipt + np) < fval(ipt + 1)) xipt = -xipt
XJPT = RHOBEG
IF (fval(jpt + np) < fval(jpt + 1)) xjpt = -xjpt
xpt(nf, ipt) = xipt
xpt(nf, jpt) = xjpt
END IF
!
! Calculate the next value of F, label 70 being reached immediately
! after this calculation. The least function value so far and its in
! are required.
!
DO j = 1, n
x(j) = xpt(nf, j) + xbase(j)
END DO
GOTO 310
70 fval(nf) = f
IF (nf == 1) THEN
fbeg = f
fopt = f
kopt = 1
ELSE IF (f < fopt) THEN
fopt = f
kopt = nf
END IF
!
! Set the nonzero initial elements of BMAT and the quadratic model i
! the cases when NF is at most 2*N+1.
!
IF (NFM <= 2*N) THEN
IF (nfm >= 1 .AND. nfm <= n) THEN
gq(nfm) = (f - fbeg)/rhobeg
IF (npt < nf + n) THEN
bmat(1, nfm) = -one/rhobeg
bmat(nf, nfm) = one/rhobeg
bmat(npt + nfm, nfm) = -half*rhosq
END IF
ELSE IF (nfm > n) THEN
bmat(nf - n, nfmm) = half/rhobeg
bmat(nf, nfmm) = -half/rhobeg
zmat(1, nfmm) = -reciq - reciq
zmat(nf - n, nfmm) = reciq
zmat(nf, nfmm) = reciq
ih = (nfmm*(nfmm + 1))/2
temp = (fbeg - f)/rhobeg
hq(ih) = (gq(nfmm) - temp)/rhobeg
gq(nfmm) = half*(gq(nfmm) + temp)
END IF
!
! Set the off-diagonal second derivatives of the Lagrange functions
! the initial quadratic model.
!
ELSE
ih = (ipt*(ipt - 1))/2 + jpt
IF (xipt < zero) ipt = ipt + n
IF (xjpt < zero) jpt = jpt + n
zmat(1, nfmm) = recip
zmat(nf, nfmm) = recip
zmat(ipt + 1, nfmm) = -recip
zmat(jpt + 1, nfmm) = -recip
hq(ih) = (fbeg - fval(ipt + 1) - fval(jpt + 1) + f)/(xipt*xjpt)
END IF
IF (nf < npt) GOTO 50
!
! Begin the iterative procedure, because the initial model is comple
!
rho = rhobeg
delta = rho
idz = 1
diffa = zero
diffb = zero
itest = 0
xoptsq = zero
DO i = 1, n
xopt(i) = xpt(kopt, i)
xoptsq = xoptsq + xopt(i)**2
END DO
90 nfsav = nf
!
! Generate the next trust region step and test its length. Set KNEW
! to -1 if the purpose of the next F will be to improve the model.
!
100 knew = 0
CALL trsapp(n, npt, xopt, xpt, gq, hq, pq, delta, d, w, w(np), w(np + n), w(np + 2*n), crvmin)
dsq = zero
DO i = 1, n
dsq = dsq + d(i)**2
END DO
dnorm = MIN(delta, SQRT(dsq))
IF (dnorm < half*rho) THEN
knew = -1
delta = tenth*delta
ratio = -1.0_dp
IF (delta <= 1.5_dp*rho) delta = rho
IF (nf <= nfsav + 2) GOTO 460
temp = 0.125_dp*crvmin*rho*rho
IF (temp <= MAX(diffa, diffb, diffc)) GOTO 460
GOTO 490
END IF
!
! Shift XBASE if XOPT may be too far from XBASE. First make the chan
! to BMAT that do not depend on ZMAT.
!
120 IF (dsq <= 1.0e-3_dp*xoptsq) THEN
tempq = 0.25_dp*xoptsq
DO k = 1, npt
sum = zero
DO i = 1, n
sum = sum + xpt(k, i)*xopt(i)
END DO
temp = pq(k)*sum
sum = sum - half*xoptsq
w(npt + k) = sum
DO i = 1, n
gq(i) = gq(i) + temp*xpt(k, i)
xpt(k, i) = xpt(k, i) - half*xopt(i)
vlag(i) = bmat(k, i)
w(i) = sum*xpt(k, i) + tempq*xopt(i)
ip = npt + i
DO j = 1, i
bmat(ip, j) = bmat(ip, j) + vlag(i)*w(j) + w(i)*vlag(j)
END DO
END DO
END DO
!
! Then the revisions of BMAT that depend on ZMAT are calculated.
!
DO k = 1, nptm
sumz = zero
DO i = 1, npt
sumz = sumz + zmat(i, k)
w(i) = w(npt + i)*zmat(i, k)
END DO
DO j = 1, n
sum = tempq*sumz*xopt(j)
DO i = 1, npt
sum = sum + w(i)*xpt(i, j)
vlag(j) = sum
IF (k < idz) sum = -sum
END DO
DO i = 1, npt
bmat(i, j) = bmat(i, j) + sum*zmat(i, k)
END DO
END DO
DO i = 1, n
ip = i + npt
temp = vlag(i)
IF (k < idz) temp = -temp
DO j = 1, i
bmat(ip, j) = bmat(ip, j) + temp*vlag(j)
END DO
END DO
END DO
!
! The following instructions complete the shift of XBASE, including
! the changes to the parameters of the quadratic model.
!
ih = 0
DO j = 1, n
w(j) = zero
DO k = 1, npt
w(j) = w(j) + pq(k)*xpt(k, j)
xpt(k, j) = xpt(k, j) - half*xopt(j)
END DO
DO i = 1, j
ih = ih + 1
IF (i < j) gq(j) = gq(j) + hq(ih)*xopt(i)
gq(i) = gq(i) + hq(ih)*xopt(j)
hq(ih) = hq(ih) + w(i)*xopt(j) + xopt(i)*w(j)
bmat(npt + i, j) = bmat(npt + j, i)
END DO
END DO
DO j = 1, n
xbase(j) = xbase(j) + xopt(j)
xopt(j) = zero
END DO
xoptsq = zero
END IF
!
! Pick the model step if KNEW is positive. A different choice of D
! may be made later, if the choice of D by BIGLAG causes substantial
! cancellation in DENOM.
!
IF (knew > 0) THEN
CALL biglag(n, npt, xopt, xpt, bmat, zmat, idz, ndim, knew, dstep, &
d, alpha, vlag, vlag(npt + 1), w, w(np), w(np + n))
END IF
!
! Calculate VLAG and BETA for the current choice of D. The first NPT
! components of W_check will be held in W.
!
DO k = 1, npt
suma = zero
sumb = zero
sum = zero
DO j = 1, n
suma = suma + xpt(k, j)*d(j)
sumb = sumb + xpt(k, j)*xopt(j)
sum = sum + bmat(k, j)*d(j)
END DO
w(k) = suma*(half*suma + sumb)
vlag(k) = sum
END DO
beta = zero
DO k = 1, nptm
sum = zero
DO i = 1, npt
sum = sum + zmat(i, k)*w(i)
END DO
IF (k < idz) THEN
beta = beta + sum*sum
sum = -sum
ELSE
beta = beta - sum*sum
END IF
DO i = 1, npt
vlag(i) = vlag(i) + sum*zmat(i, k)
END DO
END DO
bsum = zero
dx = zero
DO j = 1, n
sum = zero
DO i = 1, npt
sum = sum + w(i)*bmat(i, j)
END DO
bsum = bsum + sum*d(j)
jp = npt + j
DO k = 1, n
sum = sum + bmat(jp, k)*d(k)
END DO
vlag(jp) = sum
bsum = bsum + sum*d(j)
dx = dx + d(j)*xopt(j)
END DO
beta = dx*dx + dsq*(xoptsq + dx + dx + half*dsq) + beta - bsum
vlag(kopt) = vlag(kopt) + one
!
! If KNEW is positive and if the cancellation in DENOM is unacceptab
! then BIGDEN calculates an alternative model step, XNEW being used
! working space.
!
IF (knew > 0) THEN
temp = one + alpha*beta/vlag(knew)**2
IF (ABS(temp) <= 0.8_dp) THEN
CALL bigden(n, npt, xopt, xpt, bmat, zmat, idz, ndim, kopt, &
knew, d, w, vlag, beta, xnew, w(ndim + 1), w(6*ndim + 1))
END IF
END IF
!
! Calculate the next value of the objective function.
!
290 DO i = 1, n
xnew(i) = xopt(i) + d(i)
x(i) = xbase(i) + xnew(i)
END DO
nf = nf + 1
310 IF (nf > nftest) THEN
! return to many steps
nf = nf - 1
opt%state = 3
CALL get_state
GOTO 530
END IF
CALL get_state
opt%state = 2
RETURN
1000 CONTINUE
CALL set_state
IF (nf <= npt) GOTO 70
IF (knew == -1) THEN
opt%state = 6
CALL get_state
GOTO 530
END IF
!
! Use the quadratic model to predict the change in F due to the step
! and set DIFF to the error of this prediction.
!
vquad = zero
ih = 0
DO j = 1, n
vquad = vquad + d(j)*gq(j)
DO i = 1, j
ih = ih + 1
temp = d(i)*xnew(j) + d(j)*xopt(i)
IF (i == j) temp = half*temp
vquad = vquad + temp*hq(ih)
END DO
END DO
DO k = 1, npt
vquad = vquad + pq(k)*w(k)
END DO
diff = f - fopt - vquad
diffc = diffb
diffb = diffa
diffa = ABS(diff)
IF (dnorm > rho) nfsav = nf
!
! Update FOPT and XOPT if the new F is the least value of the object
! function so far. The branch when KNEW is positive occurs if D is n
! a trust region step.
!
fsave = fopt
IF (f < fopt) THEN
fopt = f
xoptsq = zero
DO i = 1, n
xopt(i) = xnew(i)
xoptsq = xoptsq + xopt(i)**2
END DO
END IF
ksave = knew
IF (knew > 0) GOTO 410
!
! Pick the next value of DELTA after a trust region step.
!
IF (vquad >= zero) THEN
! Return because a trust region step has failed to reduce Q
opt%state = 4
CALL get_state
GOTO 530
END IF
ratio = (f - fsave)/vquad
IF (ratio <= tenth) THEN
delta = half*dnorm
ELSE IF (ratio <= 0.7_dp) THEN
delta = MAX(half*delta, dnorm)
ELSE
delta = MAX(half*delta, dnorm + dnorm)
END IF
IF (delta <= 1.5_dp*rho) delta = rho
!
! Set KNEW to the index of the next interpolation point to be delete
!
rhosq = MAX(tenth*delta, rho)**2
ktemp = 0
detrat = zero
IF (f >= fsave) THEN
ktemp = kopt
detrat = one
END IF
DO k = 1, npt
hdiag = zero
DO j = 1, nptm
temp = one
IF (j < idz) temp = -one
hdiag = hdiag + temp*zmat(k, j)**2
END DO
temp = ABS(beta*hdiag + vlag(k)**2)
distsq = zero
DO j = 1, n
distsq = distsq + (xpt(k, j) - xopt(j))**2
END DO
IF (distsq > rhosq) temp = temp*(distsq/rhosq)**3
IF (temp > detrat .AND. k /= ktemp) THEN
detrat = temp
knew = k
END IF
END DO
IF (knew == 0) GOTO 460
!
! Update BMAT, ZMAT and IDZ, so that the KNEW-th interpolation point
! can be moved. Begin the updating of the quadratic model, starting
! with the explicit second derivative term.
!
410 CALL update(n, npt, bmat, zmat, idz, ndim, vlag, beta, knew, w)
fval(knew) = f
ih = 0
DO i = 1, n
temp = pq(knew)*xpt(knew, i)
DO j = 1, i
ih = ih + 1
hq(ih) = hq(ih) + temp*xpt(knew, j)
END DO
END DO
pq(knew) = zero
!
! Update the other second derivative parameters, and then the gradie
! vector of the model. Also include the new interpolation point.
!
DO j = 1, nptm
temp = diff*zmat(knew, j)
IF (j < idz) temp = -temp
DO k = 1, npt
pq(k) = pq(k) + temp*zmat(k, j)
END DO
END DO
gqsq = zero
DO i = 1, n
gq(i) = gq(i) + diff*bmat(knew, i)
gqsq = gqsq + gq(i)**2
xpt(knew, i) = xnew(i)
END DO
!
! If a trust region step makes a small change to the objective funct
! then calculate the gradient of the least Frobenius norm interpolan
! XBASE, and store it in W, using VLAG for a vector of right hand si
!
IF (ksave == 0 .AND. delta == rho) THEN
IF (ABS(ratio) > 1.0e-2_dp) THEN
itest = 0
ELSE
DO k = 1, npt
vlag(k) = fval(k) - fval(kopt)
END DO
gisq = zero
DO i = 1, n
sum = zero
DO k = 1, npt
sum = sum + bmat(k, i)*vlag(k)
END DO
gisq = gisq + sum*sum
w(i) = sum
END DO
!
! Test whether to replace the new quadratic model by the least Frobe
! norm interpolant, making the replacement if the test is satisfied.
!
itest = itest + 1
IF (gqsq < 1.0e2_dp*gisq) itest = 0
IF (itest >= 3) THEN
DO i = 1, n
gq(i) = w(i)
END DO
DO ih = 1, nh
hq(ih) = zero
END DO
DO j = 1, nptm
w(j) = zero
DO k = 1, npt
w(j) = w(j) + vlag(k)*zmat(k, j)
END DO
IF (j < idz) w(j) = -w(j)
END DO
DO k = 1, npt
pq(k) = zero
DO j = 1, nptm
pq(k) = pq(k) + zmat(k, j)*w(j)
END DO
END DO
itest = 0
END IF
END IF
END IF
IF (f < fsave) kopt = knew
!
! If a trust region step has provided a sufficient decrease in F, th
! branch for another trust region calculation. The case KSAVE>0 occu
! when the new function value was calculated by a model step.
!
IF (f <= fsave + tenth*vquad) GOTO 100
IF (ksave > 0) GOTO 100
!
! Alternatively, find out if the interpolation points are close enou
! to the best point so far.
!
knew = 0
460 distsq = 4.0_dp*delta*delta
DO k = 1, npt
sum = zero
DO j = 1, n
sum = sum + (xpt(k, j) - xopt(j))**2
END DO
IF (sum > distsq) THEN
knew = k
distsq = sum
END IF
END DO
!
! If KNEW is positive, then set DSTEP, and branch back for the next
! iteration, which will generate a "model step".
!
IF (knew > 0) THEN
dstep = MAX(MIN(tenth*SQRT(distsq), half*delta), rho)
dsq = dstep*dstep
GOTO 120
END IF
IF (ratio > zero) GOTO 100
IF (MAX(delta, dnorm) > rho) GOTO 100
!
! The calculations with the current value of RHO are complete. Pick
! next values of RHO and DELTA.
!
490 IF (rho > rhoend) THEN
delta = half*rho
ratio = rho/rhoend
IF (ratio <= 16.0_dp) THEN
rho = rhoend
ELSE IF (ratio <= 250.0_dp) THEN
rho = SQRT(ratio)*rhoend
ELSE
rho = tenth*rho
END IF
delta = MAX(delta, rho)
GOTO 90
END IF
!
! Return from the calculation, after another Newton-Raphson step, if
! it is too short to have been tried before.
!
IF (knew == -1) GOTO 290
opt%state = 7
CALL get_state
530 IF (fopt <= f) THEN
DO i = 1, n
x(i) = xbase(i) + xopt(i)
END DO
f = fopt
END IF
CALL get_state
!******************************************************************************
CONTAINS
!******************************************************************************
! **************************************************************************************************
!> \brief ...
! **************************************************************************************************
SUBROUTINE get_state()
opt%np = np
opt%nh = nh
opt%nptm = nptm
opt%nftest = nftest
opt%idz = idz
opt%itest = itest
opt%nf = nf
opt%nfm = nfm
opt%nfmm = nfmm
opt%nfsav = nfsav
opt%knew = knew
opt%kopt = kopt
opt%ksave = ksave
opt%ktemp = ktemp
opt%rhosq = rhosq
opt%recip = recip
opt%reciq = reciq
opt%fbeg = fbeg
opt%fopt = fopt
opt%diffa = diffa
opt%xoptsq = xoptsq
opt%rho = rho
opt%delta = delta
opt%dsq = dsq
opt%dnorm = dnorm
opt%ratio = ratio
opt%temp = temp
opt%tempq = tempq
opt%beta = beta
opt%dx = dx
opt%vquad = vquad
opt%diff = diff
opt%diffc = diffc
opt%diffb = diffb
opt%fsave = fsave
opt%detrat = detrat
opt%hdiag = hdiag
opt%distsq = distsq
opt%gisq = gisq
opt%gqsq = gqsq
opt%f = f
opt%bstep = bstep
opt%alpha = alpha
opt%dstep = dstep
END SUBROUTINE get_state
!******************************************************************************
! **************************************************************************************************
!> \brief ...
! **************************************************************************************************
SUBROUTINE set_state()
np = opt%np
nh = opt%nh
nptm = opt%nptm
nftest = opt%nftest
idz = opt%idz
itest = opt%itest
nf = opt%nf
nfm = opt%nfm
nfmm = opt%nfmm
nfsav = opt%nfsav
knew = opt%knew
kopt = opt%kopt
ksave = opt%ksave
ktemp = opt%ktemp
rhosq = opt%rhosq
recip = opt%recip
reciq = opt%reciq
fbeg = opt%fbeg
fopt = opt%fopt
diffa = opt%diffa
xoptsq = opt%xoptsq
rho = opt%rho
delta = opt%delta
dsq = opt%dsq
dnorm = opt%dnorm
ratio = opt%ratio
temp = opt%temp
tempq = opt%tempq
beta = opt%beta
dx = opt%dx
vquad = opt%vquad
diff = opt%diff
diffc = opt%diffc
diffb = opt%diffb
fsave = opt%fsave
detrat = opt%detrat