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