chiark / gitweb /
add preliminary BOBYQA (not quite working yet)
[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      "Augmented Lagrangian method (local, no-derivative)",
100      "Augmented Lagrangian method (local, derivative)",
101      "Augmented Lagrangian method for equality constraints (local, no-derivative)",
102      "Augmented Lagrangian method for equality constraints (local, derivative)",
103      "BOBYQA bound-constrained optimization via quadratic models (local, no-derivative)",
104 };
105
106 const char *nlopt_algorithm_name(nlopt_algorithm a)
107 {
108      if (a < 0 || a >= NLOPT_NUM_ALGORITHMS) return "UNKNOWN";
109      return nlopt_algorithm_names[a];
110 }
111
112 /*************************************************************************/
113
114 static int nlopt_srand_called = 0;
115 void nlopt_srand(unsigned long seed) {
116      nlopt_srand_called = 1;
117      nlopt_init_genrand(seed);
118 }
119
120 void nlopt_srand_time() {
121      nlopt_srand(nlopt_time_seed());
122 }
123
124 /*************************************************************************/
125
126 /* crude heuristics for initial step size of nonderivative algorithms */
127 static double initial_step(int n, 
128                            const double *lb, const double *ub, const double *x)
129 {
130      int i;
131      double step = HUGE_VAL;
132      for (i = 0; i < n; ++i) {
133           if (!nlopt_isinf(ub[i]) && !nlopt_isinf(lb[i])
134               && (ub[i] - lb[i]) * 0.25 < step && ub[i] > lb[i])
135                step = (ub[i] - lb[i]) * 0.25;
136           if (!nlopt_isinf(ub[i]) 
137               && ub[i] - x[i] < step && ub[i] > x[i])
138                step = ub[i] - x[i];
139           if (!nlopt_isinf(lb[i]) 
140               && x[i] - lb[i] < step && x[i] > lb[i])
141                step = x[i] - lb[i];
142      }
143      if (nlopt_isinf(step))
144           for (i = 0; i < n; ++i) {
145                if (!nlopt_isinf(ub[i]) 
146                    && ub[i] - x[i] < step && ub[i] > x[i] + 1e-4)
147                     step = ub[i] - x[i];
148                if (!nlopt_isinf(lb[i]) 
149                    && x[i] - lb[i] < step && x[i] > lb[i] + 1e-4)
150                     step = x[i] - lb[i];
151           }
152      if (nlopt_isinf(step)) {
153           step = 0;
154           for (i = 0; i < n; ++i)
155                if (fabs(x[i]) * 0.25 > step)
156                     step = fabs(x[i]) * 0.25;
157           if (step == 0)
158                step = 1;
159      }
160      return step;
161 }
162               
163 /*************************************************************************/
164
165 typedef struct {
166      nlopt_func f;
167      void *f_data;
168      const double *lb, *ub;
169 } nlopt_data;
170
171 #include "praxis.h"
172
173 static double f_bound(int n, const double *x, void *data_)
174 {
175      int i;
176      nlopt_data *data = (nlopt_data *) data_;
177      double f;
178
179      /* some methods do not support bound constraints, but support
180         discontinuous objectives so we can just return Inf for invalid x */
181      for (i = 0; i < n; ++i)
182           if (x[i] < data->lb[i] || x[i] > data->ub[i])
183                return HUGE_VAL;
184
185      f = data->f(n, x, NULL, data->f_data);
186      return (isnan(f) || nlopt_isinf(f) ? HUGE_VAL : f);
187 }
188
189 static double f_noderiv(int n, const double *x, void *data_)
190 {
191      nlopt_data *data = (nlopt_data *) data_;
192      return data->f(n, x, NULL, data->f_data);
193 }
194
195 #include "direct.h"
196
197 static double f_direct(int n, const double *x, int *undefined, void *data_)
198 {
199      nlopt_data *data = (nlopt_data *) data_;
200      double f;
201      f = data->f(n, x, NULL, data->f_data);
202      *undefined = isnan(f) || nlopt_isinf(f);
203      return f;
204 }
205
206 #ifdef WITH_CXX
207 #  include "stogo.h"
208 #endif
209
210 #include "cdirect.h"
211
212 #ifdef WITH_NOCEDAL
213 #  include "l-bfgs-b.h"
214 #endif
215
216 #include "luksan.h"
217
218 #include "crs.h"
219
220 #include "mlsl.h"
221 #include "mma.h"
222 #include "cobyla.h"
223 #include "newuoa.h"
224 #include "neldermead.h"
225 #include "auglag.h"
226 #include "bobyqa.h"
227
228 #define AUGLAG_ALG(a) ((a) == NLOPT_LN_AUGLAG ||        \
229                        (a) == NLOPT_LN_AUGLAG_EQ ||     \
230                        (a) == NLOPT_LD_AUGLAG ||        \
231                        (a) == NLOPT_LD_AUGLAG_EQ)
232
233 /*************************************************************************/
234
235 /* for "hybrid" algorithms that combine global search with some
236    local search algorithm, most of the time we anticipate that people
237    will stick with the default local search algorithm, so we
238    don't add this as a parameter to nlopt_minimize.  Instead, the user
239    can call nlopt_{set/get}_hybrid_local_algorithm to get/set the defaults. */
240
241 /* default local-search algorithms */
242 static nlopt_algorithm local_search_alg_deriv = NLOPT_LD_MMA;
243 static nlopt_algorithm local_search_alg_nonderiv = NLOPT_LN_COBYLA;
244
245 static int local_search_maxeval = -1; /* no maximum by default */
246
247 void nlopt_get_local_search_algorithm(nlopt_algorithm *deriv,
248                                       nlopt_algorithm *nonderiv,
249                                       int *maxeval)
250 {
251      *deriv = local_search_alg_deriv;
252      *nonderiv = local_search_alg_nonderiv;
253      *maxeval = local_search_maxeval;
254 }
255
256 void nlopt_set_local_search_algorithm(nlopt_algorithm deriv,
257                                       nlopt_algorithm nonderiv,
258                                       int maxeval)
259 {
260      local_search_alg_deriv = deriv;
261      local_search_alg_nonderiv = nonderiv;
262      local_search_maxeval = maxeval;
263 }
264
265 /*************************************************************************/
266
267 /* same as nlopt_minimize_econstrained, 
268    but xtol_abs is required to be non-NULL */
269 static nlopt_result nlopt_minimize_(
270      nlopt_algorithm algorithm,
271      int n, nlopt_func f, void *f_data,
272      int m, nlopt_func fc, void *fc_data, ptrdiff_t fc_datum_size,
273      int p, nlopt_func h, void *h_data, ptrdiff_t h_datum_size,
274      const double *lb, const double *ub, /* bounds */
275      double *x, /* in: initial guess, out: minimizer */
276      double *minf, /* out: minimum */
277      double minf_max, double ftol_rel, double ftol_abs,
278      double xtol_rel, const double *xtol_abs,
279      double htol_rel, double htol_abs,
280      int maxeval, double maxtime)
281 {
282      int i;
283      nlopt_data d;
284      nlopt_stopping stop;
285
286      /* some basic argument checks */
287      if (!minf || !f || n < 0 || m < 0 || p < 0 || (p > 0 && !h)
288           || (m > 0 && !fc)) return NLOPT_INVALID_ARGS;
289      if (n == 0) { /* trivial case: no degrees of freedom */
290           *minf = f(n, x, NULL, f_data);
291           return NLOPT_SUCCESS;
292      }
293      else if (n < 0 || !lb || !ub || !x)
294           return NLOPT_INVALID_ARGS;
295
296      /* nonlinear constraints are only supported with some algorithms */
297      if (m != 0 && algorithm != NLOPT_LD_MMA && algorithm != NLOPT_LN_COBYLA
298           && !AUGLAG_ALG(algorithm)) 
299           return NLOPT_INVALID_ARGS;
300
301      /* nonlinear equality constraints (h(x) = 0) only via some algorithms */
302      if (p != 0 && !AUGLAG_ALG(algorithm))
303           return NLOPT_INVALID_ARGS;
304
305      d.f = f;
306      d.f_data = f_data;
307      d.lb = lb;
308      d.ub = ub;
309
310      /* make sure rand generator is inited */
311      if (!nlopt_srand_called)
312           nlopt_srand_time(); /* default is non-deterministic */
313
314      /* check bound constraints */
315      for (i = 0; i < n; ++i)
316           if (lb[i] > ub[i] || x[i] < lb[i] || x[i] > ub[i])
317                return NLOPT_INVALID_ARGS;
318
319      stop.n = n;
320      stop.minf_max = (isnan(minf_max) 
321                       || (nlopt_isinf(minf_max) && minf_max < 0))
322           ? -HUGE_VAL : minf_max;
323      stop.ftol_rel = ftol_rel;
324      stop.ftol_abs = ftol_abs;
325      stop.xtol_rel = xtol_rel;
326      stop.xtol_abs = xtol_abs;
327      stop.htol_rel = htol_rel;
328      stop.htol_abs = htol_abs;
329      stop.nevals = 0;
330      stop.maxeval = maxeval;
331      stop.maxtime = maxtime;
332      stop.start = nlopt_seconds();
333
334      switch (algorithm) {
335          case NLOPT_GN_DIRECT:
336          case NLOPT_GN_DIRECT_L: 
337          case NLOPT_GN_DIRECT_L_RAND: 
338               return cdirect(n, f, f_data, lb, ub, x, minf, &stop, 0.0, 
339                              (algorithm != NLOPT_GN_DIRECT)
340                              + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
341                              + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
342               
343          case NLOPT_GN_DIRECT_NOSCAL:
344          case NLOPT_GN_DIRECT_L_NOSCAL: 
345          case NLOPT_GN_DIRECT_L_RAND_NOSCAL: 
346               return cdirect_unscaled(n, f, f_data, lb, ub, x, minf, 
347                                       &stop, 0.0, 
348                                       (algorithm != NLOPT_GN_DIRECT)
349                                       + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
350                                       + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
351               
352          case NLOPT_GN_ORIG_DIRECT:
353          case NLOPT_GN_ORIG_DIRECT_L: 
354               switch (direct_optimize(f_direct, &d, n, lb, ub, x, minf,
355                                       maxeval, -1, 0.0, 0.0,
356                                       pow(xtol_rel, (double) n), -1.0,
357                                       stop.minf_max, 0.0,
358                                       NULL, 
359                                       algorithm == NLOPT_GN_ORIG_DIRECT
360                                       ? DIRECT_ORIGINAL
361                                       : DIRECT_GABLONSKY)) {
362                   case DIRECT_INVALID_BOUNDS:
363                   case DIRECT_MAXFEVAL_TOOBIG:
364                   case DIRECT_INVALID_ARGS:
365                        return NLOPT_INVALID_ARGS;
366                   case DIRECT_INIT_FAILED:
367                   case DIRECT_SAMPLEPOINTS_FAILED:
368                   case DIRECT_SAMPLE_FAILED:
369                        return NLOPT_FAILURE;
370                   case DIRECT_MAXFEVAL_EXCEEDED:
371                   case DIRECT_MAXITER_EXCEEDED:
372                        return NLOPT_MAXEVAL_REACHED;
373                   case DIRECT_GLOBAL_FOUND:
374                        return NLOPT_MINF_MAX_REACHED;
375                   case DIRECT_VOLTOL:
376                   case DIRECT_SIGMATOL:
377                        return NLOPT_XTOL_REACHED;
378                   case DIRECT_OUT_OF_MEMORY:
379                        return NLOPT_OUT_OF_MEMORY;
380               break;
381          }
382
383          case NLOPT_GD_STOGO:
384          case NLOPT_GD_STOGO_RAND:
385 #ifdef WITH_CXX
386               if (!stogo_minimize(n, f, f_data, x, minf, lb, ub, &stop,
387                                   algorithm == NLOPT_GD_STOGO
388                                   ? 0 : 2*n))
389                    return NLOPT_FAILURE;
390               break;
391 #else
392               return NLOPT_FAILURE;
393 #endif
394
395 #if 0
396               /* lacking a free/open-source license, we no longer use
397                  Rowan's code, and instead use by "sbplx" re-implementation */
398          case NLOPT_LN_SUBPLEX: {
399               int iret;
400               double *scale = (double *) malloc(sizeof(double) * n);
401               if (!scale) return NLOPT_OUT_OF_MEMORY;
402               for (i = 0; i < n; ++i)
403                    scale[i] = initial_step(1, lb+i, ub+i, x+i);
404               iret = nlopt_subplex(f_bound, minf, x, n, &d, &stop, scale);
405               free(scale);
406               switch (iret) {
407                   case -2: return NLOPT_INVALID_ARGS;
408                   case -10: return NLOPT_MAXTIME_REACHED;
409                   case -1: return NLOPT_MAXEVAL_REACHED;
410                   case 0: return NLOPT_XTOL_REACHED;
411                   case 1: return NLOPT_SUCCESS;
412                   case 2: return NLOPT_MINF_MAX_REACHED;
413                   case 20: return NLOPT_FTOL_REACHED;
414                   case -200: return NLOPT_OUT_OF_MEMORY;
415                   default: return NLOPT_FAILURE; /* unknown return code */
416               }
417               break;
418          }
419 #endif
420
421          case NLOPT_LN_PRAXIS:
422               return praxis_(0.0, DBL_EPSILON, 
423                              initial_step(n, lb, ub, x), n, x, f_bound, &d,
424                              &stop, minf);
425
426 #ifdef WITH_NOCEDAL
427          case NLOPT_LD_LBFGS_NOCEDAL: {
428               int iret, *nbd = (int *) malloc(sizeof(int) * n);
429               if (!nbd) return NLOPT_OUT_OF_MEMORY;
430               for (i = 0; i < n; ++i) {
431                    int linf = nlopt_isinf(lb[i]) && lb[i] < 0;
432                    int uinf = nlopt_isinf(ub[i]) && ub[i] > 0;
433                    nbd[i] = linf && uinf ? 0 : (uinf ? 1 : (linf ? 3 : 2));
434               }
435               iret = lbfgsb_minimize(n, f, f_data, x, nbd, lb, ub,
436                                      MIN(n, 5), 0.0, ftol_rel, 
437                                      xtol_abs ? *xtol_abs : xtol_rel,
438                                      maxeval);
439               free(nbd);
440               if (iret <= 0) {
441                    switch (iret) {
442                        case -1: return NLOPT_INVALID_ARGS;
443                        case -2: default: return NLOPT_FAILURE;
444                    }
445               }
446               else {
447                    *minf = f(n, x, NULL, f_data);
448                    switch (iret) {
449                        case 5: return NLOPT_MAXEVAL_REACHED;
450                        case 2: return NLOPT_XTOL_REACHED;
451                        case 1: return NLOPT_FTOL_REACHED;
452                        default: return NLOPT_SUCCESS;
453                    }
454               }
455               break;
456          }
457 #endif
458
459          case NLOPT_LD_LBFGS: 
460               return luksan_plis(n, f, f_data, lb, ub, x, minf, &stop);
461
462          case NLOPT_LD_VAR1: 
463          case NLOPT_LD_VAR2: 
464               return luksan_plip(n, f, f_data, lb, ub, x, minf, &stop,
465                    algorithm == NLOPT_LD_VAR1 ? 1 : 2);
466
467          case NLOPT_LD_TNEWTON: 
468          case NLOPT_LD_TNEWTON_RESTART: 
469          case NLOPT_LD_TNEWTON_PRECOND: 
470          case NLOPT_LD_TNEWTON_PRECOND_RESTART: 
471               return luksan_pnet(n, f, f_data, lb, ub, x, minf, &stop,
472                                  1 + (algorithm - NLOPT_LD_TNEWTON) % 2,
473                                  1 + (algorithm - NLOPT_LD_TNEWTON) / 2);
474
475          case NLOPT_GN_CRS2_LM:
476               return crs_minimize(n, f, f_data, lb, ub, x, minf, &stop, 0);
477
478          case NLOPT_GN_MLSL:
479          case NLOPT_GD_MLSL:
480          case NLOPT_GN_MLSL_LDS:
481          case NLOPT_GD_MLSL_LDS:
482               for (i = 0; i < n && stop.xtol_abs[i] <= 0; ++i) ;
483               if (stop.ftol_rel <= 0 && stop.ftol_abs <= 0 &&
484                   stop.xtol_rel <= 0 && i == n) {
485                    /* it is not sensible to call MLSL without *some*
486                       nonzero tolerance for the local search */
487                    stop.ftol_rel = 1e-15;
488                    stop.xtol_rel = 1e-7;
489               }
490               return mlsl_minimize(n, f, f_data, lb, ub, x, minf, &stop,
491                                    (algorithm == NLOPT_GN_MLSL ||
492                                     algorithm == NLOPT_GN_MLSL_LDS)
493                                    ? local_search_alg_nonderiv
494                                    : local_search_alg_deriv,
495                                    local_search_maxeval,
496                                    algorithm >= NLOPT_GN_MLSL_LDS);
497
498          case NLOPT_LD_MMA:
499               return mma_minimize(n, f, f_data, m, fc, fc_data, fc_datum_size,
500                                   lb, ub, x, minf, &stop,
501                                   local_search_alg_deriv, 1e-12, 100000);
502
503          case NLOPT_LN_COBYLA:
504               return cobyla_minimize(n, f, f_data, 
505                                      m, fc, fc_data, fc_datum_size,
506                                      lb, ub, x, minf, &stop,
507                                      initial_step(n, lb, ub, x));
508                                      
509          case NLOPT_LN_NEWUOA:
510               return newuoa(n, 2*n+1, x, 0, 0, initial_step(n, lb, ub, x),
511                             &stop, minf, f_noderiv, &d);
512                                      
513          case NLOPT_LN_NEWUOA_BOUND:
514               return newuoa(n, 2*n+1, x, lb, ub, initial_step(n, lb, ub, x),
515                             &stop, minf, f_noderiv, &d);
516
517          case NLOPT_LN_BOBYQA:
518               return bobyqa(n, 2*n+1, x, lb, ub, initial_step(n, lb, ub, x),
519                             &stop, minf, f_noderiv, &d);
520
521          case NLOPT_LN_NELDERMEAD: 
522          case NLOPT_LN_SBPLX: 
523          {
524               nlopt_result ret;
525               double *scale = (double *) malloc(sizeof(double) * n);
526               if (!scale) return NLOPT_OUT_OF_MEMORY;
527               for (i = 0; i < n; ++i)
528                    scale[i] = initial_step(1, lb+i, ub+i, x+i);
529               if (algorithm == NLOPT_LN_NELDERMEAD)
530                    ret = nldrmd_minimize(n,f,f_data,lb,ub,x,minf,scale,&stop);
531               else
532                    ret = sbplx_minimize(n,f,f_data,lb,ub,x,minf,scale,&stop);
533               free(scale);
534               return ret;
535          }
536
537          case NLOPT_LN_AUGLAG:
538          case NLOPT_LN_AUGLAG_EQ:
539               return auglag_minimize(n, f, f_data, 
540                                      m, fc, fc_data, fc_datum_size,
541                                      p, h, h_data, h_datum_size,
542                                      lb, ub, x, minf, &stop,
543                                      local_search_alg_nonderiv,
544                                      algorithm == NLOPT_LN_AUGLAG_EQ);
545
546          case NLOPT_LD_AUGLAG:
547          case NLOPT_LD_AUGLAG_EQ:
548               return auglag_minimize(n, f, f_data, 
549                                      m, fc, fc_data, fc_datum_size,
550                                      p, h, h_data, h_datum_size,
551                                      lb, ub, x, minf, &stop,
552                                      local_search_alg_deriv,
553                                      algorithm == NLOPT_LD_AUGLAG_EQ);
554
555          default:
556               return NLOPT_INVALID_ARGS;
557      }
558
559      return NLOPT_SUCCESS;
560 }
561
562 nlopt_result nlopt_minimize_econstrained(
563      nlopt_algorithm algorithm,
564      int n, nlopt_func f, void *f_data,
565      int m, nlopt_func fc, void *fc_data, ptrdiff_t fc_datum_size,
566      int p, nlopt_func h, void *h_data, ptrdiff_t h_datum_size,
567      const double *lb, const double *ub, /* bounds */
568      double *x, /* in: initial guess, out: minimizer */
569      double *minf, /* out: minimum */
570      double minf_max, double ftol_rel, double ftol_abs,
571      double xtol_rel, const double *xtol_abs,
572      double htol_rel, double htol_abs,
573      int maxeval, double maxtime)
574 {
575      nlopt_result ret;
576      if (xtol_abs)
577           ret = nlopt_minimize_(algorithm, n, f, f_data,
578                                 m, fc, fc_data, fc_datum_size, 
579                                 p, h, h_data, h_datum_size,
580                                 lb, ub,
581                                 x, minf, minf_max, ftol_rel, ftol_abs,
582                                 xtol_rel, xtol_abs, htol_rel, htol_abs,
583                                 maxeval, maxtime);
584      else {
585           int i;
586           double *xtol = (double *) malloc(sizeof(double) * n);
587           if (!xtol) return NLOPT_OUT_OF_MEMORY;
588           for (i = 0; i < n; ++i) xtol[i] = -1;
589           ret = nlopt_minimize_(algorithm, n, f, f_data, 
590                                 m, fc, fc_data, fc_datum_size,
591                                 p, h, h_data, h_datum_size,
592                                 lb, ub,
593                                 x, minf, minf_max, ftol_rel, ftol_abs,
594                                 xtol_rel, xtol, htol_rel, htol_abs,
595                                 maxeval, maxtime);
596           free(xtol);
597      }
598      return ret;
599 }
600
601 nlopt_result nlopt_minimize_constrained(
602      nlopt_algorithm algorithm,
603      int n, nlopt_func f, void *f_data,
604      int m, nlopt_func fc, void *fc_data, ptrdiff_t fc_datum_size,
605      const double *lb, const double *ub, /* bounds */
606      double *x, /* in: initial guess, out: minimizer */
607      double *minf, /* out: minimum */
608      double minf_max, double ftol_rel, double ftol_abs,
609      double xtol_rel, const double *xtol_abs,
610      int maxeval, double maxtime)
611 {
612      return nlopt_minimize_econstrained(
613           algorithm, n, f, f_data, 
614           m, fc, fc_data, fc_datum_size, 0, NULL, NULL, 0,
615           lb, ub, x, minf, minf_max, ftol_rel, ftol_abs,
616           xtol_rel, xtol_abs, ftol_rel, ftol_abs, maxeval, maxtime);
617 }
618
619 nlopt_result nlopt_minimize(
620      nlopt_algorithm algorithm,
621      int n, nlopt_func f, void *f_data,
622      const double *lb, const double *ub, /* bounds */
623      double *x, /* in: initial guess, out: minimizer */
624      double *minf, /* out: minimum */
625      double minf_max, double ftol_rel, double ftol_abs,
626      double xtol_rel, const double *xtol_abs,
627      int maxeval, double maxtime)
628 {
629      return nlopt_minimize_constrained(
630           algorithm, n, f, f_data, 0, NULL, NULL, 0,
631           lb, ub, x, minf, minf_max, ftol_rel, ftol_abs,
632           xtol_rel, xtol_abs, maxeval, maxtime);
633 }