chiark / gitweb /
3e497532a31a2cf1743cf161e36c3ea3bc4ca752
[nlopt.git] / api / nlopt.c
1 /* Copyright (c) 2007-2008 Massachusetts Institute of Technology
2  *
3  * Permission is hereby granted, free of charge, to any person obtaining
4  * a copy of this software and associated documentation files (the
5  * "Software"), to deal in the Software without restriction, including
6  * without limitation the rights to use, copy, modify, merge, publish,
7  * distribute, sublicense, and/or sell copies of the Software, and to
8  * permit persons to whom the Software is furnished to do so, subject to
9  * the following conditions:
10  * 
11  * The above copyright notice and this permission notice shall be
12  * included in all copies or substantial portions of the Software.
13  * 
14  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND,
15  * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
16  * MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND
17  * NONINFRINGEMENT. IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE
18  * LIABLE FOR ANY CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION
19  * OF CONTRACT, TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION
20  * WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE. 
21  */
22
23 #include <stdlib.h>
24 #include <math.h>
25 #include <string.h>
26 #include <float.h>
27
28 #include "nlopt.h"
29 #include "nlopt-util.h"
30 #include "config.h"
31
32 #define MIN(a,b) ((a) < (b) ? (a) : (b))
33
34 /*************************************************************************/
35
36 #ifdef INFINITY
37 #  define MY_INF INFINITY
38 #else
39 #  define MY_INF HUGE_VAL
40 #endif
41
42 int nlopt_isinf(double x) {
43      return fabs(x) >= HUGE_VAL * 0.99
44 #ifdef HAVE_ISINF
45           || isinf(x)
46 #endif
47           ;
48 }
49
50 #ifndef HAVE_ISNAN
51 static int my_isnan(double x) { return x != x; }
52 #  define isnan my_isnan
53 #endif
54
55 /*************************************************************************/
56
57 void nlopt_version(int *major, int *minor, int *bugfix)
58 {
59      *major = MAJOR_VERSION;
60      *minor = MINOR_VERSION;
61      *bugfix = BUGFIX_VERSION;
62 }
63
64 /*************************************************************************/
65
66 static const char nlopt_algorithm_names[NLOPT_NUM_ALGORITHMS][256] = {
67      "DIRECT (global, no-derivative)",
68      "DIRECT-L (global, no-derivative)",
69      "Randomized DIRECT-L (global, no-derivative)",
70      "Unscaled DIRECT (global, no-derivative)",
71      "Unscaled DIRECT-L (global, no-derivative)",
72      "Unscaled Randomized DIRECT-L (global, no-derivative)",
73      "Original DIRECT version (global, no-derivative)",
74      "Original DIRECT-L version (global, no-derivative)",
75 #ifdef WITH_CXX
76      "StoGO (global, derivative-based)",
77      "StoGO with randomized search (global, derivative-based)",
78 #else
79      "StoGO (NOT COMPILED)",
80      "StoGO randomized (NOT COMPILED)",
81 #endif
82 #ifdef WITH_NOCEDAL_LBFGS
83      "original NON-FREE L-BFGS code by Nocedal et al. (local, deriv.-based)",
84 #else
85      "original NON-FREE L-BFGS code by Nocedal et al. (NOT COMPILED)",
86 #endif
87      "Limited-memory BFGS (L-BFGS) (local, derivative-based)",
88      "Principal-axis, praxis (local, no-derivative)",
89      "Limited-memory variable-metric, rank 1 (local, derivative-based)",
90      "Limited-memory variable-metric, rank 2 (local, derivative-based)",
91      "Truncated Newton (local, derivative-based)",
92      "Truncated Newton with restarting (local, derivative-based)",
93      "Preconditioned truncated Newton (local, derivative-based)",
94      "Preconditioned truncated Newton with restarting (local, derivative-based)",
95      "Controlled random search (CRS2) with local mutation (global, no-derivative)",
96      "Multi-level single-linkage (MLSL), random (global, no-derivative)",
97      "Multi-level single-linkage (MLSL), random (global, derivative)",
98      "Multi-level single-linkage (MLSL), quasi-random (global, no-derivative)",
99      "Multi-level single-linkage (MLSL), quasi-random (global, derivative)",
100      "Method of Moving Asymptotes (MMA) (local, derivative)",
101      "COBYLA (Constrained Optimization BY Linear Approximations) (local, no-derivative)",
102      "NEWUOA unconstrained optimization via quadratic models (local, no-derivative)",
103      "Bound-constrained optimization via NEWUOA-based quadratic models (local, no-derivative)",
104      "Nelder-Mead simplex algorithm (local, no-derivative)",
105      "Sbplx variant of Nelder-Mead (re-implementation of Rowan's Subplex) (local, no-derivative)"
106 };
107
108 const char *nlopt_algorithm_name(nlopt_algorithm a)
109 {
110      if (a < 0 || a >= NLOPT_NUM_ALGORITHMS) return "UNKNOWN";
111      return nlopt_algorithm_names[a];
112 }
113
114 /*************************************************************************/
115
116 static int nlopt_srand_called = 0;
117 void nlopt_srand(unsigned long seed) {
118      nlopt_srand_called = 1;
119      nlopt_init_genrand(seed);
120 }
121
122 void nlopt_srand_time() {
123      nlopt_srand(nlopt_time_seed());
124 }
125
126 /*************************************************************************/
127
128 /* crude heuristics for initial step size of nonderivative algorithms */
129 static double initial_step(int n, 
130                            const double *lb, const double *ub, const double *x)
131 {
132      int i;
133      double step = HUGE_VAL;
134      for (i = 0; i < n; ++i) {
135           if (!nlopt_isinf(ub[i]) && !nlopt_isinf(lb[i])
136               && (ub[i] - lb[i]) * 0.25 < step && ub[i] > lb[i])
137                step = (ub[i] - lb[i]) * 0.25;
138           if (!nlopt_isinf(ub[i]) 
139               && ub[i] - x[i] < step && ub[i] > x[i])
140                step = ub[i] - x[i];
141           if (!nlopt_isinf(lb[i]) 
142               && x[i] - lb[i] < step && x[i] > lb[i])
143                step = x[i] - lb[i];
144      }
145      if (nlopt_isinf(step))
146           for (i = 0; i < n; ++i) {
147                if (!nlopt_isinf(ub[i]) 
148                    && ub[i] - x[i] < step && ub[i] > x[i] + 1e-4)
149                     step = ub[i] - x[i];
150                if (!nlopt_isinf(lb[i]) 
151                    && x[i] - lb[i] < step && x[i] > lb[i] + 1e-4)
152                     step = x[i] - lb[i];
153           }
154      if (nlopt_isinf(step)) {
155           step = 0;
156           for (i = 0; i < n; ++i)
157                if (fabs(x[i]) * 0.25 > step)
158                     step = fabs(x[i]) * 0.25;
159           if (step == 0)
160                step = 1;
161      }
162      return step;
163 }
164               
165 /*************************************************************************/
166
167 typedef struct {
168      nlopt_func f;
169      void *f_data;
170      const double *lb, *ub;
171 } nlopt_data;
172
173 #include "praxis.h"
174
175 static double f_bound(int n, const double *x, void *data_)
176 {
177      int i;
178      nlopt_data *data = (nlopt_data *) data_;
179      double f;
180
181      /* some methods do not support bound constraints, but support
182         discontinuous objectives so we can just return Inf for invalid x */
183      for (i = 0; i < n; ++i)
184           if (x[i] < data->lb[i] || x[i] > data->ub[i])
185                return MY_INF;
186
187      f = data->f(n, x, NULL, data->f_data);
188      return (isnan(f) || nlopt_isinf(f) ? MY_INF : f);
189 }
190
191 static double f_noderiv(int n, const double *x, void *data_)
192 {
193      nlopt_data *data = (nlopt_data *) data_;
194      return data->f(n, x, NULL, data->f_data);
195 }
196
197 #include "direct.h"
198
199 static double f_direct(int n, const double *x, int *undefined, void *data_)
200 {
201      nlopt_data *data = (nlopt_data *) data_;
202      double f;
203      f = data->f(n, x, NULL, data->f_data);
204      *undefined = isnan(f) || nlopt_isinf(f);
205      return f;
206 }
207
208 #ifdef WITH_CXX
209 #  include "stogo.h"
210 #endif
211
212 #include "cdirect.h"
213
214 #ifdef WITH_NOCEDAL
215 #  include "l-bfgs-b.h"
216 #endif
217
218 #include "luksan.h"
219
220 #include "crs.h"
221
222 #include "mlsl.h"
223 #include "mma.h"
224 #include "cobyla.h"
225 #include "newuoa.h"
226 #include "neldermead.h"
227
228 /*************************************************************************/
229
230 /* for "hybrid" algorithms that combine global search with some
231    local search algorithm, most of the time we anticipate that people
232    will stick with the default local search algorithm, so we
233    don't add this as a parameter to nlopt_minimize.  Instead, the user
234    can call nlopt_{set/get}_hybrid_local_algorithm to get/set the defaults. */
235
236 /* default local-search algorithms */
237 static nlopt_algorithm local_search_alg_deriv = NLOPT_LD_MMA;
238 static nlopt_algorithm local_search_alg_nonderiv = NLOPT_LN_COBYLA;
239
240 static int local_search_maxeval = -1; /* no maximum by default */
241
242 void nlopt_get_local_search_algorithm(nlopt_algorithm *deriv,
243                                       nlopt_algorithm *nonderiv,
244                                       int *maxeval)
245 {
246      *deriv = local_search_alg_deriv;
247      *nonderiv = local_search_alg_nonderiv;
248      *maxeval = local_search_maxeval;
249 }
250
251 void nlopt_set_local_search_algorithm(nlopt_algorithm deriv,
252                                       nlopt_algorithm nonderiv,
253                                       int maxeval)
254 {
255      local_search_alg_deriv = deriv;
256      local_search_alg_nonderiv = nonderiv;
257      local_search_maxeval = maxeval;
258 }
259
260 /*************************************************************************/
261
262 /* same as nlopt_minimize, but xtol_abs is required to be non-NULL */
263 static nlopt_result nlopt_minimize_(
264      nlopt_algorithm algorithm,
265      int n, nlopt_func f, void *f_data,
266      int m, nlopt_func fc, void *fc_data, ptrdiff_t fc_datum_size,
267      const double *lb, const double *ub, /* bounds */
268      double *x, /* in: initial guess, out: minimizer */
269      double *minf, /* out: minimum */
270      double minf_max, double ftol_rel, double ftol_abs,
271      double xtol_rel, const double *xtol_abs,
272      int maxeval, double maxtime)
273 {
274      int i;
275      nlopt_data d;
276      nlopt_stopping stop;
277
278      /* some basic argument checks */
279      if (!minf || !f || n < 0 || m < 0
280           || (m > 0 && !fc)) return NLOPT_INVALID_ARGS;
281      if (n == 0) { /* trivial case: no degrees of freedom */
282           *minf = f(n, x, NULL, f_data);
283           return NLOPT_SUCCESS;
284      }
285      else if (n < 0 || !lb || !ub || !x)
286           return NLOPT_INVALID_ARGS;
287
288      /* nonlinear constraints are only supported with MMA or COBYLA */
289      if (m != 0 && algorithm != NLOPT_LD_MMA && algorithm != NLOPT_LN_COBYLA) 
290           return NLOPT_INVALID_ARGS;
291
292      d.f = f;
293      d.f_data = f_data;
294      d.lb = lb;
295      d.ub = ub;
296
297      /* make sure rand generator is inited */
298      if (!nlopt_srand_called)
299           nlopt_srand_time(); /* default is non-deterministic */
300
301      /* check bound constraints */
302      for (i = 0; i < n; ++i)
303           if (lb[i] > ub[i] || x[i] < lb[i] || x[i] > ub[i])
304                return NLOPT_INVALID_ARGS;
305
306      stop.n = n;
307      stop.minf_max = (isnan(minf_max) 
308                       || (nlopt_isinf(minf_max) && minf_max < 0))
309           ? -MY_INF : minf_max;
310      stop.ftol_rel = ftol_rel;
311      stop.ftol_abs = ftol_abs;
312      stop.xtol_rel = xtol_rel;
313      stop.xtol_abs = xtol_abs;
314      stop.nevals = 0;
315      stop.maxeval = maxeval;
316      stop.maxtime = maxtime;
317      stop.start = nlopt_seconds();
318
319      switch (algorithm) {
320          case NLOPT_GN_DIRECT:
321          case NLOPT_GN_DIRECT_L: 
322          case NLOPT_GN_DIRECT_L_RAND: 
323               return cdirect(n, f, f_data, lb, ub, x, minf, &stop, 0.0, 
324                              (algorithm != NLOPT_GN_DIRECT)
325                              + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
326                              + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
327               
328          case NLOPT_GN_DIRECT_NOSCAL:
329          case NLOPT_GN_DIRECT_L_NOSCAL: 
330          case NLOPT_GN_DIRECT_L_RAND_NOSCAL: 
331               return cdirect_unscaled(n, f, f_data, lb, ub, x, minf, 
332                                       &stop, 0.0, 
333                                       (algorithm != NLOPT_GN_DIRECT)
334                                       + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
335                                       + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
336               
337          case NLOPT_GN_ORIG_DIRECT:
338          case NLOPT_GN_ORIG_DIRECT_L: 
339               switch (direct_optimize(f_direct, &d, n, lb, ub, x, minf,
340                                       maxeval, -1, 0.0, 0.0,
341                                       pow(xtol_rel, (double) n), -1.0,
342                                       stop.minf_max, 0.0,
343                                       NULL, 
344                                       algorithm == NLOPT_GN_ORIG_DIRECT
345                                       ? DIRECT_ORIGINAL
346                                       : DIRECT_GABLONSKY)) {
347                   case DIRECT_INVALID_BOUNDS:
348                   case DIRECT_MAXFEVAL_TOOBIG:
349                   case DIRECT_INVALID_ARGS:
350                        return NLOPT_INVALID_ARGS;
351                   case DIRECT_INIT_FAILED:
352                   case DIRECT_SAMPLEPOINTS_FAILED:
353                   case DIRECT_SAMPLE_FAILED:
354                        return NLOPT_FAILURE;
355                   case DIRECT_MAXFEVAL_EXCEEDED:
356                   case DIRECT_MAXITER_EXCEEDED:
357                        return NLOPT_MAXEVAL_REACHED;
358                   case DIRECT_GLOBAL_FOUND:
359                        return NLOPT_MINF_MAX_REACHED;
360                   case DIRECT_VOLTOL:
361                   case DIRECT_SIGMATOL:
362                        return NLOPT_XTOL_REACHED;
363                   case DIRECT_OUT_OF_MEMORY:
364                        return NLOPT_OUT_OF_MEMORY;
365               break;
366          }
367
368          case NLOPT_GD_STOGO:
369          case NLOPT_GD_STOGO_RAND:
370 #ifdef WITH_CXX
371               if (!stogo_minimize(n, f, f_data, x, minf, lb, ub, &stop,
372                                   algorithm == NLOPT_GD_STOGO
373                                   ? 0 : 2*n))
374                    return NLOPT_FAILURE;
375               break;
376 #else
377               return NLOPT_FAILURE;
378 #endif
379
380 #if 0
381               /* lacking a free/open-source license, we no longer use
382                  Rowan's code, and instead use by "sbplx" re-implementation */
383          case NLOPT_LN_SUBPLEX: {
384               int iret;
385               double *scale = (double *) malloc(sizeof(double) * n);
386               if (!scale) return NLOPT_OUT_OF_MEMORY;
387               for (i = 0; i < n; ++i)
388                    scale[i] = initial_step(1, lb+i, ub+i, x+i);
389               iret = nlopt_subplex(f_bound, minf, x, n, &d, &stop, scale);
390               free(scale);
391               switch (iret) {
392                   case -2: return NLOPT_INVALID_ARGS;
393                   case -10: return NLOPT_MAXTIME_REACHED;
394                   case -1: return NLOPT_MAXEVAL_REACHED;
395                   case 0: return NLOPT_XTOL_REACHED;
396                   case 1: return NLOPT_SUCCESS;
397                   case 2: return NLOPT_MINF_MAX_REACHED;
398                   case 20: return NLOPT_FTOL_REACHED;
399                   case -200: return NLOPT_OUT_OF_MEMORY;
400                   default: return NLOPT_FAILURE; /* unknown return code */
401               }
402               break;
403          }
404 #endif
405
406          case NLOPT_LN_PRAXIS:
407               return praxis_(0.0, DBL_EPSILON, 
408                              initial_step(n, lb, ub, x), n, x, f_bound, &d,
409                              &stop, minf);
410
411 #ifdef WITH_NOCEDAL
412          case NLOPT_LD_LBFGS_NOCEDAL: {
413               int iret, *nbd = (int *) malloc(sizeof(int) * n);
414               if (!nbd) return NLOPT_OUT_OF_MEMORY;
415               for (i = 0; i < n; ++i) {
416                    int linf = nlopt_isinf(lb[i]) && lb[i] < 0;
417                    int uinf = nlopt_isinf(ub[i]) && ub[i] > 0;
418                    nbd[i] = linf && uinf ? 0 : (uinf ? 1 : (linf ? 3 : 2));
419               }
420               iret = lbfgsb_minimize(n, f, f_data, x, nbd, lb, ub,
421                                      MIN(n, 5), 0.0, ftol_rel, 
422                                      xtol_abs ? *xtol_abs : xtol_rel,
423                                      maxeval);
424               free(nbd);
425               if (iret <= 0) {
426                    switch (iret) {
427                        case -1: return NLOPT_INVALID_ARGS;
428                        case -2: default: return NLOPT_FAILURE;
429                    }
430               }
431               else {
432                    *minf = f(n, x, NULL, f_data);
433                    switch (iret) {
434                        case 5: return NLOPT_MAXEVAL_REACHED;
435                        case 2: return NLOPT_XTOL_REACHED;
436                        case 1: return NLOPT_FTOL_REACHED;
437                        default: return NLOPT_SUCCESS;
438                    }
439               }
440               break;
441          }
442 #endif
443
444          case NLOPT_LD_LBFGS: 
445               return luksan_plis(n, f, f_data, lb, ub, x, minf, &stop);
446
447          case NLOPT_LD_VAR1: 
448          case NLOPT_LD_VAR2: 
449               return luksan_plip(n, f, f_data, lb, ub, x, minf, &stop,
450                    algorithm == NLOPT_LD_VAR1 ? 1 : 2);
451
452          case NLOPT_LD_TNEWTON: 
453          case NLOPT_LD_TNEWTON_RESTART: 
454          case NLOPT_LD_TNEWTON_PRECOND: 
455          case NLOPT_LD_TNEWTON_PRECOND_RESTART: 
456               return luksan_pnet(n, f, f_data, lb, ub, x, minf, &stop,
457                                  1 + (algorithm - NLOPT_LD_TNEWTON) % 2,
458                                  1 + (algorithm - NLOPT_LD_TNEWTON) / 2);
459
460          case NLOPT_GN_CRS2_LM:
461               return crs_minimize(n, f, f_data, lb, ub, x, minf, &stop, 0);
462
463          case NLOPT_GN_MLSL:
464          case NLOPT_GD_MLSL:
465          case NLOPT_GN_MLSL_LDS:
466          case NLOPT_GD_MLSL_LDS:
467               for (i = 0; i < n && stop.xtol_abs[i] <= 0; ++i) ;
468               if (stop.ftol_rel <= 0 && stop.ftol_abs <= 0 &&
469                   stop.xtol_rel <= 0 && i == n) {
470                    /* it is not sensible to call MLSL without *some*
471                       nonzero tolerance for the local search */
472                    stop.ftol_rel = 1e-15;
473                    stop.xtol_rel = 1e-7;
474               }
475               return mlsl_minimize(n, f, f_data, lb, ub, x, minf, &stop,
476                                    (algorithm == NLOPT_GN_MLSL ||
477                                     algorithm == NLOPT_GN_MLSL_LDS)
478                                    ? local_search_alg_nonderiv
479                                    : local_search_alg_deriv,
480                                    local_search_maxeval,
481                                    algorithm >= NLOPT_GN_MLSL_LDS);
482
483          case NLOPT_LD_MMA:
484               return mma_minimize(n, f, f_data, m, fc, fc_data, fc_datum_size,
485                                   lb, ub, x, minf, &stop,
486                                   local_search_alg_deriv, 1e-12, 100000);
487
488          case NLOPT_LN_COBYLA:
489               return cobyla_minimize(n, f, f_data, 
490                                      m, fc, fc_data, fc_datum_size,
491                                      lb, ub, x, minf, &stop,
492                                      initial_step(n, lb, ub, x));
493                                      
494          case NLOPT_LN_NEWUOA:
495               return newuoa(n, 2*n+1, x, 0, 0, initial_step(n, lb, ub, x),
496                             &stop, minf, f_noderiv, &d);
497                                      
498          case NLOPT_LN_NEWUOA_BOUND:
499               return newuoa(n, 2*n+1, x, lb, ub, initial_step(n, lb, ub, x),
500                             &stop, minf, f_noderiv, &d);
501
502          case NLOPT_LN_NELDERMEAD: 
503          case NLOPT_LN_SBPLX: 
504          {
505               nlopt_result ret;
506               double *scale = (double *) malloc(sizeof(double) * n);
507               if (!scale) return NLOPT_OUT_OF_MEMORY;
508               for (i = 0; i < n; ++i)
509                    scale[i] = initial_step(1, lb+i, ub+i, x+i);
510               if (algorithm == NLOPT_LN_NELDERMEAD)
511                    ret = nldrmd_minimize(n,f,f_data,lb,ub,x,minf,scale,&stop);
512               else
513                    ret = sbplx_minimize(n,f,f_data,lb,ub,x,minf,scale,&stop);
514               free(scale);
515               return ret;
516          }
517
518          default:
519               return NLOPT_INVALID_ARGS;
520      }
521
522      return NLOPT_SUCCESS;
523 }
524
525 nlopt_result nlopt_minimize_constrained(
526      nlopt_algorithm algorithm,
527      int n, nlopt_func f, void *f_data,
528      int m, nlopt_func fc, void *fc_data, ptrdiff_t fc_datum_size,
529      const double *lb, const double *ub, /* bounds */
530      double *x, /* in: initial guess, out: minimizer */
531      double *minf, /* out: minimum */
532      double minf_max, double ftol_rel, double ftol_abs,
533      double xtol_rel, const double *xtol_abs,
534      int maxeval, double maxtime)
535 {
536      nlopt_result ret;
537      if (xtol_abs)
538           ret = nlopt_minimize_(algorithm, n, f, f_data,
539                                 m, fc, fc_data, fc_datum_size, lb, ub,
540                                 x, minf, minf_max, ftol_rel, ftol_abs,
541                                 xtol_rel, xtol_abs, maxeval, maxtime);
542      else {
543           int i;
544           double *xtol = (double *) malloc(sizeof(double) * n);
545           if (!xtol) return NLOPT_OUT_OF_MEMORY;
546           for (i = 0; i < n; ++i) xtol[i] = -1;
547           ret = nlopt_minimize_(algorithm, n, f, f_data, 
548                                 m, fc, fc_data, fc_datum_size, lb, ub,
549                                 x, minf, minf_max, ftol_rel, ftol_abs,
550                                 xtol_rel, xtol, maxeval, maxtime);
551           free(xtol);
552      }
553      return ret;
554 }
555
556 nlopt_result nlopt_minimize(
557      nlopt_algorithm algorithm,
558      int n, nlopt_func f, void *f_data,
559      const double *lb, const double *ub, /* bounds */
560      double *x, /* in: initial guess, out: minimizer */
561      double *minf, /* out: minimum */
562      double minf_max, double ftol_rel, double ftol_abs,
563      double xtol_rel, const double *xtol_abs,
564      int maxeval, double maxtime)
565 {
566      return nlopt_minimize_constrained(
567           algorithm, n, f, f_data, 0, NULL, NULL, 0,
568           lb, ub, x, minf, minf_max, ftol_rel, ftol_abs,
569           xtol_rel, xtol_abs, maxeval, maxtime);
570 }