chiark / gitweb /
added test function with minimum at side (i.e. one active bound constraint), which...
[nlopt.git] / test / testopt.cpp
1 #include <stdlib.h>
2 #include <stdio.h>
3 #include <math.h>
4 #include <string.h>
5
6 #include "config.h"
7
8 #ifdef HAVE_UNISTD_H
9 #  include <unistd.h>
10 #endif
11 #ifdef HAVE_GETOPT_H
12 #  include <getopt.h>
13 #endif
14
15 #include "nlopt.h"
16 #include "nlopt-util.h"
17 #include "testfuncs.h"
18
19 static nlopt_algorithm algorithm = NLOPT_GN_DIRECT_L;
20 static double ftol_rel = 0, ftol_abs = 0, xtol_rel = 0, xtol_abs = 0, minf_max_delta = -HUGE_VAL;
21 static int maxeval = 1000, iterations = 1, center_start = 0;
22 static double maxtime = 0.0;
23 static double xinit_tol = -1;
24 static int force_constraints = 0;
25
26 static void listalgs(FILE *f)
27 {
28   int i;
29   fprintf(f, "Available algorithms:\n");
30   for (i = 0; i < NLOPT_NUM_ALGORITHMS; ++i)
31     fprintf(f, "  %2d: %s\n", i, nlopt_algorithm_name((nlopt_algorithm) i));
32 }
33
34 static void listfuncs(FILE *f)
35 {
36   int i;
37   fprintf(f, "Available objective functions:\n");
38   for (i = 0; i < NTESTFUNCS; ++i)
39     fprintf(f, "  %2d: %s (%d dims)\n", i, testfuncs[i].name, testfuncs[i].n);
40 }
41
42 typedef struct {
43   const double *lb, *ub;
44   nlopt_func f;
45   void *f_data;
46 } bounds_wrap_data;
47
48 static double bounds_wrap_func(int n, const double *x, double *grad, void *d_)
49 {
50   bounds_wrap_data *d = (bounds_wrap_data *) d_;
51   int i;
52   for (i = 0; i < n; ++i)
53     if (x[i] < d->lb[i] || x[i] > d->ub[i])
54       break;
55   if (i < n)
56     fprintf(stderr, "WARNING: bounds violated by x[%d] = %g\n", i, x[i]);
57   return d->f(n, x, grad, d->f_data);
58 }
59
60 static int test_function(int ifunc)
61 {
62   testfunc func;
63   int i, iter;
64   double *x, minf, minf_max, f0, *xtabs, *lb, *ub;
65   nlopt_result ret;
66   double start = nlopt_seconds();
67   int total_count = 0, max_count = 0, min_count = 1<<30;
68   double total_err = 0, max_err = 0;
69   bounds_wrap_data bw;
70   
71   if (ifunc < 0 || ifunc >= NTESTFUNCS) {
72     fprintf(stderr, "testopt: invalid function %d\n", ifunc);
73     listfuncs(stderr);
74     return 0;
75   }
76   func = testfuncs[ifunc];
77   x = (double *) malloc(sizeof(double) * func.n * 5);
78   if (!x) { fprintf(stderr, "testopt: Out of memory!\n"); return 0; }
79
80   lb = x + func.n * 3;
81   ub = lb + func.n;
82   xtabs = x + func.n * 2;
83   bw.lb = lb;
84   bw.ub = ub;
85   bw.f = func.f;
86   bw.f_data = func.f_data;
87
88   for (i = 0; i < func.n; ++i) xtabs[i] = xtol_abs;
89   minf_max = minf_max_delta > (-HUGE_VAL) ? minf_max_delta + func.minf : (-HUGE_VAL);
90   
91   printf("-----------------------------------------------------------\n");
92   printf("Optimizing %s (%d dims) using %s algorithm\n",
93          func.name, func.n, nlopt_algorithm_name(algorithm));
94   printf("lower bounds at lb = [");
95   for (i = 0; i < func.n; ++i) printf(" %g", func.lb[i]);
96   printf("]\n");
97   printf("upper bounds at ub = [");
98   for (i = 0; i < func.n; ++i) printf(" %g", func.ub[i]);
99   printf("]\n");
100   memcpy(lb, func.lb, func.n * sizeof(double));
101   memcpy(ub, func.ub, func.n * sizeof(double));
102   if (force_constraints) {
103     for (i = 0; i < func.n; ++i) {
104       if (nlopt_iurand(2) == 0)
105         ub[i] = nlopt_urand(lb[i], func.xmin[i]);
106       else
107         lb[i] = nlopt_urand(func.xmin[i], ub[i]);
108     }
109     printf("adjusted lower bounds at lb = [");
110     for (i = 0; i < func.n; ++i) printf(" %g", lb[i]);
111     printf("]\n");
112     printf("adjusted upper bounds at ub = [");
113     for (i = 0; i < func.n; ++i) printf(" %g", ub[i]);
114     printf("]\n");
115   }
116
117   if (fabs(func.f(func.n, func.xmin, 0, func.f_data) - func.minf) > 1e-8) {
118     fprintf(stderr, "BUG: function does not achieve given lower bound!\n");
119     fprintf(stderr, "f(%g", func.xmin[0]);
120     for (i = 1; i < func.n; ++i) fprintf(stderr, ", %g", func.xmin[i]);
121     fprintf(stderr, ") = %0.16g instead of %0.16g, |diff| = %g\n", 
122             func.f(func.n, func.xmin, 0, func.f_data), func.minf,
123             fabs(func.f(func.n, func.xmin, 0, func.f_data) - func.minf));
124     return 0;
125   }
126
127   for (iter = 0; iter < iterations; ++iter) {
128     double val;
129     testfuncs_counter = 0;
130
131     printf("Starting guess x = [");
132     for (i = 0; i < func.n; ++i) {
133       if (center_start)
134         x[i] = (ub[i] + lb[i]) * 0.5;
135       else if (xinit_tol < 0) { /* random starting point near center of box */
136         double dx = (ub[i] - lb[i]) * 0.25;
137         double xm = 0.5 * (ub[i] + lb[i]);
138         x[i] = nlopt_urand(xm - dx, xm + dx);
139       }
140       else {
141         x[i] = nlopt_urand(-xinit_tol, xinit_tol)
142           + (1 + nlopt_urand(-xinit_tol, xinit_tol)) * func.xmin[i];
143         if (x[i] > ub[i]) x[i] = ub[i];
144         else if (x[i] < lb[i]) x[i] = lb[i];
145       }
146       printf(" %g", x[i]);
147     }
148     printf("]\n");
149     f0 = func.f(func.n, x, x + func.n, func.f_data);
150     printf("Starting function value = %g\n", f0);
151     
152     if (iter == 0 && testfuncs_verbose && func.has_gradient) {
153       printf("checking gradient:\n");
154       for (i = 0; i < func.n; ++i) {
155         double f;
156         x[i] *= 1 + 1e-6;
157         f = func.f(func.n, x, NULL, func.f_data);
158         x[i] /= 1 + 1e-6;
159         printf("  grad[%d] = %g vs. numerical derivative %g\n",
160                i, x[i + func.n], (f - f0) / (x[i] * 1e-6));
161       }
162     }
163     
164     ret = nlopt_minimize(algorithm,
165                          func.n, bounds_wrap_func, &bw,
166                          lb, ub,
167                          x, &minf,
168                          minf_max, ftol_rel, ftol_abs, xtol_rel, xtabs,
169                          maxeval, maxtime);
170     printf("finished after %g seconds.\n", nlopt_seconds() - start);
171     printf("return code %d from nlopt_minimize\n", ret);
172     if (ret < 0) {
173       fprintf(stderr, "testopt: error in nlopt_minimize\n");
174       return 0;
175     }
176     printf("Found minimum f = %g after %d evaluations.\n", 
177            minf, testfuncs_counter);
178     total_count += testfuncs_counter;
179     if (testfuncs_counter > max_count) max_count = testfuncs_counter;
180     if (testfuncs_counter < min_count) min_count = testfuncs_counter;
181     printf("Minimum at x = [");
182     for (i = 0; i < func.n; ++i) printf(" %g", x[i]);
183     printf("]\n");
184     printf("|f - minf| = %g, |f - minf| / |minf| = %e\n",
185            fabs(minf - func.minf), fabs(minf - func.minf) / fabs(func.minf));
186     total_err += fabs(minf - func.minf);
187     if (fabs(minf - func.minf) > max_err)
188       max_err = fabs(minf - func.minf);
189     printf("vs. global minimum f = %g at x = [", func.minf);
190     for (i = 0; i < func.n; ++i) printf(" %g", func.xmin[i]);
191     printf("]\n");
192
193     val = func.f(func.n, x, NULL, func.f_data);
194     if (val != minf) {
195       fprintf(stderr, "Mismatch %g between returned minf=%g and f(x) = %g\n", 
196               minf - val, minf, val);
197       return 0;
198     }
199   }
200   if (iterations > 1)
201     printf("average #evaluations = %g (%d-%d)\naverage |f-minf| = %g, max |f-minf| = %g\n", total_count * 1.0 / iterations, min_count, max_count, total_err / iterations, max_err);
202
203   free(x);
204   return 1;
205 }
206
207 static void usage(FILE *f)
208 {
209   fprintf(f, "Usage: testopt [OPTIONS]\n"
210           "Options:\n"
211           "     -h : print this help\n"
212           "     -L : list available algorithms and objective functions\n"
213           "     -v : verbose mode\n"
214           " -a <n> : use optimization algorithm <n>\n"
215           " -o <n> : use objective function <n>\n"
216           " -0 <x> : starting guess within <x> + (1+<x>) * optimum\n"
217           "     -c : starting guess at center of cell\n"
218           "     -C : put optimum outside of bound constraints\n"
219           " -e <n> : use at most <n> evals (default: %d, 0 to disable)\n"
220           " -t <t> : use at most <t> seconds (default: disabled)\n"
221           " -x <t> : relative tolerance <t> on x (default: disabled)\n"
222           " -X <t> : absolute tolerance <t> on x (default: disabled)\n"
223           " -f <t> : relative tolerance <t> on f (default: disabled)\n"
224           " -F <t> : absolute tolerance <t> on f (default: disabled)\n"
225           " -m <m> : stop when minf+<m> is reached (default: disabled)\n"
226           " -i <n> : iterate optimization <n> times (default: 1)\n"
227           " -r <s> : use random seed <s> for starting guesses\n"
228           , maxeval);
229 }
230
231 int main(int argc, char **argv)
232 {
233   int c;
234   
235   nlopt_srand_time();
236   testfuncs_verbose = 0;
237   
238   if (argc <= 1)
239     usage(stdout);
240   
241   while ((c = getopt(argc, argv, "hLvCc0:r:a:o:i:e:t:x:X:f:F:m:")) != -1)
242     switch (c) {
243     case 'h':
244       usage(stdout);
245       return EXIT_SUCCESS;
246     case 'L':
247       listalgs(stdout);
248       listfuncs(stdout);
249       return EXIT_SUCCESS;
250     case 'v':
251       testfuncs_verbose = 1;
252       break;
253     case 'C':
254       force_constraints = 1;
255       break;
256     case 'r':
257       nlopt_srand((unsigned long) atoi(optarg));
258       break;
259     case 'a':
260       c = atoi(optarg);
261       if (c < 0 || c >= NLOPT_NUM_ALGORITHMS) {
262         fprintf(stderr, "testopt: invalid algorithm %d\n", c);
263         listalgs(stderr);
264         return EXIT_FAILURE;
265       }
266       algorithm = (nlopt_algorithm) c;
267       break;
268     case 'o':
269       if (!test_function(atoi(optarg)))
270         return EXIT_FAILURE;
271       break;
272     case 'e':
273       maxeval = atoi(optarg);
274       break;
275     case 'i':
276       iterations = atoi(optarg);
277       break;
278     case 't':
279       maxtime = atof(optarg);
280       break;
281     case 'x':
282       xtol_rel = atof(optarg);
283       break;
284     case 'X':
285       xtol_abs = atof(optarg);
286       break;
287     case 'f':
288       ftol_rel = atof(optarg);
289       break;
290     case 'F':
291       ftol_abs = atof(optarg);
292       break;
293     case 'm':
294       minf_max_delta = atof(optarg);
295       break;
296     case 'c':
297       center_start = 1;
298       break;
299     case '0':
300       center_start = 0;
301       xinit_tol = atof(optarg);
302       break;
303     default:
304       fprintf(stderr, "harminv: invalid argument -%c\n", c);
305       usage(stderr);
306       return EXIT_FAILURE;
307     }
308   
309   return EXIT_SUCCESS;
310 }