-
Notifications
You must be signed in to change notification settings - Fork 5
/
newuoa.h
1725 lines (1695 loc) · 57.2 KB
/
newuoa.h
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
/*
This is NEWUOA for unconstrained minimization. The codes were written
by Powell in Fortran and then translated to C with f2c. I further
modified the code to make it independent of libf2c and f2c.h. Please
refer to "The NEWUOA software for unconstrained optimization without
derivatives", which is available at www.damtp.cam.ac.uk, for more
information.
*/
/*
The original fortran codes are distributed without restrictions. The
C++ codes are distributed under MIT license.
*/
/* The MIT License
Copyright (c) 2004, by M.J.D. Powell <mjdp@cam.ac.uk>
2008, by Attractive Chaos <attractivechaos@aol.co.uk>
Permission is hereby granted, free of charge, to any person obtaining
a copy of this software and associated documentation files (the
"Software"), to deal in the Software without restriction, including
without limitation the rights to use, copy, modify, merge, publish,
distribute, sublicense, and/or sell copies of the Software, and to
permit persons to whom the Software is furnished to do so, subject to
the following conditions:
The above copyright notice and this permission notice shall be
included in all copies or substantial portions of the Software.
THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND
NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS
BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN
ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN
CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
SOFTWARE.
*/
/*
XCSoar notes:
- this has not been tested with fixed math type
- error conditions are not checked (printf disabled)
*/
#ifndef AC_NEWUOA_HH_
#define AC_NEWUOA_HH_
#include <math.h>
#include <algorithm>
#define M_PI 3.14159265358979323846
using namespace std;
template<class TYPE, class Func>
TYPE min_newuoa(int n, TYPE *x, Func &func, TYPE r_start=1e7, TYPE tol=1e-8, int max_iter=5000);
template<class TYPE, class Func>
static int biglag_(int n, int npt, TYPE *xopt, TYPE *xpt, TYPE *bmat, TYPE *zmat, int *idz,
int *ndim, int *knew, TYPE *delta, TYPE *d__, TYPE *alpha, TYPE *hcol, TYPE *gc,
TYPE *gd, TYPE *s, TYPE *w, Func &func)
{
/* N is the number of variables. NPT is the number of interpolation
* equations. XOPT is the best interpolation point so far. XPT
* contains the coordinates of the current interpolation
* points. BMAT provides the last N columns of H. ZMAT and IDZ give
* a factorization of the first NPT by NPT submatrix of H. NDIM is
* the first dimension of BMAT and has the value NPT+N. KNEW is the
* index of the interpolation point that is going to be moved. DEBLLTA
* is the current trust region bound. D will be set to the step from
* XOPT to the new point. ABLLPHA will be set to the KNEW-th diagonal
* element of the H matrix. HCOBLL, GC, GD, S and W will be used for
* working space. */
/* The step D is calculated in a way that attempts to maximize the
* modulus of BLLFUNC(XOPT+D), subject to the bound ||D|| .BLLE. DEBLLTA,
* where BLLFUNC is the KNEW-th BLLagrange function. */
int xpt_dim1, xpt_offset, bmat_dim1, bmat_offset, zmat_dim1, zmat_offset,
i__1, i__2, i__, j, k, iu, nptm, iterc, isave;
TYPE sp, ss, cf1, cf2, cf3, cf4, cf5, dhd, cth, tau, sth, sum, temp, step,
angle, scale, denom, delsq, tempa, tempb, twopi, taubeg, tauold, taumax,
d__1, dd, gg;
/* Parameter adjustments */
tempa = tempb = 0.0;
zmat_dim1 = npt;
zmat_offset = 1 + zmat_dim1;
zmat -= zmat_offset;
xpt_dim1 = npt;
xpt_offset = 1 + xpt_dim1;
xpt -= xpt_offset;
--xopt;
bmat_dim1 = *ndim;
bmat_offset = 1 + bmat_dim1;
bmat -= bmat_offset;
--d__; --hcol; --gc; --gd; --s; --w;
/* Function Body */
twopi = 2.0 * M_PI;
delsq = *delta * *delta;
nptm = npt - n - 1;
/* Set the first NPT components of HCOBLL to the leading elements of
* the KNEW-th column of H. */
iterc = 0;
i__1 = npt;
for (k = 1; k <= i__1; ++k) hcol[k] = 0;
i__1 = nptm;
for (j = 1; j <= i__1; ++j) {
temp = zmat[*knew + j * zmat_dim1];
if (j < *idz) temp = -temp;
i__2 = npt;
for (k = 1; k <= i__2; ++k)
hcol[k] += temp * zmat[k + j * zmat_dim1];
}
*alpha = hcol[*knew];
/* Set the unscaled initial direction D. Form the gradient of BLLFUNC
* atXOPT, and multiply D by the second derivative matrix of
* BLLFUNC. */
dd = 0;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__) {
d__[i__] = xpt[*knew + i__ * xpt_dim1] - xopt[i__];
gc[i__] = bmat[*knew + i__ * bmat_dim1];
gd[i__] = 0;
/* Computing 2nd power */
d__1 = d__[i__];
dd += d__1 * d__1;
}
i__2 = npt;
for (k = 1; k <= i__2; ++k) {
temp = 0;
sum = 0;
i__1 = n;
for (j = 1; j <= i__1; ++j) {
temp += xpt[k + j * xpt_dim1] * xopt[j];
sum += xpt[k + j * xpt_dim1] * d__[j];
}
temp = hcol[k] * temp;
sum = hcol[k] * sum;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
gc[i__] += temp * xpt[k + i__ * xpt_dim1];
gd[i__] += sum * xpt[k + i__ * xpt_dim1];
}
}
/* Scale D and GD, with a sign change if required. Set S to another
* vector in the initial two dimensional subspace. */
gg = sp = dhd = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
/* Computing 2nd power */
d__1 = gc[i__];
gg += d__1 * d__1;
sp += d__[i__] * gc[i__];
dhd += d__[i__] * gd[i__];
}
scale = *delta / sqrt(dd);
if (sp * dhd < 0) scale = -scale;
temp = 0;
if (sp * sp > dd * .99 * gg) temp = 1.0;
tau = scale * (fabs(sp) + 0.5 * scale * fabs(dhd));
if (gg * delsq < tau * .01 * tau) temp = 1.0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
d__[i__] = scale * d__[i__];
gd[i__] = scale * gd[i__];
s[i__] = gc[i__] + temp * gd[i__];
}
/* Begin the iteration by overwriting S with a vector that has the
* required length and direction, except that termination occurs if
* the given D and S are nearly parallel. */
for (iterc = 0; iterc != n; ++iterc) {
dd = sp = ss = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
/* Computing 2nd power */
d__1 = d__[i__];
dd += d__1 * d__1;
sp += d__[i__] * s[i__];
/* Computing 2nd power */
d__1 = s[i__];
ss += d__1 * d__1;
}
temp = dd * ss - sp * sp;
if (temp <= dd * 1e-8 * ss) return 0;
denom = sqrt(temp);
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
s[i__] = (dd * s[i__] - sp * d__[i__]) / denom;
w[i__] = 0;
}
/* Calculate the coefficients of the objective function on the
* circle, beginning with the multiplication of S by the second
* derivative matrix. */
i__1 = npt;
for (k = 1; k <= i__1; ++k) {
sum = 0;
i__2 = n;
for (j = 1; j <= i__2; ++j)
sum += xpt[k + j * xpt_dim1] * s[j];
sum = hcol[k] * sum;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__)
w[i__] += sum * xpt[k + i__ * xpt_dim1];
}
cf1 = cf2 = cf3 = cf4 = cf5 = 0;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__) {
cf1 += s[i__] * w[i__];
cf2 += d__[i__] * gc[i__];
cf3 += s[i__] * gc[i__];
cf4 += d__[i__] * gd[i__];
cf5 += s[i__] * gd[i__];
}
cf1 = 0.5 * cf1;
cf4 = 0.5 * cf4 - cf1;
/* Seek the value of the angle that maximizes the modulus of TAU. */
taubeg = cf1 + cf2 + cf4;
taumax = tauold = taubeg;
isave = 0;
iu = 49;
temp = twopi / (TYPE) (iu + 1);
i__2 = iu;
for (i__ = 1; i__ <= i__2; ++i__) {
angle = (TYPE) i__ *temp;
cth = cos(angle);
sth = sin(angle);
tau = cf1 + (cf2 + cf4 * cth) * cth + (cf3 + cf5 * cth) * sth;
if (fabs(tau) > fabs(taumax)) {
taumax = tau;
isave = i__;
tempa = tauold;
} else if (i__ == isave + 1) tempb = tau;
tauold = tau;
}
if (isave == 0) tempa = tau;
if (isave == iu) tempb = taubeg;
step = 0;
if (tempa != tempb) {
tempa -= taumax;
tempb -= taumax;
step = 0.5 * (tempa - tempb) / (tempa + tempb);
}
angle = temp * ((TYPE) isave + step);
/* Calculate the new D and GD. Then test for convergence. */
cth = cos(angle);
sth = sin(angle);
tau = cf1 + (cf2 + cf4 * cth) * cth + (cf3 + cf5 * cth) * sth;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__) {
d__[i__] = cth * d__[i__] + sth * s[i__];
gd[i__] = cth * gd[i__] + sth * w[i__];
s[i__] = gc[i__] + gd[i__];
}
if (fabs(tau) <= fabs(taubeg) * 1.1) return 0;
}
return 0;
}
template<class TYPE>
static int bigden_(int n, int npt, TYPE *xopt, TYPE *xpt, TYPE *bmat, TYPE *zmat, int *idz,
int *ndim, int *kopt, int *knew, TYPE *d__, TYPE *w, TYPE *vlag, TYPE *beta,
TYPE *s, TYPE *wvec, TYPE *prod)
{
/* N is the number of variables.
* NPT is the number of interpolation equations.
* XOPT is the best interpolation point so far.
* XPT contains the coordinates of the current interpolation points.
* BMAT provides the last N columns of H.
* ZMAT and IDZ give a factorization of the first NPT by NPT
submatrix of H.
* NDIM is the first dimension of BMAT and has the value NPT+N.
* KOPT is the index of the optimal interpolation point.
* KNEW is the index of the interpolation point that is going to be
moved.
* D will be set to the step from XOPT to the new point, and on
entry it should be the D that was calculated by the last call of
BIGBDLAG. The length of the initial D provides a trust region bound
on the final D.
* W will be set to Wcheck for the final choice of D.
* VBDLAG will be set to Theta*Wcheck+e_b for the final choice of D.
* BETA will be set to the value that will occur in the updating
formula when the KNEW-th interpolation point is moved to its new
position.
* S, WVEC, PROD and the private arrays DEN, DENEX and PAR will be
used for working space.
* D is calculated in a way that should provide a denominator with a
large modulus in the updating formula when the KNEW-th
interpolation point is shifted to the new position XOPT+D. */
int xpt_dim1, xpt_offset, bmat_dim1, bmat_offset, zmat_dim1, zmat_offset,
wvec_dim1, wvec_offset, prod_dim1, prod_offset, i__1, i__2, i__, j, k,
isave, iterc, jc, ip, iu, nw, ksav, nptm;
TYPE dd, d__1, ds, ss, den[9], par[9], tau, sum, diff, temp, step,
alpha, angle, denex[9], tempa, tempb, tempc, ssden, dtest, xoptd,
twopi, xopts, denold, denmax, densav, dstemp, sumold, sstemp, xoptsq;
/* Parameter adjustments */
zmat_dim1 = npt;
zmat_offset = 1 + zmat_dim1;
zmat -= zmat_offset;
xpt_dim1 = npt;
xpt_offset = 1 + xpt_dim1;
xpt -= xpt_offset;
--xopt;
prod_dim1 = *ndim;
prod_offset = 1 + prod_dim1;
prod -= prod_offset;
wvec_dim1 = *ndim;
wvec_offset = 1 + wvec_dim1;
wvec -= wvec_offset;
bmat_dim1 = *ndim;
bmat_offset = 1 + bmat_dim1;
bmat -= bmat_offset;
--d__; --w; --vlag; --s;
/* Function Body */
twopi = atan(1.0) * 8.;
nptm = npt - n - 1;
/* Store the first NPT elements of the KNEW-th column of H in W(N+1)
* to W(N+NPT). */
i__1 = npt;
for (k = 1; k <= i__1; ++k) w[n + k] = 0;
i__1 = nptm;
for (j = 1; j <= i__1; ++j) {
temp = zmat[*knew + j * zmat_dim1];
if (j < *idz) temp = -temp;
i__2 = npt;
for (k = 1; k <= i__2; ++k)
w[n + k] += temp * zmat[k + j * zmat_dim1];
}
alpha = w[n + *knew];
/* The initial search direction D is taken from the last call of
* BIGBDLAG, and the initial S is set below, usually to the direction
* from X_OPT to X_KNEW, but a different direction to an
* interpolation point may be chosen, in order to prevent S from
* being nearly parallel to D. */
dd = ds = ss = xoptsq = 0;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__) {
/* Computing 2nd power */
d__1 = d__[i__];
dd += d__1 * d__1;
s[i__] = xpt[*knew + i__ * xpt_dim1] - xopt[i__];
ds += d__[i__] * s[i__];
/* Computing 2nd power */
d__1 = s[i__];
ss += d__1 * d__1;
/* Computing 2nd power */
d__1 = xopt[i__];
xoptsq += d__1 * d__1;
}
if (ds * ds > dd * .99 * ss) {
ksav = *knew;
dtest = ds * ds / ss;
i__2 = npt;
for (k = 1; k <= i__2; ++k) {
if (k != *kopt) {
dstemp = 0;
sstemp = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
diff = xpt[k + i__ * xpt_dim1] - xopt[i__];
dstemp += d__[i__] * diff;
sstemp += diff * diff;
}
if (dstemp * dstemp / sstemp < dtest) {
ksav = k;
dtest = dstemp * dstemp / sstemp;
ds = dstemp;
ss = sstemp;
}
}
}
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__)
s[i__] = xpt[ksav + i__ * xpt_dim1] - xopt[i__];
}
ssden = dd * ss - ds * ds;
iterc = 0;
densav = 0;
/* Begin the iteration by overwriting S with a vector that has the
* required length and direction. */
BDL70:
++iterc;
temp = 1.0 / sqrt(ssden);
xoptd = xopts = 0;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__) {
s[i__] = temp * (dd * s[i__] - ds * d__[i__]);
xoptd += xopt[i__] * d__[i__];
xopts += xopt[i__] * s[i__];
}
/* Set the coefficients of the first 2.0 terms of BETA. */
tempa = 0.5 * xoptd * xoptd;
tempb = 0.5 * xopts * xopts;
den[0] = dd * (xoptsq + 0.5 * dd) + tempa + tempb;
den[1] = 2.0 * xoptd * dd;
den[2] = 2.0 * xopts * dd;
den[3] = tempa - tempb;
den[4] = xoptd * xopts;
for (i__ = 6; i__ <= 9; ++i__) den[i__ - 1] = 0;
/* Put the coefficients of Wcheck in WVEC. */
i__2 = npt;
for (k = 1; k <= i__2; ++k) {
tempa = tempb = tempc = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
tempa += xpt[k + i__ * xpt_dim1] * d__[i__];
tempb += xpt[k + i__ * xpt_dim1] * s[i__];
tempc += xpt[k + i__ * xpt_dim1] * xopt[i__];
}
wvec[k + wvec_dim1] = 0.25 * (tempa * tempa + tempb * tempb);
wvec[k + (wvec_dim1 << 1)] = tempa * tempc;
wvec[k + wvec_dim1 * 3] = tempb * tempc;
wvec[k + (wvec_dim1 << 2)] = 0.25 * (tempa * tempa - tempb * tempb);
wvec[k + wvec_dim1 * 5] = 0.5 * tempa * tempb;
}
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__) {
ip = i__ + npt;
wvec[ip + wvec_dim1] = 0;
wvec[ip + (wvec_dim1 << 1)] = d__[i__];
wvec[ip + wvec_dim1 * 3] = s[i__];
wvec[ip + (wvec_dim1 << 2)] = 0;
wvec[ip + wvec_dim1 * 5] = 0;
}
/* Put the coefficents of THETA*Wcheck in PROD. */
for (jc = 1; jc <= 5; ++jc) {
nw = npt;
if (jc == 2 || jc == 3) nw = *ndim;
i__2 = npt;
for (k = 1; k <= i__2; ++k) prod[k + jc * prod_dim1] = 0;
i__2 = nptm;
for (j = 1; j <= i__2; ++j) {
sum = 0;
i__1 = npt;
for (k = 1; k <= i__1; ++k) sum += zmat[k + j * zmat_dim1] * wvec[k + jc * wvec_dim1];
if (j < *idz) sum = -sum;
i__1 = npt;
for (k = 1; k <= i__1; ++k)
prod[k + jc * prod_dim1] += sum * zmat[k + j * zmat_dim1];
}
if (nw == *ndim) {
i__1 = npt;
for (k = 1; k <= i__1; ++k) {
sum = 0;
i__2 = n;
for (j = 1; j <= i__2; ++j)
sum += bmat[k + j * bmat_dim1] * wvec[npt + j + jc * wvec_dim1];
prod[k + jc * prod_dim1] += sum;
}
}
i__1 = n;
for (j = 1; j <= i__1; ++j) {
sum = 0;
i__2 = nw;
for (i__ = 1; i__ <= i__2; ++i__)
sum += bmat[i__ + j * bmat_dim1] * wvec[i__ + jc * wvec_dim1];
prod[npt + j + jc * prod_dim1] = sum;
}
}
/* Include in DEN the part of BETA that depends on THETA. */
i__1 = *ndim;
for (k = 1; k <= i__1; ++k) {
sum = 0;
for (i__ = 1; i__ <= 5; ++i__) {
par[i__ - 1] = 0.5 * prod[k + i__ * prod_dim1] * wvec[k + i__ * wvec_dim1];
sum += par[i__ - 1];
}
den[0] = den[0] - par[0] - sum;
tempa = prod[k + prod_dim1] * wvec[k + (wvec_dim1 << 1)] + prod[k + (
prod_dim1 << 1)] * wvec[k + wvec_dim1];
tempb = prod[k + (prod_dim1 << 1)] * wvec[k + (wvec_dim1 << 2)] +
prod[k + (prod_dim1 << 2)] * wvec[k + (wvec_dim1 << 1)];
tempc = prod[k + prod_dim1 * 3] * wvec[k + wvec_dim1 * 5] + prod[k +
prod_dim1 * 5] * wvec[k + wvec_dim1 * 3];
den[1] = den[1] - tempa - 0.5 * (tempb + tempc);
den[5] -= 0.5 * (tempb - tempc);
tempa = prod[k + prod_dim1] * wvec[k + wvec_dim1 * 3] + prod[k +
prod_dim1 * 3] * wvec[k + wvec_dim1];
tempb = prod[k + (prod_dim1 << 1)] * wvec[k + wvec_dim1 * 5] + prod[k
+ prod_dim1 * 5] * wvec[k + (wvec_dim1 << 1)];
tempc = prod[k + prod_dim1 * 3] * wvec[k + (wvec_dim1 << 2)] + prod[k
+ (prod_dim1 << 2)] * wvec[k + wvec_dim1 * 3];
den[2] = den[2] - tempa - 0.5 * (tempb - tempc);
den[6] -= 0.5 * (tempb + tempc);
tempa = prod[k + prod_dim1] * wvec[k + (wvec_dim1 << 2)] + prod[k + (
prod_dim1 << 2)] * wvec[k + wvec_dim1];
den[3] = den[3] - tempa - par[1] + par[2];
tempa = prod[k + prod_dim1] * wvec[k + wvec_dim1 * 5] + prod[k +
prod_dim1 * 5] * wvec[k + wvec_dim1];
tempb = prod[k + (prod_dim1 << 1)] * wvec[k + wvec_dim1 * 3] + prod[k
+ prod_dim1 * 3] * wvec[k + (wvec_dim1 << 1)];
den[4] = den[4] - tempa - 0.5 * tempb;
den[7] = den[7] - par[3] + par[4];
tempa = prod[k + (prod_dim1 << 2)] * wvec[k + wvec_dim1 * 5] + prod[k
+ prod_dim1 * 5] * wvec[k + (wvec_dim1 << 2)];
den[8] -= 0.5 * tempa;
}
/* Extend DEN so that it holds all the coefficients of DENOM. */
sum = 0;
for (i__ = 1; i__ <= 5; ++i__) {
/* Computing 2nd power */
d__1 = prod[*knew + i__ * prod_dim1];
par[i__ - 1] = 0.5 * (d__1 * d__1);
sum += par[i__ - 1];
}
denex[0] = alpha * den[0] + par[0] + sum;
tempa = 2.0 * prod[*knew + prod_dim1] * prod[*knew + (prod_dim1 << 1)];
tempb = prod[*knew + (prod_dim1 << 1)] * prod[*knew + (prod_dim1 << 2)];
tempc = prod[*knew + prod_dim1 * 3] * prod[*knew + prod_dim1 * 5];
denex[1] = alpha * den[1] + tempa + tempb + tempc;
denex[5] = alpha * den[5] + tempb - tempc;
tempa = 2.0 * prod[*knew + prod_dim1] * prod[*knew + prod_dim1 * 3];
tempb = prod[*knew + (prod_dim1 << 1)] * prod[*knew + prod_dim1 * 5];
tempc = prod[*knew + prod_dim1 * 3] * prod[*knew + (prod_dim1 << 2)];
denex[2] = alpha * den[2] + tempa + tempb - tempc;
denex[6] = alpha * den[6] + tempb + tempc;
tempa = 2.0 * prod[*knew + prod_dim1] * prod[*knew + (prod_dim1 << 2)];
denex[3] = alpha * den[3] + tempa + par[1] - par[2];
tempa = 2.0 * prod[*knew + prod_dim1] * prod[*knew + prod_dim1 * 5];
denex[4] = alpha * den[4] + tempa + prod[*knew + (prod_dim1 << 1)] * prod[
*knew + prod_dim1 * 3];
denex[7] = alpha * den[7] + par[3] - par[4];
denex[8] = alpha * den[8] + prod[*knew + (prod_dim1 << 2)] * prod[*knew +
prod_dim1 * 5];
/* Seek the value of the angle that maximizes the modulus of DENOM. */
sum = denex[0] + denex[1] + denex[3] + denex[5] + denex[7];
denold = denmax = sum;
isave = 0;
iu = 49;
temp = twopi / (TYPE) (iu + 1);
par[0] = 1.0;
i__1 = iu;
for (i__ = 1; i__ <= i__1; ++i__) {
angle = (TYPE) i__ *temp;
par[1] = cos(angle);
par[2] = sin(angle);
for (j = 4; j <= 8; j += 2) {
par[j - 1] = par[1] * par[j - 3] - par[2] * par[j - 2];
par[j] = par[1] * par[j - 2] + par[2] * par[j - 3];
}
sumold = sum;
sum = 0;
for (j = 1; j <= 9; ++j)
sum += denex[j - 1] * par[j - 1];
if (fabs(sum) > fabs(denmax)) {
denmax = sum;
isave = i__;
tempa = sumold;
} else if (i__ == isave + 1) {
tempb = sum;
}
}
if (isave == 0) tempa = sum;
if (isave == iu) tempb = denold;
step = 0;
if (tempa != tempb) {
tempa -= denmax;
tempb -= denmax;
step = 0.5 * (tempa - tempb) / (tempa + tempb);
}
angle = temp * ((TYPE) isave + step);
/* Calculate the new parameters of the denominator, the new VBDLAG
* vector and the new D. Then test for convergence. */
par[1] = cos(angle);
par[2] = sin(angle);
for (j = 4; j <= 8; j += 2) {
par[j - 1] = par[1] * par[j - 3] - par[2] * par[j - 2];
par[j] = par[1] * par[j - 2] + par[2] * par[j - 3];
}
*beta = 0;
denmax = 0;
for (j = 1; j <= 9; ++j) {
*beta += den[j - 1] * par[j - 1];
denmax += denex[j - 1] * par[j - 1];
}
i__1 = *ndim;
for (k = 1; k <= i__1; ++k) {
vlag[k] = 0;
for (j = 1; j <= 5; ++j)
vlag[k] += prod[k + j * prod_dim1] * par[j - 1];
}
tau = vlag[*knew];
dd = tempa = tempb = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
d__[i__] = par[1] * d__[i__] + par[2] * s[i__];
w[i__] = xopt[i__] + d__[i__];
/* Computing 2nd power */
d__1 = d__[i__];
dd += d__1 * d__1;
tempa += d__[i__] * w[i__];
tempb += w[i__] * w[i__];
}
if (iterc >= n) goto BDL340;
if (iterc > 1) densav = std::max(densav, denold);
if (fabs(denmax) <= fabs(densav) * 1.1) goto BDL340;
densav = denmax;
/* Set S to 0.5 the gradient of the denominator with respect to
* D. Then branch for the next iteration. */
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
temp = tempa * xopt[i__] + tempb * d__[i__] - vlag[npt + i__];
s[i__] = tau * bmat[*knew + i__ * bmat_dim1] + alpha * temp;
}
i__1 = npt;
for (k = 1; k <= i__1; ++k) {
sum = 0;
i__2 = n;
for (j = 1; j <= i__2; ++j)
sum += xpt[k + j * xpt_dim1] * w[j];
temp = (tau * w[n + k] - alpha * vlag[k]) * sum;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__)
s[i__] += temp * xpt[k + i__ * xpt_dim1];
}
ss = 0;
ds = 0;
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__) {
/* Computing 2nd power */
d__1 = s[i__];
ss += d__1 * d__1;
ds += d__[i__] * s[i__];
}
ssden = dd * ss - ds * ds;
if (ssden >= dd * 1e-8 * ss) goto BDL70;
/* Set the vector W before the RETURN from the subroutine. */
BDL340:
i__2 = *ndim;
for (k = 1; k <= i__2; ++k) {
w[k] = 0;
for (j = 1; j <= 5; ++j) w[k] += wvec[k + j * wvec_dim1] * par[j - 1];
}
vlag[*kopt] += 1.0;
return 0;
}
template<class TYPE>
int trsapp_(int n, int npt, TYPE * xopt, TYPE * xpt, TYPE * gq, TYPE * hq, TYPE * pq,
TYPE * delta, TYPE * step, TYPE * d__, TYPE * g, TYPE * hd, TYPE * hs, TYPE * crvmin)
{
/* The arguments NPT, XOPT, XPT, GQ, HQ and PQ have their usual
* meanings, in order to define the current quadratic model Q.
* DETRLTA is the trust region radius, and has to be positive. STEP
* will be set to the calculated trial step. The arrays D, G, HD and
* HS will be used for working space. CRVMIN will be set to the
* least curvature of H aint the conjugate directions that occur,
* except that it is set to 0 if STEP goes all the way to the trust
* region boundary. The calculation of STEP begins with the
* truncated conjugate gradient method. If the boundary of the trust
* region is reached, then further changes to STEP may be made, each
* one being in the 2D space spanned by the current STEP and the
* corresponding gradient of Q. Thus STEP should provide a
* substantial reduction to Q within the trust region. */
int xpt_dim1, xpt_offset, i__1, i__2, i__, j, k, ih, iu, iterc,
isave, itersw, itermax;
TYPE d__1, d__2, dd, cf, dg, gg, ds, sg, ss, dhd, dhs,
cth, sgk, shs, sth, qadd, qbeg, qred, qmin, temp,
qsav, qnew, ggbeg, alpha, angle, reduc, ggsav, delsq,
tempa, tempb, bstep, ratio, twopi, angtest;
/* Parameter adjustments */
tempa = tempb = shs = sg = bstep = ggbeg = gg = qred = dd = 0.0;
xpt_dim1 = npt;
xpt_offset = 1 + xpt_dim1;
xpt -= xpt_offset;
--xopt; --gq; --hq; --pq; --step; --d__; --g; --hd; --hs;
/* Function Body */
twopi = 2.0 * M_PI;
delsq = *delta * *delta;
iterc = 0;
itermax = n;
itersw = itermax;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) d__[i__] = xopt[i__];
goto TRL170;
/* Prepare for the first line search. */
TRL20:
qred = dd = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
step[i__] = 0;
hs[i__] = 0;
g[i__] = gq[i__] + hd[i__];
d__[i__] = -g[i__];
/* Computing 2nd power */
d__1 = d__[i__];
dd += d__1 * d__1;
}
*crvmin = 0;
if (dd == 0) goto TRL160;
ds = ss = 0;
gg = dd;
ggbeg = gg;
/* Calculate the step to the trust region boundary and the product
* HD. */
TRL40:
++iterc;
temp = delsq - ss;
bstep = temp / (ds + sqrt(ds * ds + dd * temp));
goto TRL170;
TRL50:
dhd = 0;
i__1 = n;
for (j = 1; j <= i__1; ++j) dhd += d__[j] * hd[j];
/* Update CRVMIN and set the step-length ATRLPHA. */
alpha = bstep;
if (dhd > 0) {
temp = dhd / dd;
if (iterc == 1) *crvmin = temp;
*crvmin = min(*crvmin, temp);
/* Computing MIN */
d__1 = alpha, d__2 = gg / dhd;
alpha = min(d__1, d__2);
}
qadd = alpha * (gg - 0.5 * alpha * dhd);
qred += qadd;
/* Update STEP and HS. */
ggsav = gg;
gg = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
step[i__] += alpha * d__[i__];
hs[i__] += alpha * hd[i__];
/* Computing 2nd power */
d__1 = g[i__] + hs[i__];
gg += d__1 * d__1;
}
/* Begin another conjugate direction iteration if required. */
if (alpha < bstep) {
if (qadd <= qred * .01 || gg <= ggbeg * 1e-4 || iterc == itermax) goto TRL160;
temp = gg / ggsav;
dd = ds = ss = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
d__[i__] = temp * d__[i__] - g[i__] - hs[i__];
/* Computing 2nd power */
d__1 = d__[i__];
dd += d__1 * d__1;
ds += d__[i__] * step[i__];
/* Computing 2nd power */
d__1 = step[i__];
ss += d__1 * d__1;
}
if (ds <= 0) goto TRL160;
if (ss < delsq) goto TRL40;
}
*crvmin = 0;
itersw = iterc;
/* Test whether an alternative iteration is required. */
TRL90:
if (gg <= ggbeg * 1e-4) goto TRL160;
sg = 0;
shs = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
sg += step[i__] * g[i__];
shs += step[i__] * hs[i__];
}
sgk = sg + shs;
angtest = sgk / sqrt(gg * delsq);
if (angtest <= -.99) goto TRL160;
/* Begin the alternative iteration by calculating D and HD and some
* scalar products. */
++iterc;
temp = sqrt(delsq * gg - sgk * sgk);
tempa = delsq / temp;
tempb = sgk / temp;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__)
d__[i__] = tempa * (g[i__] + hs[i__]) - tempb * step[i__];
goto TRL170;
TRL120:
dg = dhd = dhs = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
dg += d__[i__] * g[i__];
dhd += hd[i__] * d__[i__];
dhs += hd[i__] * step[i__];
}
/* Seek the value of the angle that minimizes Q. */
cf = 0.5 * (shs - dhd);
qbeg = sg + cf;
qsav = qmin = qbeg;
isave = 0;
iu = 49;
temp = twopi / (TYPE) (iu + 1);
i__1 = iu;
for (i__ = 1; i__ <= i__1; ++i__) {
angle = (TYPE) i__ *temp;
cth = cos(angle);
sth = sin(angle);
qnew = (sg + cf * cth) * cth + (dg + dhs * cth) * sth;
if (qnew < qmin) {
qmin = qnew;
isave = i__;
tempa = qsav;
} else if (i__ == isave + 1) tempb = qnew;
qsav = qnew;
}
if ((TYPE) isave == 0) tempa = qnew;
if (isave == iu) tempb = qbeg;
angle = 0;
if (tempa != tempb) {
tempa -= qmin;
tempb -= qmin;
angle = 0.5 * (tempa - tempb) / (tempa + tempb);
}
angle = temp * ((TYPE) isave + angle);
/* Calculate the new STEP and HS. Then test for convergence. */
cth = cos(angle);
sth = sin(angle);
reduc = qbeg - (sg + cf * cth) * cth - (dg + dhs * cth) * sth;
gg = 0;
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) {
step[i__] = cth * step[i__] + sth * d__[i__];
hs[i__] = cth * hs[i__] + sth * hd[i__];
/* Computing 2nd power */
d__1 = g[i__] + hs[i__];
gg += d__1 * d__1;
}
qred += reduc;
ratio = reduc / qred;
if (iterc < itermax && ratio > .01) goto TRL90;
TRL160:
return 0;
/* The following instructions act as a subroutine for setting the
* vector HD to the vector D multiplied by the second derivative
* matrix of Q. They are called from three different places, which
* are distinguished by the value of ITERC. */
TRL170:
i__1 = n;
for (i__ = 1; i__ <= i__1; ++i__) hd[i__] = 0;
i__1 = npt;
for (k = 1; k <= i__1; ++k) {
temp = 0;
i__2 = n;
for (j = 1; j <= i__2; ++j)
temp += xpt[k + j * xpt_dim1] * d__[j];
temp *= pq[k];
i__2 = n;
for (i__ = 1; i__ <= i__2; ++i__)
hd[i__] += temp * xpt[k + i__ * xpt_dim1];
}
ih = 0;
i__2 = n;
for (j = 1; j <= i__2; ++j) {
i__1 = j;
for (i__ = 1; i__ <= i__1; ++i__) {
++ih;
if (i__ < j) hd[j] += hq[ih] * d__[i__];
hd[i__] += hq[ih] * d__[j];
}
}
if (iterc == 0) goto TRL20;
if (iterc <= itersw) goto TRL50;
goto TRL120;
}
template<class TYPE>
static int update_(int n, int npt, TYPE *bmat, TYPE *zmat, int *idz, int *ndim, TYPE *vlag,
TYPE *beta, int *knew, TYPE *w)
{
/* The arrays BMAT and ZMAT with IDZ are updated, in order to shift
* the interpolation point that has index KNEW. On entry, VLAG
* contains the components of the vector Theta*Wcheck+e_b of the
* updating formula (6.11), and BETA holds the value of the
* parameter that has this name. The vector W is used for working
* space. */
int bmat_dim1, bmat_offset, zmat_dim1, zmat_offset, i__1, i__2, i__,
j, ja, jb, jl, jp, nptm, iflag;
TYPE d__1, d__2, tau, temp, scala, scalb, alpha, denom, tempa, tempb, tausq;
/* Parameter adjustments */
tempb = 0.0;
zmat_dim1 = npt;
zmat_offset = 1 + zmat_dim1;
zmat -= zmat_offset;
bmat_dim1 = *ndim;
bmat_offset = 1 + bmat_dim1;
bmat -= bmat_offset;
--vlag;
--w;
/* Function Body */
nptm = npt - n - 1;
/* Apply the rotations that put zeros in the KNEW-th row of ZMAT. */
jl = 1;
i__1 = nptm;
for (j = 2; j <= i__1; ++j) {
if (j == *idz) {
jl = *idz;
} else if (zmat[*knew + j * zmat_dim1] != 0) {
/* Computing 2nd power */
d__1 = zmat[*knew + jl * zmat_dim1];
/* Computing 2nd power */
d__2 = zmat[*knew + j * zmat_dim1];
temp = sqrt(d__1 * d__1 + d__2 * d__2);
tempa = zmat[*knew + jl * zmat_dim1] / temp;
tempb = zmat[*knew + j * zmat_dim1] / temp;
i__2 = npt;
for (i__ = 1; i__ <= i__2; ++i__) {
temp = tempa * zmat[i__ + jl * zmat_dim1] + tempb * zmat[i__
+ j * zmat_dim1];
zmat[i__ + j * zmat_dim1] = tempa * zmat[i__ + j * zmat_dim1]
- tempb * zmat[i__ + jl * zmat_dim1];
zmat[i__ + jl * zmat_dim1] = temp;
}
zmat[*knew + j * zmat_dim1] = 0;
}
}
/* Put the first NPT components of the KNEW-th column of HLAG into
* W, and calculate the parameters of the updating formula. */
tempa = zmat[*knew + zmat_dim1];
if (*idz >= 2) tempa = -tempa;
if (jl > 1) tempb = zmat[*knew + jl * zmat_dim1];
i__1 = npt;
for (i__ = 1; i__ <= i__1; ++i__) {
w[i__] = tempa * zmat[i__ + zmat_dim1];
if (jl > 1) w[i__] += tempb * zmat[i__ + jl * zmat_dim1];
}
alpha = w[*knew];
tau = vlag[*knew];
tausq = tau * tau;
denom = alpha * *beta + tausq;
vlag[*knew] -= 1.0;
/* Complete the updating of ZMAT when there is only 1.0 nonzero
* element in the KNEW-th row of the new matrix ZMAT, but, if IFLAG
* is set to 1.0, then the first column of ZMAT will be exchanged
* with another 1.0 later. */
iflag = 0;
if (jl == 1) {
temp = sqrt((fabs(denom)));
tempb = tempa / temp;
tempa = tau / temp;
i__1 = npt;
for (i__ = 1; i__ <= i__1; ++i__)
zmat[i__ + zmat_dim1] = tempa * zmat[i__ + zmat_dim1] - tempb *
vlag[i__];
if (*idz == 1 && temp < 0) *idz = 2;
if (*idz >= 2 && temp >= 0) iflag = 1;
} else {
/* Complete the updating of ZMAT in the alternative case. */
ja = 1;
if (*beta >= 0) {
ja = jl;
}
jb = jl + 1 - ja;
temp = zmat[*knew + jb * zmat_dim1] / denom;
tempa = temp * *beta;
tempb = temp * tau;
temp = zmat[*knew + ja * zmat_dim1];
scala = 1.0 / sqrt(fabs(*beta) * temp * temp + tausq);
scalb = scala * sqrt((fabs(denom)));
i__1 = npt;
for (i__ = 1; i__ <= i__1; ++i__) {
zmat[i__ + ja * zmat_dim1] = scala * (tau * zmat[i__ + ja *
zmat_dim1] - temp * vlag[i__]);
zmat[i__ + jb * zmat_dim1] = scalb * (zmat[i__ + jb * zmat_dim1]
- tempa * w[i__] - tempb * vlag[i__]);
}
if (denom <= 0) {
if (*beta < 0) ++(*idz);
if (*beta >= 0) iflag = 1;
}
}
/* IDZ is reduced in the following case, and usually the first
* column of ZMAT is exchanged with a later 1.0. */
if (iflag == 1) {
--(*idz);
i__1 = npt;
for (i__ = 1; i__ <= i__1; ++i__) {
temp = zmat[i__ + zmat_dim1];
zmat[i__ + zmat_dim1] = zmat[i__ + *idz * zmat_dim1];
zmat[i__ + *idz * zmat_dim1] = temp;
}
}
/* Finally, update the matrix BMAT. */
i__1 = n;
for (j = 1; j <= i__1; ++j) {
jp = npt + j;
w[jp] = bmat[*knew + j * bmat_dim1];
tempa = (alpha * vlag[jp] - tau * w[jp]) / denom;
tempb = (-(*beta) * w[jp] - tau * vlag[jp]) / denom;
i__2 = jp;
for (i__ = 1; i__ <= i__2; ++i__) {
bmat[i__ + j * bmat_dim1] = bmat[i__ + j * bmat_dim1] + tempa *
vlag[i__] + tempb * w[i__];
if (i__ > npt) {
bmat[jp + (i__ - npt) * bmat_dim1] = bmat[i__ + j *
bmat_dim1];
}
}
}
return 0;
}
template<class TYPE, class Func>
static TYPE newuob_(int n, int npt, TYPE *x,