-
Notifications
You must be signed in to change notification settings - Fork 14
/
cobyla.c
1948 lines (1802 loc) · 50.9 KB
/
cobyla.c
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
/*
* cobyla.c -
*
* Implementation (in C) of Mike Powell's COBYLA algorithm for minimizing a
* function of a few variables. The method is "derivatives free" (only the
* function values are needed) and accounts for constraints on the variables.
* The algorithm is described in:
*
* M.J.D. Powell, "A direct search optimization method that models the
* objective and constraint functions by linear interpolation," in Advances
* in Optimization and Numerical Analysis Mathematics and Its Applications,
* vol. 275 (eds. Susana Gomez and Jean-Pierre Hennart), Kluwer Academic
* Publishers, pp. 51-67 (1994).
*
* The present code is based on the original FORTRAN version written by Mike
* Powell who kindly provides his code on demand and has been converted to C by
* É. Thiébaut.
*
* Copyright (c) 1992, Mike Powell (FORTRAN version).
* Copyright (c) 2015, Éric Thiébaut (C version).
*/
/* To re-use as much as the code for the reverse communication routines, we use
a trick which consists in "self-including" this file with different macros
(_COBYLA_PART1, _COBYLA_PART2, etc.) defined so as to skip or modify certain
parts of the source file. */
#ifndef _COBYLA_PART1
#define _COBYLA_PART1 1
#include <stdlib.h>
#include <string.h>
#include <errno.h>
#include <math.h>
#include <stdio.h>
#include "cobyla.h"
/* Macros to deal with single/double precision. */
#undef REAL
#ifdef SINGLE_PRECISION
# define FLT(x) x ## f /* floating point literal constant */
# define REAL float
# define ABS(x) fabsf(x)
# define SQRT(x) sqrtf(x)
# define HYPOT(x) hypotf(x)
# define LOG(x) logf(x)
# define EXP(x) expf(x)
# define SIN(x) sinf(x)
# define COS(x) cosf(x)
# define TAN(x) tanf(x)
# define ASIN(x) asinf(x)
# define ACOS(x) acosf(x)
# define ATAN(x) atanf(x)
#else
# define FLT(x) x /* floating point literal constant */
# define REAL double
# define ABS(x) fabs(x)
# define SQRT(x) sqrt(x)
# define HYPOT(x) hypot(x)
# define LOG(x) log(x)
# define EXP(x) exp(x)
# define SIN(x) sin(x)
# define COS(x) cos(x)
# define TAN(x) tan(x)
# define ASIN(x) asin(x)
# define ACOS(x) acos(x)
# define ATAN(x) atan(x)
#endif
/* Helper macro for simple FORTRAN-like loops. */
#define LOOP(var,num) for (var = 1; var <= num; ++var)
#define MAX(a, b) ((a) >= (b) ? (a) : (b))
#define MIN(a, b) ((a) <= (b) ? (a) : (b))
/* Macro `HOW_MANY(a,b)` yields the minimal number of chunks with `b` elements
needed to store `a` elements. Both `a` and `b` must be integers. */
#define HOW_MANY(a,b) ((((b) - 1) + (a))/(b))
/* Macro `ROUND_UP(a,b)` yields the integer `a` rounded up to a multiple of
integer `b`. */
#define ROUND_UP(a,b) (HOW_MANY(a,b)*(b))
/* Macro `ADDRESS(type,base,offset)` yields a `type*` address at `offset` (in
bytes) from `base` address. */
#define ADDRESS(type, base, offset) ((type*)((char*)(base) + (offset)))
/* Macro `OFFSET(type, field)` yields the offset (in bytes) of member `field`
in structure/class `type`. */
#ifdef offsetof
# define OFFSET(type, field) offsetof(type, field)
#else
# define OFFSET(type, field) ((char*)&((type*)0)->field - (char*)0)
#endif
/*---------------------------------------------------------------------------*/
/* DECLARATIONS OF PRIVATE FUNCTIONS */
static int
cobylb(const INTEGER n, const INTEGER m,
cobyla_calcfc* calcfc, void* calcfc_data, REAL x[],
const REAL rhobeg, const REAL rhoend, const INTEGER iprint,
INTEGER* maxfun, REAL con[], REAL sim[], REAL simi[],
REAL datmat[], REAL a[], REAL vsig[], REAL veta[],
REAL sigbar[], REAL dx[], REAL w[], INTEGER iact[]);
static void
trstlp(const INTEGER n, const INTEGER m, const REAL a[], const REAL b[],
const REAL rho, REAL dx[], INTEGER* ifull, INTEGER iact[], REAL z[],
REAL zdota[], REAL vmultc[], REAL sdirn[], REAL dxnew[], REAL vmultd[]);
static void
print_calcfc(FILE* output, INTEGER n, INTEGER nfvals,
REAL f, REAL maxcv, const REAL x[]);
int
cobyla(INTEGER n, INTEGER m, cobyla_calcfc* calcfc, void* calcfc_data,
REAL x[], REAL rhobeg, REAL rhoend, INTEGER iprint,
INTEGER* maxfun, REAL w[], INTEGER iact[])
{
/* Partition the working space array W to provide the storage that is needed
for the main calculation. */
INTEGER mpp = m + 2;
REAL* con = w;
REAL* sim = con + mpp;
REAL* simi = sim + n*n + n;
REAL* datmat = simi + n*n;
REAL* a = datmat + n*mpp + mpp;
REAL* vsig = a + m*n + n;
REAL* veta = vsig + n;
REAL* sigbar = veta + n;
REAL* dx = sigbar + n;
REAL* work = dx + n;
return cobylb(n, m, calcfc, calcfc_data, x, rhobeg, rhoend, iprint, maxfun,
con, sim, simi, datmat, a, vsig, veta, sigbar, dx, work, iact);
}
/*---------------------------------------------------------------------------*/
/* FORTRAN support. */
#ifdef FORTRAN_NAME
/* CALCFC_WRAPPER is a wrapper function, it assumes a CALCFC subroutine written
in FORTRAN. */
static REAL
calcfc_wrapper(INTEGER n, INTEGER m, const REAL x[], REAL con[], void* data)
{
# define fc FORTRAN_NAME(calcfc,CALCFC)
REAL f;
extern void fc(INTEGER* n, INTEGER* m, const REAL x[], REAL* f, REAL con[]);
fc(&n, &m, x, &f, con);
return f;
# undef fc
}
int
FORTRAN_NAME(cobyla,COBYLA)(INTEGER* n, INTEGER* m, REAL x[], REAL* rhobeg,
REAL* rhoend, INTEGER* iprint, INTEGER* maxfun,
REAL w[], INTEGER iact[])
{
(void)cobyla(*n, *m, calcfc_wrapper, NULL, x, *rhobeg, *rhoend,
*iprint, maxfun, w, iact);
return 0;
}
#endif /* FORTRAN_NAME */
/*---------------------------------------------------------------------------*/
/* Reverse communication version. */
struct _cobyla_context {
INTEGER n; /* number of variables */
INTEGER m; /* number of constraints */
INTEGER iprint;
INTEGER maxfun;
INTEGER nfvals;
REAL rhobeg;
REAL rhoend;
/* Workspace addresses. */
INTEGER* iact;
REAL* con;
REAL* sim;
REAL* simi;
REAL* datmat;
REAL* a;
REAL* vsig;
REAL* veta;
REAL* sigbar;
REAL* dx;
REAL* w;
REAL parmu, parsig, prerec, prerem, rho, f;
INTEGER ibrnch, iflag, ifull, jdrop;
int status;
};
cobyla_context_t*
cobyla_create(INTEGER n, INTEGER m, REAL rhobeg, REAL rhoend,
INTEGER iprint, INTEGER maxfun)
{
cobyla_context_t* ctx;
long size, offset1, offset2;
INTEGER mpp;
/* Check arguments. */
if (n < 1 || m < 0 || rhobeg < rhoend || rhoend <= 0 || maxfun < 1) {
errno = EINVAL;
return NULL;
}
/* Allocate memory. */
size = sizeof(cobyla_context_t);
offset1 = ROUND_UP(size, sizeof(INTEGER));
size = offset1 + (m + 1)*sizeof(INTEGER);
offset2 = ROUND_UP(size, sizeof(REAL));
size = offset2 + (n*(3*n + 2*m + 11) + 4*m + 6)*sizeof(REAL);
ctx = (cobyla_context_t*)malloc(size);
if (ctx == NULL) {
return NULL;
}
memset(ctx, 0, size);
ctx->n = n;
ctx->m = m;
ctx->nfvals = 0; /* indicate a fresh start */
ctx->status = COBYLA_ITERATE;
ctx->iprint = iprint;
ctx->maxfun = maxfun;
ctx->rhobeg = rhobeg;
ctx->rhoend = rhoend;
/* Partition the working space array W to provide the storage that is needed
for the main calculation. */
mpp = m + 2;
ctx->iact = ADDRESS(INTEGER, ctx, offset1);
ctx->con = ADDRESS(REAL, ctx, offset2);
ctx->sim = ctx->con + mpp;
ctx->simi = ctx->sim + n*n + n;
ctx->datmat = ctx->simi + n*n;
ctx->a = ctx->datmat + n*mpp + mpp;
ctx->vsig = ctx->a + m*n + n;
ctx->veta = ctx->vsig + n;
ctx->sigbar = ctx->veta + n;
ctx->dx = ctx->sigbar + n;
ctx->w = ctx->dx + n;
return ctx;
}
void
cobyla_delete(cobyla_context_t* ctx)
{
if (ctx != NULL) {
free((void*)ctx);
}
}
int
cobyla_restart(cobyla_context_t* ctx)
{
if (ctx == NULL) {
errno = EFAULT;
return COBYLA_BAD_ADDRESS;
}
ctx->nfvals = 0;
ctx->status = COBYLA_ITERATE;
return ctx->status;
}
int
cobyla_get_status(const cobyla_context_t* ctx)
{
if (ctx == NULL) {
errno = EFAULT;
return COBYLA_BAD_ADDRESS;
}
return ctx->status;
}
INTEGER
cobyla_get_nevals(const cobyla_context_t* ctx)
{
if (ctx == NULL) {
errno = EFAULT;
return -1;
}
return ctx->nfvals;
}
REAL
cobyla_get_rho(const cobyla_context_t* ctx)
{
if (ctx == NULL) {
errno = EFAULT;
return -1;
}
return ctx->rho;
}
REAL
cobyla_get_last_f(const cobyla_context_t* ctx)
{
if (ctx == NULL) {
errno = EFAULT;
return -1;
}
return ctx->f;
}
/* Include this file with the macro _COBYLA_REVCOM defined to
generate the code `newuoa_iterate` in place of `newuob`. */
#define _COBYLA_REVCOM 1
#include __FILE__
#undef _COBYLA_REVCOM
#endif /* _COBYLA_PART1 */
/* Define macros mimicking FORTRAN indexing. */
#define A(a1,a2) a[a1 - 1 + n*(a2 - 1)]
#define CON(a1) con[a1 - 1]
#define DATMAT(a1,a2) datmat[a1 - 1 + mpp*(a2 - 1)]
#define DX(a1) dx[a1 - 1]
#define IACT(a1) iact[a1 - 1]
#define SIGBAR(a1) sigbar[a1 - 1]
#define SIM(a1,a2) sim[a1 - 1 + n*(a2 - 1)]
#define SIMI(a1,a2) simi[a1 - 1 + n*(a2 - 1)]
#define VETA(a1) veta[a1 - 1]
#define VSIG(a1) vsig[a1 - 1]
#define W(a1) w[a1 - 1]
#define X(a1) x[a1 - 1]
#ifdef _COBYLA_REVCOM
# define RESTORE(var) var = ctx->var
# define SAVE(var) ctx->var = var
#endif
#ifdef _COBYLA_REVCOM
int
cobyla_iterate(cobyla_context_t* ctx, REAL f, REAL x[], REAL c[])
#else
static int
cobylb(const INTEGER n, const INTEGER m,
cobyla_calcfc* calcfc, void* calcfc_data, REAL x[],
const REAL rhobeg, const REAL rhoend, const INTEGER iprint,
INTEGER* _maxfun, REAL con[], REAL sim[], REAL simi[],
REAL datmat[], REAL a[], REAL vsig[], REAL veta[],
REAL sigbar[], REAL dx[], REAL w[], INTEGER iact[])
#endif
{
/* Constants. */
const REAL zero = FLT(0.0);
const REAL one = FLT(1.0);
const REAL alpha = FLT(0.25);
const REAL beta = FLT(2.1);
const REAL gamma = FLT(0.5);
const REAL delta = FLT(1.1);
/* Local variables (for each type, real or integer, the first group of
variables are those which must be preserved between calls in a reverse
communication version of the algorithm). */
REAL parmu, parsig, prerec, prerem, rho;
REAL barmu, cvmaxm, cvmaxp, dxsign, edgmax, error, pareta, phi, phimin,
ratio, resmax, resnew, sum, temp, tempa, trured, vmnew, vmold;
INTEGER ibrnch, iflag, ifull, jdrop, nfvals, maxfun;
INTEGER i, j, k, l, mp, mpp, np, nbest;
#ifdef _COBYLA_REVCOM
INTEGER n, m, iprint;
REAL rhobeg, rhoend;
REAL *con, *sim, *simi, *datmat, *a, *vsig, *veta, *sigbar, *dx, *w;
INTEGER *iact;
#else
REAL f;
#endif
int status;
#ifdef _COBYLA_REVCOM
/* Minimal checking and restore initial set of variables. */
if (ctx == NULL) {
errno = EFAULT;
return COBYLA_BAD_ADDRESS;
}
RESTORE(n);
RESTORE(m);
RESTORE(iprint);
RESTORE(maxfun);
RESTORE(nfvals);
RESTORE(rhobeg);
RESTORE(rhoend);
RESTORE(iact);
RESTORE(con);
RESTORE(sim);
RESTORE(simi);
RESTORE(datmat);
RESTORE(a);
RESTORE(vsig);
RESTORE(veta);
RESTORE(sigbar);
RESTORE(dx);
RESTORE(w);
RESTORE(status);
if (x == NULL || (c == NULL && m > 0)) {
errno = EFAULT;
ctx->status = COBYLA_BAD_ADDRESS;
return COBYLA_BAD_ADDRESS;
}
if (status != COBYLA_ITERATE || nfvals < 0) {
errno = EINVAL;
ctx->status = COBYLA_CORRUPTED;
return ctx->status;
}
#endif
/* Set the initial values of some parameters. The last column of SIM holds
the optimal vertex of the current simplex, and the preceding N columns
hold the displacements from the optimal vertex to the other vertices.
Further, SIMI holds the inverse of the matrix that is contained in the
first N columns of SIM. */
np = n + 1;
mp = m + 1;
mpp = m + 2;
#ifdef _COBYLA_REVCOM
if (nfvals == 0) {
/* This is the first function evaluation. Proceed with initialization. */
status = COBYLA_INITIAL_ITERATE;
prerec = zero; /* avoids warnings */
prerem = zero; /* avoids warnings */
parsig = zero; /* avoids warnings */
iflag = 0; /* avoids warnings */
#else
maxfun = *_maxfun;
nfvals = 0;
#endif
rho = rhobeg;
parmu = zero;
jdrop = np;
ibrnch = 0;
temp = one/rho;
LOOP(i,n) {
SIM(i,np) = X(i);
LOOP(j,n) {
SIM(i,j) = zero;
SIMI(i,j) = zero;
}
SIM(i,i) = rho;
SIMI(i,i) = temp;
}
if (iprint >= 2) {
fprintf(stdout, "\n The initial value of RHO is%13.6E"
" and PARMU is set to zero.\n", (double)rho);
}
#ifdef _COBYLA_REVCOM
} else {
/* This is not the first function evaluation. Restore other local
variables and jump to the place where function value is expected. */
RESTORE(parmu);
RESTORE(parsig);
RESTORE(prerec);
RESTORE(prerem);
RESTORE(rho);
RESTORE(ibrnch);
RESTORE(iflag);
RESTORE(ifull);
RESTORE(jdrop);
goto new_eval;
}
#endif
/* Make the next call of the user-supplied subroutine CALCFC. These
instructions are also used for calling CALCFC during the iterations
of the algorithm. */
call_calcfc:
if (nfvals >= maxfun && nfvals > 0) {
status = COBYLA_TOO_MANY_EVALUATIONS;
fprintf(stderr, "Return from subroutine COBYLA because %s.\n",
cobyla_reason(status));
goto L_600;
}
#ifdef _COBYLA_REVCOM
if (status == COBYLA_INITIAL_ITERATE) {
/* We already know the functiuon value. */
status = COBYLA_ITERATE;
} else {
/* Request a new function evaluation from the caller and branch to save
local variables in the context structure. */
status = COBYLA_ITERATE;
goto save;
}
/* We arrive here when caller has computed a new function value. */
new_eval:
#else
f = calcfc(n, m, x, con, calcfc_data);
#endif
++nfvals;
/* Estimate the worst constraint RESMAX. */
resmax = zero;
#ifdef _COBYLA_REVCOM
for (k = 0; k < m; ++k) {
temp = c[k];
con[k] = temp;
if (resmax < -temp) {
resmax = -temp;
}
}
#else
for (k = 0; k < m; ++k) {
if (resmax < (temp = -con[k])) {
resmax = temp;
}
}
#endif
if (nfvals == iprint - 1 || iprint == 3) {
print_calcfc(stdout, n, nfvals, f, resmax, &X(1));
}
CON(mp) = f;
CON(mpp) = resmax;
if (ibrnch == 1) goto L_440;
/* Set the recently calculated function values in a column of DATMAT. This
array has a column for each vertex of the current simplex, the entries of
each column being the values of the constraint functions (if any) followed
by the objective function and the greatest constraint violation at the
vertex. */
LOOP(k,mpp) {
DATMAT(k,jdrop) = CON(k);
}
if (nfvals > np) goto L_130;
/* Exchange the new vertex of the initial simplex with the optimal vertex if
necessary. Then, if the initial simplex is not complete, pick its next
vertex and calculate the function values there. */
if (jdrop <= n) {
if (DATMAT(mp,np) <= f) {
X(jdrop) = SIM(jdrop,np);
} else {
SIM(jdrop,np) = X(jdrop);
LOOP(k,mpp) {
DATMAT(k,jdrop) = DATMAT(k,np);
DATMAT(k,np) = CON(k);
}
LOOP(k,jdrop) {
temp = zero;
SIM(jdrop,k) = -rho;
for (i = k; i <= jdrop; ++i) {
temp -= SIMI(i,k);
}
SIMI(jdrop,k) = temp;
}
}
}
if (nfvals <= n) {
jdrop = nfvals;
X(jdrop) += rho;
goto call_calcfc;
}
L_130:
ibrnch = 1;
/* Identify the optimal vertex of the current simplex. */
L_140:
phimin = DATMAT(mp,np) + parmu*DATMAT(mpp,np);
nbest = np;
LOOP(j,n) {
temp = DATMAT(mp,j) + parmu*DATMAT(mpp,j);
if (temp < phimin) {
nbest = j;
phimin = temp;
} else if (temp == phimin && parmu == zero) {
if (DATMAT(mpp,j) < DATMAT(mpp,nbest)) nbest = j;
}
}
/* Switch the best vertex into pole position if it is not there already,
and also update SIM, SIMI and DATMAT. */
if (nbest <= n) {
LOOP(i,mpp) {
temp = DATMAT(i,np);
DATMAT(i,np) = DATMAT(i,nbest);
DATMAT(i,nbest) = temp;
}
LOOP(i,n) {
temp = SIM(i,nbest);
SIM(i,nbest) = zero;
SIM(i,np) += temp;
tempa = zero;
LOOP(k,n) {
SIM(i,k) -= temp;
tempa -= SIMI(k,i);
}
SIMI(nbest,i) = tempa;
}
}
/* Make an error return if SIGI is a poor approximation to the inverse of
the leading N by N submatrix of SIG. */
error = zero;
LOOP(i,n) {
LOOP(j,n) {
temp = (i == j ? -one : zero);
LOOP(k,n) {
temp += SIMI(i,k)*SIM(k,j);
}
temp = ABS(temp);
error = MAX(error, temp);
}
}
if (error > FLT(0.1)) {
status = COBYLA_ROUNDING_ERRORS;
if (iprint >= 1) {
fprintf(stderr, "Return from subroutine COBYLA because %s.\n",
cobyla_reason(status));
}
goto L_600;
}
/* Calculate the coefficients of the linear approximations to the objective
and constraint functions, placing minus the objective function gradient
after the constraint gradients in the array A. The vector W is used for
working space. */
LOOP(k,mp) {
CON(k) = -DATMAT(k,np);
LOOP(j,n) {
W(j) = DATMAT(k,j) + CON(k);
}
LOOP(i,n) {
temp = zero;
LOOP(j,n) {
temp += W(j)*SIMI(j,i);
}
A(i,k) = (k == mp ? -temp : temp);
}
}
/* Calculate the values of SIGMA and ETA, and set IFLAG=0 if the current
simplex is not acceptable. */
iflag = 1;
parsig = alpha*rho;
pareta = beta*rho;
LOOP(j,n) {
REAL wsig = zero;
REAL weta = zero;
LOOP(i,n) {
wsig += SIMI(j,i)*SIMI(j,i);
weta += SIM(i,j)*SIM(i,j);
}
VSIG(j) = one/SQRT(wsig);
VETA(j) = SQRT(weta);
if (VSIG(j) < parsig || VETA(j) > pareta) iflag = 0;
}
/* If a new vertex is needed to improve acceptability, then decide which
vertex to drop from the simplex. */
if (ibrnch == 1 || iflag == 1) goto L_370;
jdrop = 0;
temp = pareta;
LOOP(j,n) {
if (VETA(j) > temp) {
jdrop = j;
temp = VETA(j);
}
}
if (jdrop == 0) {
LOOP(j,n) {
if (VSIG(j) < temp) {
jdrop = j;
temp = VSIG(j);
}
}
}
/* Calculate the step to the new vertex and its sign. */
temp = gamma*rho*VSIG(jdrop);
LOOP(i,n) {
DX(i) = temp*SIMI(jdrop,i);
}
cvmaxp = zero;
cvmaxm = zero;
sum = zero; /* FIXME: `sum` was uninitialized */
LOOP(k,mp) {
sum = zero;
LOOP(i,n) {
sum = sum + A(i,k)*DX(i);
}
if (k < mp) {
temp = DATMAT(k,np);
cvmaxp = MAX(cvmaxp, -sum - temp);
cvmaxm = MAX(cvmaxm, sum - temp);
}
}
if (parmu*(cvmaxp - cvmaxm) > sum + sum) {
dxsign = -one;
} else {
dxsign = one;
}
/* Update the elements of SIM and SIMI, and set the next X. */
temp = zero;
LOOP(i,n) {
DX(i) *= dxsign;
SIM(i,jdrop) = DX(i);
temp += SIMI(jdrop,i)*DX(i);
}
LOOP(i,n) {
SIMI(jdrop,i) /= temp;
}
LOOP(j,n) {
if (j != jdrop) {
temp = zero;
LOOP(i,n) {
temp += SIMI(j,i)*DX(i);
}
LOOP(i,n) {
SIMI(j,i) -= temp*SIMI(jdrop,i);
}
}
X(j) = SIM(j,np) + DX(j);
}
goto call_calcfc;
/* Calculate DX = X(*) - X(0). Branch if the length of DX is less than
0.5*RHO. */
L_370:
{
REAL* z = w;
REAL* zdota = z + n*n;
REAL* vmc = zdota + n;
REAL* sdirn = vmc + mp;
REAL* dxnew = sdirn + n;
REAL* vmd = dxnew + n;
trstlp(n, m, a, con, rho, dx, &ifull, iact, z, zdota,
vmc, sdirn, dxnew, vmd);
}
if (ifull == 0) {
temp = zero;
LOOP(i,n) {
temp += DX(i)*DX(i);
}
if (temp < FLT(0.25)*rho*rho) {
ibrnch = 1;
goto L_550;
}
}
/* Predict the change to F and the new maximum constraint violation if the
variables are altered from X(0) to X(0)+DX. */
resnew = zero;
CON(mp) = zero;
sum = zero; /* FIXME: `sum` was uninitialized */
LOOP(k,mp) {
sum = CON(k);
LOOP(i,n) {
sum -= A(i,k)*DX(i);
}
if (k < mp) resnew = MAX(resnew, sum);
}
/* Increase PARMU if necessary and branch back if this change alters the
optimal vertex. Otherwise PREREM and PREREC will be set to the predicted
reductions in the merit function and the maximum constraint violation
respectively. */
prerec = DATMAT(mpp,np) - resnew;
barmu = (prerec > zero ? sum/prerec : zero);
if (parmu < FLT(1.5)*barmu) {
parmu = FLT(2.0)*barmu;
if (iprint >= 2) {
fprintf(stdout, "\n Increase in PARMU to%13.6E\n", (double)parmu);
}
phi = DATMAT(mp,np) + parmu*DATMAT(mpp,np);
LOOP(j,n) {
temp = DATMAT(mp,j) + parmu*DATMAT(mpp,j);
if (temp < phi) goto L_140;
if (temp == phi && parmu == zero) {
if (DATMAT(mpp,j) < DATMAT(mpp,np)) goto L_140;
}
}
}
prerem = parmu*prerec - sum;
/* Calculate the constraint and objective functions at X(*). Then find the
actual reduction in the merit function. */
LOOP(i,n) {
X(i) = SIM(i,np) + DX(i);
}
ibrnch = 1;
goto call_calcfc;
L_440:
vmold = DATMAT(mp,np) + parmu*DATMAT(mpp,np);
vmnew = f + parmu*resmax;
trured = vmold - vmnew;
if (parmu == zero && f == DATMAT(mp,np)) {
prerem = prerec;
trured = DATMAT(mpp,np) - resmax;
}
/* Begin the operations that decide whether X(*) should replace one of the
vertices of the current simplex, the change being mandatory if TRURED is
positive. Firstly, JDROP is set to the index of the vertex that is to be
replaced. */
ratio = (trured <= zero ? one : zero);
jdrop = 0;
LOOP(j,n) {
temp = zero;
LOOP(i,n) {
temp += SIMI(j,i)*DX(i);
}
temp = ABS(temp);
if (temp > ratio) {
jdrop = j;
ratio = temp;
}
SIGBAR(j) = temp*VSIG(j);
}
/* Calculate the value of ell. */
edgmax = delta*rho;
l = 0;
LOOP(j,n) {
if (SIGBAR(j) >= parsig || SIGBAR(j) >= VSIG(j)) {
temp = VETA(j);
if (trured > zero) {
temp = zero;
LOOP(i,n) {
REAL tempb = DX(i) - SIM(i,j);
temp += tempb*tempb;
}
temp = SQRT(temp);
}
if (temp > edgmax) {
l = j;
edgmax = temp;
}
}
}
if (l > 0) jdrop = l;
if (jdrop == 0) goto L_550;
/* Revise the simplex by updating the elements of SIM, SIMI and DATMAT. */
temp = zero;
LOOP(i,n) {
SIM(i,jdrop) = DX(i);
temp += SIMI(jdrop,i)*DX(i);
}
LOOP(i,n) {
SIMI(jdrop,i) = SIMI(jdrop,i)/temp;
}
LOOP(j,n) {
if (j != jdrop) {
temp = zero;
LOOP(i,n) {
temp += SIMI(j,i)*DX(i);
}
LOOP(i,n) {
SIMI(j,i) -= temp*SIMI(jdrop,i);
}
}
}
LOOP(k,mpp) {
DATMAT(k,jdrop) = CON(k);
}
/* Branch back for further iterations with the current RHO. */
if (trured > zero && trured >= FLT(0.1)*prerem) goto L_140;
L_550:
if (iflag == 0) {
ibrnch = 0;
goto L_140;
}
/* Otherwise reduce RHO if it is not at its least value and reset PARMU. */
if (rho > rhoend) {
rho = FLT(0.5)*rho;
if (rho <= FLT(1.5)*rhoend) rho = rhoend;
if (parmu > zero) {
REAL cmin, cmax, denom;
denom = zero;
LOOP(k,mp) {
cmax = cmin = DATMAT(k,np);
LOOP(i,n) {
temp = DATMAT(k,i);
cmin = MIN(cmin, temp);
cmax = MAX(cmax, temp);
}
if (k <= m && cmin < FLT(0.5)*cmax) {
temp = MAX(cmax, zero) - cmin;
if (denom <= zero) {
denom = temp;
} else {
denom = MIN(denom, temp);
}
}
}
if (denom == zero) {
parmu = zero;
} else if (cmax - cmin < parmu*denom) {
parmu = (cmax - cmin)/denom;
}
}
if (iprint >= 2) {
fprintf(stdout, "\n Reduction in RHO to%13.6E and PARMU =%13.6E\n",
(double)rho, (double)parmu);
if (iprint == 2) {
print_calcfc(stdout, n, nfvals, DATMAT(mp,np), DATMAT(mpp,np),
&SIM(1,np));
}
}
goto L_140;
}
/* Return the best calculated values of the variables. */
if (iprint >= 1) {
fprintf(stdout, "\n Normal return from subroutine COBYLA\n");
}
status = COBYLA_SUCCESS;
if (ifull == 1) goto L_620;
L_600:
LOOP(i,n) {
X(i) = SIM(i,np);
}
f = DATMAT(mp,np);
resmax = DATMAT(mpp,np);
L_620:
if (iprint >= 1) {
print_calcfc(stdout, n, nfvals, f, resmax, &X(1));
}
#ifdef _COBYLA_REVCOM
/* Save local variables and return status. */
save:
SAVE(nfvals);
SAVE(parmu);
SAVE(parsig);
SAVE(prerec);
SAVE(prerem);
SAVE(rho);
SAVE(f);
SAVE(ibrnch);
SAVE(iflag);
SAVE(ifull);
SAVE(jdrop);
SAVE(status);
#else
/* Save number of function evaluations and return status. */
*_maxfun = nfvals;
#endif
return status;
}
/* Undefine macros mimicking FORTRAN indexing. */
#undef A
#undef CON
#undef DATMAT
#undef DX
#undef IACT
#undef SIGBAR
#undef SIM
#undef SIMI
#undef VETA
#undef VSIG
#undef W
#undef X
#ifdef _COBYLA_REVCOM
# undef RESTORE
# undef SAVE
#endif
#ifndef _COBYLA_PART2
#define _COBYLA_PART2 1
/*
* This subroutine calculates an N-component vector DX by applying the
* following two stages. In the first stage, DX is set to the shortest vector
* that minimizes the greatest violation of the constraints:
*
* A(1,K)*DX(1)+A(2,K)*DX(2)+...+A(N,K)*DX(N) >= B(K), K=2,3,...,M,
*
* subject to the Euclidean length of DX being at most RHO. If its length is
* strictly less than RHO, then we use the resultant freedom in DX to minimize
* the objective function:
*
* -A(1,M+1)*DX(1)-A(2,M+1)*DX(2)-...-A(N,M+1)*DX(N)
*
* subject to no increase in any greatest constraint violation. This notation
* allows the gradient of the objective function to be regarded as the gradient
* of a constraint. Therefore the two stages are distinguished by MCON = M and
* MCON > M respectively. It is possible that a degeneracy may prevent DX from
* attaining the target length RHO. Then the value IFULL=0 would be set, but
* usually IFULL=1 on return.
*
* In general NACT is the number of constraints in the active set and
* IACT(1),...,IACT(NACT) are their indices, while the remainder of IACT
* contains a permutation of the remaining constraint indices. Further, Z is
* an orthogonal matrix whose first NACT columns can be regarded as the result
* of Gram-Schmidt applied to the active constraint gradients. For
* J=1,2,...,NACT, the number ZDOTA(J) is the scalar product of the J-th column
* of Z with the gradient of the J-th active constraint. DX is the current
* vector of variables and here the residuals of the active constraints should