-
Notifications
You must be signed in to change notification settings - Fork 9
/
Copy pathminimize.f90
2818 lines (2553 loc) · 90.8 KB
/
minimize.f90
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
module minimize
!This module contains all the minimization routines: amoeba and bobyq
use nrtype
REAL(DP) :: eps=epsilon(1.0_DP)
contains
SUBROUTINE amoeba(p,y,ftol,func,iter,ITMAX_USER,global_fval)
USE utilities, ONLY : assert_eq,imaxloc,iminloc,nrerror,swap
IMPLICIT NONE
REAL(DP), DIMENSION(:,:), INTENT(INOUT) :: p ! vertices. If we have n vertices, then we must be
! in n-1 dimensional space (we need one extra vertex
! than dimensions. For each row, the n-1 vector
! specifies the vertex
REAL(DP), DIMENSION(:), INTENT(INOUT) :: y ! "func" evaluated at the n vertices provided in "p"
REAL(DP), INTENT(IN) :: ftol ! tolerance for amoeba
INTERFACE
REAL(DP) FUNCTION func(x)
use genericParams
IMPLICIT NONE
REAL(DP), DIMENSION(:), INTENT(IN) :: x
! REAL(DP) :: func
END FUNCTION func
END INTERFACE
INTEGER(I4B), INTENT(OUT) :: iter ! the number of iterations required
REAL(DP), INTENT(IN), OPTIONAL :: global_fval ! the best global objective value so far
INTEGER(I4B), INTENT(IN), OPTIONAL :: ITMAX_USER
INTEGER(I4B) :: ITMAX=5000
REAL(dp), PARAMETER :: TINY=1.0e-10
INTEGER(I4B) :: ihi,ndim
REAL(dp), DIMENSION(size(p,2)) :: psum
if(PRESENT(ITMAX_USER)) THEN
ITMAX=ITMAX_USER
else
ITMAX=5000
ENDIF
call amoeba_private
CONTAINS
!BL
SUBROUTINE amoeba_private
IMPLICIT NONE
INTEGER(I4B) :: i,ilo,inhi,j
REAL(dp) :: rtol,ysave,ytry,ytmp
REAL(DP) :: ylo_pre,improv
INTEGER :: iter_pre
ndim=assert_eq(size(p,2),size(p,1)-1,size(y)-1,'amoeba')
iter=0
psum(:)=sum(p(:,:),dim=1)
ihi=imaxloc(y(:))
ylo_pre=y(ihi)
iter_pre=0
do
ilo=iminloc(y(:))
ihi=imaxloc(y(:))
ytmp=y(ihi)
y(ihi)=y(ilo)
inhi=imaxloc(y(:))
y(ihi)=ytmp
rtol=2.0_DP*dabs(y(ihi)-y(ilo))/(dabs(y(ihi))+dabs(y(ilo))+TINY)
if(mod(iter,100)==0) THEN
improv=ylo_pre-y(ilo)
ylo_pre=y(ilo)
WRITE (*,*) 'in ameoba, iteration =', iter
100 FORMAT(i2,2i10,5i12)
WRITE(*,*) ' ya(ilo) ya(ihi) rtol improvement'
WRITE(6,101) y(ilo),y(ihi),rtol,improv
101 FORMAT(5f12.6)
IF(PRESENT(global_fval) .AND. iter>=400 .AND. improv>0.0_DP) THEN
IF((y(ilo)-improv*DBLE(ITMAX-iter)/DBLE(iter-iter_pre)) > 1.02_DP*global_fval) THEN
PRINT*,'NOT ENOUGH IMPROVEMENT'
RETURN
ENDIF
ENDIF
iter_pre=iter
ENDIF
ilo=iminloc(y(:))
ihi=imaxloc(y(:))
ytmp=y(ihi)
y(ihi)=y(ilo)
inhi=imaxloc(y(:))
y(ihi)=ytmp
rtol=2.0_dp*abs(y(ihi)-y(ilo))/(abs(y(ihi))+abs(y(ilo))+TINY)
if (rtol < ftol) then
call swap(y(1),y(ilo))
call swap(p(1,:),p(ilo,:))
RETURN
end if
!if (iter >= ITMAX) call nrerror('ITMAX exceeded in amoeba')
if (iter >= ITMAX) then
write(*,*) 'reached max number of evaluations in amoeba; returning best result so far'
call swap(y(1),y(ilo))
call swap(p(1,:),p(ilo,:))
RETURN
end if
ytry=amotry(-1.0_dp)
iter=iter+1
if (ytry <= y(ilo)) then
ytry=amotry(2.0_dp)
iter=iter+1
else if (ytry >= y(inhi)) then
ysave=y(ihi)
ytry=amotry(0.5_dp)
iter=iter+1
if (ytry >= ysave) then
p(:,:)=0.5_dp*(p(:,:)+spread(p(ilo,:),1,size(p,1)))
do i=1,ndim+1
if (i /= ilo) y(i)=func(p(i,:))
end do
iter=iter+ndim
psum(:)=sum(p(:,:),dim=1)
end if
end if
end do
END SUBROUTINE amoeba_private
!BL
FUNCTION amotry(fac)
IMPLICIT NONE
REAL(dp), INTENT(IN) :: fac
REAL(dp) :: amotry
REAL(dp) :: fac1,fac2,ytry
REAL(dp), DIMENSION(size(p,2)) :: ptry
fac1=(1.0_dp-fac)/ndim
fac2=fac1-fac
ptry(:)=psum(:)*fac1-p(ihi,:)*fac2
ytry=func(ptry)
if (ytry < y(ihi)) then
y(ihi)=ytry
psum(:)=psum(:)-p(ihi,:)+ptry(:)
p(ihi,:)=ptry(:)
end if
amotry=ytry
END FUNCTION amotry
END SUBROUTINE amoeba
SUBROUTINE amebsa(p,y,pb,yb,ftol,func,iter,temptr)
USE utilities, ONLY : assert_eq,imaxloc,iminloc,swap
IMPLICIT NONE
INTEGER(I4B), INTENT(INOUT) :: iter
REAL(DP), INTENT(INOUT) :: yb
REAL(DP), INTENT(IN) :: ftol,temptr
REAL(DP), DIMENSION(:), INTENT(INOUT) :: y,pb
REAL(DP), DIMENSION(:,:), INTENT(INOUT) :: p
INTERFACE
FUNCTION func(x)
USE nrtype
IMPLICIT NONE
REAL(DP), DIMENSION(:), INTENT(IN) :: x
REAL(DP) :: func
END FUNCTION func
END INTERFACE
INTEGER(I4B), PARAMETER :: NMAX=200
INTEGER(I4B) :: ihi,ndim
REAL(DP) :: yhi
REAL(DP), DIMENSION(size(p,2)) :: psum
call amebsa_private
CONTAINS
!BL
SUBROUTINE amebsa_private
INTEGER(I4B) :: i,ilo,inhi
REAL(DP) :: rtol,ylo,ynhi,ysave,ytry
REAL(DP), DIMENSION(size(y)) :: yt,harvest
ndim=assert_eq(size(p,2),size(p,1)-1,size(y)-1,size(pb),'amebsa')
psum(:)=sum(p(:,:),dim=1)
do
call random_number(harvest)
yt(:)=y(:)-temptr*log(harvest)
ilo=iminloc(yt(:))
ylo=yt(ilo)
ihi=imaxloc(yt(:))
yhi=yt(ihi)
yt(ihi)=ylo
inhi=imaxloc(yt(:))
ynhi=yt(inhi)
rtol=2.0_dp*abs(yhi-ylo)/(abs(yhi)+abs(ylo))
if (rtol < ftol .or. iter < 0) then
call swap(y(1),y(ilo))
call swap(p(1,:),p(ilo,:))
RETURN
end if
ytry=amotsa(-1.0_dp)
iter=iter-1
if (ytry <= ylo) then
ytry=amotsa(2.0_dp)
iter=iter-1
else if (ytry >= ynhi) then
ysave=yhi
ytry=amotsa(0.5_dp)
iter=iter-1
if (ytry >= ysave) then
p(:,:)=0.5_dp*(p(:,:)+spread(p(ilo,:),1,size(p,1)))
do i=1,ndim+1
if (i /= ilo) y(i)=func(p(i,:))
end do
iter=iter-ndim
psum(:)=sum(p(:,:),dim=1)
end if
end if
end do
END SUBROUTINE amebsa_private
!BL
FUNCTION amotsa(fac)
IMPLICIT NONE
REAL(DP), INTENT(IN) :: fac
REAL(DP) :: amotsa
REAL(DP) :: fac1,fac2,yflu,ytry,harv
REAL(DP), DIMENSION(size(p,2)) :: ptry
fac1=(1.0_dp-fac)/ndim
fac2=fac1-fac
ptry(:)=psum(:)*fac1-p(ihi,:)*fac2
ytry=func(ptry)
if (ytry <= yb) then
pb(:)=ptry(:)
yb=ytry
end if
call random_number(harv)
yflu=ytry+temptr*log(harv)
if (yflu < yhi) then
y(ihi)=ytry
yhi=yflu
psum(:)=psum(:)-p(ihi,:)+ptry(:)
p(ihi,:)=ptry(:)
end if
amotsa=yflu
END FUNCTION amotsa
END SUBROUTINE amebsa
! Important Notice:
! These algorithms are modifications and based on the Software BOBYQA, authored by M. J. D. Powell,
! to minimize sum of squares with bound constraints by taking advantage of the problem structure.
! i.e. Min F(x) := Sum_{i=1}^{mv} v_err_i(x)^2, s.t. xl <= x <= xu, x \in R^n,
! where v_err(x) : R^n \to R^{mv} is a vector function.
! This subroutine seeks the least value of sum of the squres of the components of v_err(x)
! by combing trust region method and Levenberg-Marquardt method
!
! References:
!
! 1. M. J. D. Powell, The NEWUOA software for unconstrained optimization without derivatives,
! DAMTP 2004/ NA 05
! 2. M. J. D. Powell, The BOBYQA algorithm for bound constrained optimization without derivatives,
! DAMTP 2009/ NA 06
! 3. H. Zhang, A. R. CONN, AND K. SCHEINBERG, A DERIVATIVE-FREE ALGORITHM FOR THE LEAST-SQUARE
! MINIMIZATION, technical report, 2009
!
! -----------------------------------------------------------------
! | This program is free software; you can redistribute it and/or |
! |modify it under the terms of the GNU General Public License as |
! |published by the Free Software Foundation; either version 2 of |
! |the License, or (at your option) any later version. |
! |This program is distributed in the hope that it will be useful, |
! |but WITHOUT ANY WARRANTY; without even the implied warranty of |
! |MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the |
! |GNU General Public License for more details. |
! | |
! |You should have received a copy of the GNU General Publi! |
! |License along with this program; if not, write to the Free |
! |Software Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, |
! |MA 02110-1301 USA |
! -----------------------------------------------------------------|
SUBROUTINE BOBYQA_H (N,NPT,X,XL,XU,RHOBEG,RHOEND,IPRINT, MAXFUN,W,mv)
IMPLICIT double precision (A-H,O-Z)
DIMENSION X(*),XL(*),XU(*),W(*)
! N must be set to the number of variables and must be at least two.
! mv must be set to the lengh of the vector function v_err(x): R^n \to R^{mv}.
! The maximum number variables in this codes: nmax = 100
! The maximum lengh of the vector function v_err(x): mmax = 400
! If n > 100 or m > 400, the parameter nmax and mmax need to be creased in
! subroutine BOBYQB_H, PRELIM_H, RESCUE_H and TRSBOX_H
! NPT is the number of interpolation conditions. Its value must be in the
! interval [N+2,2N+1]. Recommended: NPT = 2*N+1
!
! Initial values of the variables must be set in X(1),X(2),...,X(N). They
! will be changed to the values that give the least calculated F= Sum_{i=1}^{mv} v_err_i(x)^2..
!
! For I=1,2,...,N, XL(I) and XU(I) must provide the lower and upper
! bounds, respectively, on X(I). The construction of quadratic models
! requires XL(I) to be strictly less than XU(I) for each I. Further,
! the contribution to a model from changes to the I-th variable is
! damaged severely by rounding errors if XU(I)-XL(I) is too small.
! RHOBEG and RHOEND must be set to the initial and final values of a trust
! region radius, so both must be positive with RHOEND no greater than
! RHOBEG. Typically, RHOBEG should be about one tenth of the greatest
! expected change to a variable, while RHOEND should indicate the
! accuracy that is required in the final values of the variables. An
! error return occurs if any of the differences XU(I)-XL(I), I=1,...,N,
! is less than 2*RHOBEG. Default: RHOBEG = 1.0, RHOEND = 10^{-8}
! The value of IPRINT should be set to 0, 1, 2 or 3, which controls the
! amount of printing. Specifically, there is no output if IPRINT=0 and
! there is output only at the return if IPRINT=1. Otherwise, each new
! value of RHO is printed, with the best vector of variables so far and
! the corresponding value of the objective function. Further, each new
! 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 subroutine
! dfovec(n, mv, x, v_err) which provides the values of the vector function v_err(x).
! Here: n, mv, x \in R^n are input, v_err \in R^{mv} are output.
! Default: MAXFUN= 400(n+1), i.e 400 (simplex) gradients for reasonable accuracy.
! MAXFUN= infinity, to let the algorithm explore the lowest function value
! as much as it could.
! The array W will be used for working space. Its length must be at least
! (NPT+5)*(NPT+N)+3*N*(N+5)/2.
! SUBROUTINE dfovec(n, mv, x, v_err) must be provided by the user.
! It must provide the values of the vector function v_err(x) : R^n to R^{mv}
! at the variables X(1),X(2),...,X(N), which are generated automatically in
! a way that satisfies the bounds given in XL and XU.
! Return if the value of NPT is unacceptable.
NP=N+1
IF (NPT < N+2 .OR. NPT > 2*N+1) THEN
PRINT 10
10 FORMAT (/4X,'Return from NEWUOA because NPT is not in','[N+2, 2N+1]')
GO TO 40
END IF
! Partition the working space array, so that different parts of it can
! be treated separately during the calculation of BOBYQB. The partition
! requires the first (NPT+2)*(NPT+N)+3*N*(N+5)/2 elements of W plus the
! space that is taken by the last array in the argument list of BOBYQB.
NDIM=NPT+N
IXB=1
IXP=IXB+N
IFV=IXP+N*NPT
IXO=IFV+NPT
IGO=IXO+N
IHQ=IGO+N
IPQ=IHQ+(N*NP)/2
IBMAT=IPQ+NPT
IZMAT=IBMAT+NDIM*N
ISL=IZMAT+NPT*(NPT-NP)
ISU=ISL+N
IXN=ISU+N
IXA=IXN+N
ID=IXA+N
IVL=ID+N
IW=IVL+NDIM
! Return if there is insufficient space between the bounds. Modify the
! initial X if necessary in order to avoid conflicts between the bounds
! and the construction of the first quadratic model. The lower and upper
! bounds on moves from the updated X are set now, in the ISL and ISU
! partitions of W, in order to provide useful and exact information about
! components of X that become within distance RHOBEG from their bounds.
ZERO=0.0_DP
DO 30 J=1,N
TEMP=XU(J)-XL(J)
IF (TEMP < RHOBEG+RHOBEG) THEN
PRINT 20
20 FORMAT (/4X,'Return from BOBYQA_H because one of the',&
& ' differences XU(I)-XL(I)'/6X,' is less than 2*RHOBEG.')
GO TO 40
END IF
JSL=ISL+J-1
JSU=JSL+N
W(JSL)=XL(J)-X(J)
W(JSU)=XU(J)-X(J)
if (w(jsl) >= -rhobeg) then
if (w(jsl) >= zero) then
x(j)=xl(j)
w(jsl)=zero
w(jsu)=temp
else
x(j)=xl(j)+rhobeg
w(jsl)=-rhobeg
w(jsu)=dmax1(xu(j)-x(j),rhobeg)
end if
else if (w(jsu) <= rhobeg) then
if (w(jsu) <= zero) then
x(j)=xu(j)
w(jsl)=-temp
w(jsu)=zero
else
x(j)=xu(j)-rhobeg
w(jsl)=dmin1(xl(j)-x(j),-rhobeg)
w(jsu)=rhobeg
end if
end if
30 CONTINUE
! Make the call of BOBYQB_H.
CALL BOBYQB_H (N,NPT,X,XL,XU,RHOBEG,RHOEND,IPRINT,MAXFUN,W(IXB),&
W(IXP),W(IFV),W(IXO),W(IGO),W(IHQ),W(IPQ),W(IBMAT),W(IZMAT),&
NDIM,W(ISL),W(ISU),W(IXN),W(IXA),W(ID),W(IVL),W(IW),mv)
40 RETURN
END subroutine
!
! Important Notice:
! This NEWUOB_H are modifications and based on the subroutine BOBYQB in the software
! BOBYQA, authored by M. J. D. Powell.
!
SUBROUTINE BOBYQB_H (N,NPT,X,XL,XU,RHOBEG,RHOEND,IPRINT,&
MAXFUN,XBASE,XPT,FVAL,XOPT,GOPT,HQ,PQ,BMAT,ZMAT,NDIM,&
SL,SU,XNEW,XALT,D,VLAG,W,mv)
use OBJECTIVE, only: dfovec, getPenalty
IMPLICIT double precision (A-H,O-Z)
DIMENSION X(*),XL(*),XU(*),XBASE(*),XPT(NPT,*),FVAL(*),&
XOPT(*),GOPT(*),HQ(*),PQ(*),BMAT(NDIM,*),ZMAT(NPT,*),&
SL(*),SU(*),XNEW(*),XALT(*),D(*),VLAG(*),W(*)
integer nmax, mmax, nptmax
parameter (nmax = 100, mmax=3000, nptmax=2*nmax+1)
dimension GQV(mmax,nmax),HQV(mmax,(nmax+1)*nmax/2),&
& PQV(mmax,nptmax), PQVold(mmax),&
& FVAL_V(mmax,nptmax),v_sum(nptmax),&
& WV(mmax,nmax), v_sumpq(mmax), &
& v_err(mmax), v_beg(mmax), v_temp(mmax),&
& v_diff(mmax), v_vquad(mmax),&
& HD1(nmax)
logical model_update, opt_update
integer v_itest(mmax)
! external dfovec
if (n>nmax) then
print *, "in bobyqb_h.f increase the dimension &
& nmax to be at least", n
stop
endif
if (mv>mmax) then
print *, "in bobyqb_h.f increase the dimension &
& mmax to be at least", mv
stop
endif
! Set some constants.
model_update = .true.
opt_update = .true.
HALF=0.5_DP
ONE=1.0_DP
TEN=10.0_DP
TENTH=0.1_DP
TWO=2.0_DP
ZERO=0.0_DP
NP=N+1
NPTM=NPT-NP
NH=(N*NP)/2
! The call of PRELIM sets the elements of XBASE, XPT, FVAL, GOPT, HQ, PQ,
! BMAT and ZMAT for the first iteration, with the corresponding values of
! of NF and KOPT, which are the number of calls of CALFUN so far and the
! index of the interpolation point at the trust region centre. Then the
! initial XOPT is set too. The branch to label 720 occurs if MAXFUN is
! less than NPT. GOPT will be updated if KOPT is different from KBASE.
CALL PRELIM_H (N,NPT,X,XL,XU,RHOBEG,IPRINT,MAXFUN,XBASE,XPT,&
FVAL,GOPT,HQ,PQ,BMAT,ZMAT,NDIM,SL,SU,NF,KOPT,&
mv,HQV,GQV,PQV,FVAL_V,v_err,v_beg,FBEG)
XOPTSQ=ZERO
DO 10 I=1,N
XOPT(I)=XPT(KOPT,I)
10 XOPTSQ=XOPTSQ+XOPT(I)**2
FSAVE=FVAL(1)
IF (NF < NPT) THEN
IF (IPRINT > 0) PRINT 390
GOTO 720
END IF
KBASE=1
! Complete the settings that are required for the iterative procedure.
RHO=RHOBEG
DELTA=RHO
NRESC=NF
NTRITS=0
DIFFA=ZERO
DIFFB=ZERO
do 15 m1=1,mv
15 v_itest(m1)=0
NFSAV=NF
! Update GOPT if necessary before the first iteration and after each
! call of RESCUE that makes a call of CALFUN.
20 IF (KOPT .NE. KBASE) THEN
model_update = .true.
IH=0
DO 30 J=1,N
DO 30 I=1,J
IH=IH+1
IF (I < J) then
do 24 m1=1,mv
GQV(m1,J)=GQV(m1,J)+HQV(m1,IH)*XOPT(I)
24 continue
endif
do 26 m1=1,mv
GQV(m1,I)=GQV(m1,I)+HQV(m1,IH)*XOPT(J)
26 continue
30 continue
IF (NF > NPT) THEN
DO 50 K=1,NPT
TEMP=ZERO
DO 40 J=1,N
40 TEMP=TEMP+XPT(K,J)*XOPT(J)
do 45 m1=1,mv
t=PQV(m1,k)*TEMP
do 45 i=1,n
45 GQV(m1,i)=GQV(m1,i)+t*XPT(k,i)
50 continue
END IF
END IF
! Generate the next point in the trust region that provides a small value
! of the quadratic model subject to the constraints on the variables.
! The integer NTRITS is set to the number "trust region" iterations that
! have occurred since the last "alternative" iteration. If the length
! of XNEW-XOPT is less than HALF*RHO, however, then there is a branch to
! label 650 or 680 with NTRITS=-1, instead of calculating F at XNEW.
60 CALL TRSBOX_H (N,NPT,XPT,XOPT,GOPT,HQ,PQ,SL,SU,DELTA,XNEW,D,&
W,W(NP),W(NP+N),W(NP+2*N),W(NP+3*N),DSQ,CRVMIN,&
mv,HQV,GQV,PQV,FVAL_V,KOPT,HD1,v_temp,vquad1,&
model_update,opt_update)
DNORM=DMIN1(DELTA,DSQRT(DSQ))
IF (DNORM < HALF*RHO) THEN
NTRITS=-1
DISTSQ=(TEN*RHO)**2
IF (NF <= NFSAV+2) GOTO 650
! The following choice between labels 650 and 680 depends on whether or
! not our work with the current RHO seems to be complete. Either RHO is
! decreased or termination occurs if the errors in the quadratic model at
! the last three interpolation points compare favourably with predictions
! of likely improvements to the model within distance HALF*RHO of XOPT.
ERRBIG=DMAX1(DIFFA,DIFFB,DIFFC)
FRHOSQ=0.125_DP*RHO*RHO
IF (CRVMIN > ZERO .AND. ERRBIG > FRHOSQ*CRVMIN) GOTO 650
BDTOL=ERRBIG/RHO
DO 80 J=1,N
BDTEST=BDTOL
IF (dabs(XNEW(J)-SL(J)) <= eps) BDTEST=W(J)
IF (dabs(XNEW(J)-SU(J)) <= eps) BDTEST=-W(J)
IF (BDTEST < BDTOL) THEN
CURV=HQ((J+J*J)/2)
BDTEST=BDTEST+HALF*CURV*RHO
IF (BDTEST < BDTOL) GOTO 650
END IF
80 CONTINUE
GOTO 680
END IF
NTRITS=NTRITS+1
!
! Severe cancellation is likely to occur if XOPT is too far from XBASE.
! If the following test holds, then XBASE is shifted so that XOPT becomes
! zero. The appropriate changes are made to BMAT and to the second
! derivatives of the current model, beginning with the changes to BMAT
! that do not depend on ZMAT. VLAG is used temporarily for working space.
! 90 IF (DSQ <= 1.0D-3*XOPTSQ) THEN
90 IF (DSQ <= 1.0D-1*XOPTSQ) THEN
model_update = .true.
FRACSQ=0.25_DP*XOPTSQ
do 95 m1=1,mv
95 v_sumpq(m1)=ZERO
DO 110 K=1,NPT
do 98 m1=1,mv
98 v_sumpq(m1)= v_sumpq(m1)+PQV(m1,K)
SUM=-HALF*XOPTSQ
DO 100 I=1,N
100 SUM=SUM+XPT(K,I)*XOPT(I)
W(NPT+K)=SUM
TEMP=FRACSQ-HALF*SUM
DO 110 I=1,N
W(I)=BMAT(K,I)
VLAG(I)=SUM*XPT(K,I)+TEMP*XOPT(I)
IP=NPT+I
DO 110 J=1,I
110 BMAT(IP,J)=BMAT(IP,J)+W(I)*VLAG(J)+VLAG(I)*W(J)
! Then the revisions of BMAT that depend on ZMAT are calculated.
DO 150 JJ=1,NPTM
SUMZ=ZERO
SUMW=ZERO
DO 120 K=1,NPT
SUMZ=SUMZ+ZMAT(K,JJ)
VLAG(K)=W(NPT+K)*ZMAT(K,JJ)
120 SUMW=SUMW+VLAG(K)
DO 140 J=1,N
SUM=(FRACSQ*SUMZ-HALF*SUMW)*XOPT(J)
DO 130 K=1,NPT
130 SUM=SUM+VLAG(K)*XPT(K,J)
W(J)=SUM
DO 140 K=1,NPT
140 BMAT(K,J)=BMAT(K,J)+SUM*ZMAT(K,JJ)
DO 150 I=1,N
IP=I+NPT
TEMP=W(I)
DO 150 J=1,I
150 BMAT(IP,J)=BMAT(IP,J)+TEMP*W(J)
! The following instructions complete the shift, including the changes
! to the second derivative parameters of the quadratic model.
IH=0
DO 170 J=1,N
do 152 m1=1,mv
152 WV(m1,J)=-HALF*v_sumpq(m1)*XOPT(J)
DO 160 K=1,NPT
do 154 m1=1,mv
154 WV(m1,J)=WV(m1,J)+PQV(m1,K)*XPT(K,J)
160 XPT(K,J)=XPT(K,J)-XOPT(J)
DO 170 I=1,J
IH=IH+1
do 165 m1=1,mv
HQV(m1,IH)=HQV(m1,IH)+WV(m1,I)*XOPT(J) &
+XOPT(I)*WV(m1,J)
165 continue
170 BMAT(NPT+I,J)=BMAT(NPT+J,I)
DO 180 I=1,N
XBASE(I)=XBASE(I)+XOPT(I)
XNEW(I)=XNEW(I)-XOPT(I)
SL(I)=SL(I)-XOPT(I)
SU(I)=SU(I)-XOPT(I)
180 XOPT(I)=ZERO
XOPTSQ=ZERO
END IF
IF (NTRITS .EQ. 0) GOTO 210
GOTO 230
! XBASE is also moved to XOPT by a call of RESCUE. This calculation is
! more expensive than the previous shift, because new matrices BMAT and
! ZMAT are generated from scratch, which may include the replacement of
! interpolation points whose positions seem to be causing near linear
! dependence in the interpolation conditions. Therefore RESCUE is called
! only if rounding errors have reduced by at least a factor of two the
! denominator of the formula for updating the H matrix. It provides a
! useful safeguard, but is not invoked in most applications of BOBYQA.
190 NFSAV=NF
KBASE=KOPT
CALL RESCUE_H (N,NPT,XL,XU,IPRINT,MAXFUN,XBASE,XPT,FVAL,&
XOPT,GOPT,HQ,PQ,BMAT,ZMAT,NDIM,SL,SU,NF,DELTA,KOPT,&
VLAG,W,W(N+NP),W(NDIM+NP),&
mv,HQV,GQV,PQV,FVAL_V,WV,v_err,v_sumpq,v_vquad,v_diff,v_temp)
model_update = .true.
! XOPT is updated now in case the branch below to label 720 is taken.
! Any updating of GOPT occurs after the branch below to label 20, which
! leads to a trust region iteration as does the branch to label 60.
XOPTSQ=ZERO
IF (KOPT .NE. KBASE) THEN
DO 200 I=1,N
XOPT(I)=XPT(KOPT,I)
200 XOPTSQ=XOPTSQ+XOPT(I)**2
END IF
IF (NF < 0) THEN
NF=MAXFUN
IF (IPRINT > 0) PRINT 390
GOTO 720
END IF
NRESC=NF
IF (NFSAV < NF) THEN
NFSAV=NF
GOTO 20
END IF
IF (NTRITS > 0) GOTO 60
! Pick two alternative vectors of variables, relative to XBASE, that
! are suitable as new positions of the KNEW-th interpolation point.
! Firstly, XNEW is set to the point on a line through XOPT and another
! interpolation point that minimizes the predicted value of the next
! denominator, subject to ||XNEW - XOPT|| .LEQ. ADELT and to the SL
! and SU bounds. Secondly, XALT is set to the best feasible point on
! a constrained version of the Cauchy step of the KNEW-th Lagrange
! function, the corresponding value of the square of this function
! being returned in CAUCHY. The choice between these alternatives is
! going to be made when the denominator is calculated.
210 CALL ALTMOV (N,NPT,XPT,XOPT,BMAT,ZMAT,NDIM,SL,SU,KOPT,&
KNEW,ADELT,XNEW,XALT,ALPHA,CAUCHY,W,W(NP),W(NDIM+1))
DO 220 I=1,N
220 D(I)=XNEW(I)-XOPT(I)
! Calculate VLAG and BETA for the current choice of D. The scalar
! product of D with XPT(K,.) is going to be held in W(NPT+K) for
! use when VQUAD is calculated.
230 DO 250 K=1,NPT
SUMA=ZERO
SUMB=ZERO
SUM=ZERO
DO 240 J=1,N
SUMA=SUMA+XPT(K,J)*D(J)
SUMB=SUMB+XPT(K,J)*XOPT(J)
240 SUM=SUM+BMAT(K,J)*D(J)
W(K)=SUMA*(HALF*SUMA+SUMB)
VLAG(K)=SUM
250 W(NPT+K)=SUMA
BETA=ZERO
DO 270 JJ=1,NPTM
SUM=ZERO
DO 260 K=1,NPT
260 SUM=SUM+ZMAT(K,JJ)*W(K)
BETA=BETA-SUM*SUM
DO 270 K=1,NPT
270 VLAG(K)=VLAG(K)+SUM*ZMAT(K,JJ)
DSQ=ZERO
BSUM=ZERO
DX=ZERO
DO 300 J=1,N
DSQ=DSQ+D(J)**2
SUM=ZERO
DO 280 K=1,NPT
280 SUM=SUM+W(K)*BMAT(K,J)
BSUM=BSUM+SUM*D(J)
JP=NPT+J
DO 290 I=1,N
290 SUM=SUM+BMAT(JP,I)*D(I)
VLAG(JP)=SUM
BSUM=BSUM+SUM*D(J)
300 DX=DX+D(J)*XOPT(J)
BETA=DX*DX+DSQ*(XOPTSQ+DX+DX+HALF*DSQ)+BETA-BSUM
VLAG(KOPT)=VLAG(KOPT)+ONE
! If NTRITS is zero, the denominator may be increased by replacing
! the step D of ALTMOV by a Cauchy step. Then RESCUE may be called if
! rounding errors have damaged the chosen denominator.
IF (NTRITS .EQ. 0) THEN
DENOM=VLAG(KNEW)**2+ALPHA*BETA
IF (DENOM < CAUCHY .AND. CAUCHY > ZERO) THEN
DO 310 I=1,N
XNEW(I)=XALT(I)
310 D(I)=XNEW(I)-XOPT(I)
CAUCHY=ZERO
GO TO 230
END IF
IF (DENOM <= HALF*VLAG(KNEW)**2) THEN
IF (NF > NRESC) GOTO 190
IF (IPRINT > 0) PRINT 320
320 FORMAT (/5X,'Return from BOBYQA_H because of much',&
' cancellation in a denominator.')
GOTO 720
END IF
! Alternatively, if NTRITS is positive, then set KNEW to the index of
! the next interpolation point to be deleted to make room for a trust
! region step. Again RESCUE may be called if rounding errors have damaged
! the chosen denominator, which is the reason for attempting to select
! KNEW before calculating the next value of the objective function.
ELSE
DELSQ=DELTA*DELTA
SCADEN=ZERO
BIGLSQ=ZERO
KNEW=0
DO 350 K=1,NPT
IF (K .EQ. KOPT) GOTO 350
HDIAG=ZERO
DO 330 JJ=1,NPTM
330 HDIAG=HDIAG+ZMAT(K,JJ)**2
DEN=BETA*HDIAG+VLAG(K)**2
DISTSQ=ZERO
DO 340 J=1,N
340 DISTSQ=DISTSQ+(XPT(K,J)-XOPT(J))**2
TEMP=DMAX1(ONE,(DISTSQ/DELSQ)**2)
IF (TEMP*DEN > SCADEN) THEN
SCADEN=TEMP*DEN
KNEW=K
DENOM=DEN
END IF
BIGLSQ=DMAX1(BIGLSQ,TEMP*VLAG(K)**2)
350 CONTINUE
IF (SCADEN <= HALF*BIGLSQ) THEN
IF (NF > NRESC) GOTO 190
IF (IPRINT > 0) PRINT 320
GOTO 720
END IF
END IF
! Put the variables for the next calculation of the objective function
! in XNEW, with any adjustments for the bounds.
! Calculate the value of the objective function at XBASE+XNEW, unless
! the limit on the number of calculations of F has been reached.
360 DO 380 I=1,N
X(I)=DMIN1(DMAX1(XL(I),XBASE(I)+XNEW(I)),XU(I))
IF (dabs(XNEW(I)-SL(I)) <= eps) X(I)=XL(I)
IF (dabs(XNEW(I)-SU(I)) <= eps) X(I)=XU(I)
380 CONTINUE
IF (NF >= MAXFUN) THEN
IF (IPRINT > 0) PRINT 390
390 FORMAT (/4X,'Return from BOBYQA_H because CALFUN has been',&
' called MAXFUN times.')
GOTO 720
END IF
NF=NF+1
IF (IPRINT .EQ. 3) THEN
PRINT 398, (X(I),I=1,N)
398 FORMAT (/4X,'Before the call to CALFUN',&
' The corresponding X is:'/(2X,5D15.6))
END IF
! dfovec(n, mv, x, v_err) provides the values of the vector function v_err(x): R^n \to R^{mv}.
! Here: n, mv, x \in R^n are input, v_err \in R^{mv} are output.
call dfovec(n, mv, x, v_err)
! f_value(mv,v_err,F) provides the value of the sum of the squres of the components of v_err(x)
! i.e. F = sum_{i=1}^{mv} v_err_i (x)^2
call f_value(mv,v_err,F)
!LS: Add penalty to result
F = F + getPenalty(x)
IF (IPRINT .EQ. 3) THEN
PRINT 400, NF,F,(X(I),I=1,N)
400 FORMAT (/4X,'Function number',I6,' F =',1PD18.10,&
' The corresponding X is:'/(2X,5D15.6))
END IF
if(F<=dmax1(1.d-12,1.d-20*FBEG)) then
print *, ""
print *, " Return: F<=dmax1(1.d-12,1.d-20*FBEG)"
go to 720
endif
IF (NTRITS .EQ. -1) THEN
FSAVE=F
GOTO 720
END IF
! Use the quadratic model to predict the change in F due to the step D,
! and set DIFF to the error of this prediction.
! VQUAD=ZERO
do 405 m1=1,mv
405 v_vquad(m1)=ZERO
IH=0
DO 410 J=1,N
do 408 m1=1,mv
408 v_vquad(m1)=v_vquad(m1)+D(J)*GQV(m1,J)
DO 410 I=1,J
IH=IH+1
TEMP=D(I)*D(J)
IF (I .EQ. J) TEMP=HALF*TEMP
do 410 m1=1,mv
410 v_vquad(m1)=v_vquad(m1)+HQV(m1,IH)*TEMP
DO 420 K=1,NPT
do 420 m1=1,mv
420 v_vquad(m1)=v_vquad(m1)+HALF*PQV(m1,K)*W(NPT+K)**2
do 425 m1=1,mv
425 v_diff(m1)=v_err(m1)-FVAL_V(m1,KOPT)-v_vquad(m1)
if (NTRITS <= 0) then
DO 426 I=1,N
426 HD1(I) = zero
IH=0
DO 427 J=1,N
DO 427 I=1,J
IH=IH+1
IF (I < J) HD1(J)=HD1(J)+HQ(IH)*D(I)
427 HD1(I)=HD1(I)+HQ(IH)*D(J)
vquad1 = zero
do 428 i=1,n
428 vquad1 = vquad1 + D(i)*(GOPT(i)+HALF*HD1(i))
endif
FOPT=FVAL(KOPT)
DIFF=F-FOPT-VQUAD1
DIFFC=DIFFB
DIFFB=DIFFA
DIFFA=DABS(DIFF)
IF (DNORM > RHO) NFSAV=NF
! Pick the next value of DELTA after a trust region step.
IF (NTRITS > 0) THEN
IF (VQUAD1 >= ZERO) THEN
IF (IPRINT > 0) PRINT 430
430 FORMAT (/4X,'Return from BOBYQA_H because a trust',&
' region step has failed to reduce Q.')
GOTO 720
END IF
RATIO=(F-FOPT)/VQUAD1
IF (RATIO <= TENTH) THEN
DELTA=DMIN1(HALF*DELTA,DNORM)
ELSE IF (RATIO <= 0.7_DP ) THEN
DELTA=DMAX1(HALF*DELTA,DNORM)
ELSE
! DELTA=DMAX1(HALF*DELTA,DNORM+DNORM)
DELTA=DMIN1(DMAX1(2._DP*DELTA,4._DP*DNORM),1.d10)
END IF
IF (DELTA <= 1.5_DP*RHO) DELTA=RHO
! Recalculate KNEW and DENOM if the new F is less than FOPT.
IF (F < FOPT) THEN
KSAV=KNEW
DENSAV=DENOM
DELSQ=DELTA*DELTA
SCADEN=ZERO
BIGLSQ=ZERO
KNEW=0
DO 460 K=1,NPT
HDIAG=ZERO
DO 440 JJ=1,NPTM
440 HDIAG=HDIAG+ZMAT(K,JJ)**2
DEN=BETA*HDIAG+VLAG(K)**2
DISTSQ=ZERO
DO 450 J=1,N
450 DISTSQ=DISTSQ+(XPT(K,J)-XNEW(J))**2
TEMP=DMAX1(ONE,(DISTSQ/DELSQ)**2)
IF (TEMP*DEN > SCADEN) THEN
SCADEN=TEMP*DEN
KNEW=K
DENOM=DEN
END IF
460 BIGLSQ=DMAX1(BIGLSQ,TEMP*VLAG(K)**2)
IF (SCADEN <= HALF*BIGLSQ) THEN
KNEW=KSAV
DENOM=DENSAV
END IF
END IF
END IF
! Update BMAT and ZMAT, so that the KNEW-th interpolation point can be
! moved. Also update the second derivative terms of the model.
CALL UPDATE (N,NPT,BMAT,ZMAT,NDIM,VLAG,BETA,DENOM,KNEW,W)
model_update = .true.
IH=0
do 465 m1=1,mv
PQVold(m1) = PQV(m1,KNEW)
PQV(m1,KNEW)=ZERO
465 continue
DO 470 I=1,N
do 468 m1=1,mv
v_temp(m1)=PQVold(m1)*XPT(KNEW,I)
468 continue
DO 470 J=1,I
IH=IH+1
do 470 m1=1,mv
HQV(m1,IH)=HQV(m1,IH)+v_temp(m1)*XPT(KNEW,J)
470 continue
DO 480 JJ=1,NPTM
do 475 m1=1,mv
v_temp(m1)=v_diff(m1)*ZMAT(KNEW,JJ)
475 continue
DO 480 K=1,NPT
do 480 m1=1,mv
PQV(m1,K)=PQV(m1,K)+v_temp(m1)*ZMAT(K,JJ)
480 continue
! Include the new interpolation point, and make the changes to GOPT at
! the old XOPT that are caused by the updating of the quadratic model.
FVAL(KNEW)=F
do 485 m1=1,mv
FVAL_V(m1,KNEW)=v_err(m1)
485 continue
DO 490 I=1,N
XPT(KNEW,I)=XNEW(I)
490 W(I)=BMAT(KNEW,I)
DO 520 K=1,NPT
SUMA=ZERO
DO 500 JJ=1,NPTM
500 SUMA=SUMA+ZMAT(KNEW,JJ)*ZMAT(K,JJ)
SUMB=ZERO
DO 510 J=1,N
510 SUMB=SUMB+XPT(K,J)*XOPT(J)
TEMP=SUMA*SUMB
DO 520 I=1,N
520 W(I)=W(I)+TEMP*XPT(K,I)
DO 530 I=1,N
do 530 m1=1,mv
GQV(m1,I)=GQV(m1,I)+v_diff(m1)*W(I)
530 continue
! Update XOPT, GOPT and KOPT if the new calculated F is less than FOPT.
IF (F < FOPT) THEN
opt_update = .true.
model_update = .true.
KOPT=KNEW
XOPTSQ=ZERO
IH=0
DO 540 J=1,N
XOPT(J)=XNEW(J)
XOPTSQ=XOPTSQ+XOPT(J)**2
DO 540 I=1,J
IH=IH+1
IF (I < J) then