chiark / gitweb /
add mconstraint to C API (not yet tested), thanks to Dmitrey of OpenOpt for the sugge...
[nlopt.git] / api / optimize.c
1 /* Copyright (c) 2007-2010 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 <float.h>
26
27 #include "nlopt-internal.h"
28
29 /*********************************************************************/
30
31 #ifndef HAVE_ISNAN
32 static int my_isnan(double x) { return x != x; }
33 #  define isnan my_isnan
34 #endif
35
36 /*********************************************************************/
37
38 #include "praxis.h"
39 #include "direct.h"
40
41 #ifdef WITH_CXX
42 #  include "stogo.h"
43 #endif
44
45 #include "cdirect.h"
46
47 #ifdef WITH_NOCEDAL
48 #  include "l-bfgs-b.h"
49 #endif
50
51 #include "luksan.h"
52
53 #include "crs.h"
54
55 #include "mlsl.h"
56 #include "mma.h"
57 #include "cobyla.h"
58 #include "newuoa.h"
59 #include "neldermead.h"
60 #include "auglag.h"
61 #include "bobyqa.h"
62 #include "isres.h"
63
64 /*********************************************************************/
65
66 static double f_bound(int n, const double *x, void *data_)
67 {
68      int i;
69      nlopt_opt data = (nlopt_opt) data_;
70      double f;
71
72      /* some methods do not support bound constraints, but support
73         discontinuous objectives so we can just return Inf for invalid x */
74      for (i = 0; i < n; ++i)
75           if (x[i] < data->lb[i] || x[i] > data->ub[i])
76                return HUGE_VAL;
77
78      f = data->f((unsigned) n, x, NULL, data->f_data);
79      return (isnan(f) || nlopt_isinf(f) ? HUGE_VAL : f);
80 }
81
82 static double f_noderiv(int n, const double *x, void *data_)
83 {
84      nlopt_opt data = (nlopt_opt) data_;
85      return data->f((unsigned) n, x, NULL, data->f_data);
86 }
87
88 static double f_direct(int n, const double *x, int *undefined, void *data_)
89 {
90      nlopt_opt data = (nlopt_opt) data_;
91      double f;
92      unsigned i, j;
93      f = data->f((unsigned) n, x, NULL, data->f_data);
94      *undefined = isnan(f) || nlopt_isinf(f);
95      for (i = 0; i < data->m && !*undefined; ++i) {
96           nlopt_eval_constraint(data->work, NULL, data->fc+i, (unsigned) n, x);
97           for (j = 0; j < data->fc[i].m; ++j)
98                if (data->work[j] > 0)
99                     *undefined = 1;
100      }
101      return f;
102 }
103
104 /*********************************************************************/
105
106 /* get min(dx) for algorithms requiring a scalar initial step size */
107 static nlopt_result initial_step(nlopt_opt opt, const double *x, double *step)
108 {
109      unsigned freedx = 0, i;
110
111      if (!opt->dx) {
112           freedx = 1;
113           if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
114                return NLOPT_OUT_OF_MEMORY;
115      }
116
117      *step = HUGE_VAL;
118      for (i = 0; i < opt->n; ++i)
119           if (*step > fabs(opt->dx[i]))
120                *step = fabs(opt->dx[i]);
121
122      if (freedx) { free(opt->dx); opt->dx = NULL; }
123      return NLOPT_SUCCESS;
124 }
125
126 /*********************************************************************/
127
128 /* return true if [lb,ub] is finite in every dimension (n dimensions) */
129 static int finite_domain(unsigned n, const double *lb, const double *ub)
130 {
131      unsigned i;
132      for (i = 0; i < n; ++i)
133           if (nlopt_isinf(ub[i] - lb[i])) return 0;
134      return 1;
135 }
136
137 /*********************************************************************/
138
139 #define POP(defaultpop) (opt->stochastic_population > 0 ?               \
140                          opt->stochastic_population :                   \
141                          (nlopt_stochastic_population > 0 ?             \
142                           nlopt_stochastic_population : (defaultpop)))
143
144 /* unlike nlopt_optimize() below, only handles minimization case */
145 static nlopt_result nlopt_optimize_(nlopt_opt opt, double *x, double *minf)
146 {
147      const double *lb, *ub;
148      nlopt_algorithm algorithm;
149      nlopt_func f; void *f_data;
150      unsigned n, i;
151      int ni;
152      nlopt_stopping stop;
153
154      if (!opt || !x || !minf || !opt->f
155          || opt->maximize) return NLOPT_INVALID_ARGS;
156
157      /* reset stopping flag */
158      nlopt_set_force_stop(opt, 0);
159      opt->force_stop_child = NULL;
160      
161      /* copy a few params to local vars for convenience */
162      n = opt->n;
163      ni = (int) n; /* most of the subroutines take "int" arg */
164      lb = opt->lb; ub = opt->ub;
165      algorithm = opt->algorithm;
166      f = opt->f; f_data = opt->f_data;
167
168      if (n == 0) { /* trivial case: no degrees of freedom */
169           *minf = opt->f(n, x, NULL, opt->f_data);
170           return NLOPT_SUCCESS;
171      }
172
173      *minf = HUGE_VAL;
174      
175      /* make sure rand generator is inited */
176      nlopt_srand_time_default(); /* default is non-deterministic */
177
178      /* check bound constraints */
179      for (i = 0; i < n; ++i)
180           if (lb[i] > ub[i] || x[i] < lb[i] || x[i] > ub[i])
181                return NLOPT_INVALID_ARGS;
182
183      stop.n = n;
184      stop.minf_max = opt->stopval;
185      stop.ftol_rel = opt->ftol_rel;
186      stop.ftol_abs = opt->ftol_abs;
187      stop.xtol_rel = opt->xtol_rel;
188      stop.xtol_abs = opt->xtol_abs;
189      stop.nevals = 0;
190      stop.maxeval = opt->maxeval;
191      stop.maxtime = opt->maxtime;
192      stop.start = nlopt_seconds();
193      stop.force_stop = &(opt->force_stop);
194
195      switch (algorithm) {
196          case NLOPT_GN_DIRECT:
197          case NLOPT_GN_DIRECT_L: 
198          case NLOPT_GN_DIRECT_L_RAND: 
199               if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
200               return cdirect(ni, f, f_data, 
201                              lb, ub, x, minf, &stop, 0.0, 
202                              (algorithm != NLOPT_GN_DIRECT)
203                              + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND 
204                                     ? 2 : (algorithm != NLOPT_GN_DIRECT))
205                              + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND 
206                                     ? 1 : (algorithm != NLOPT_GN_DIRECT)));
207               
208          case NLOPT_GN_DIRECT_NOSCAL:
209          case NLOPT_GN_DIRECT_L_NOSCAL: 
210          case NLOPT_GN_DIRECT_L_RAND_NOSCAL: 
211               if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
212               return cdirect_unscaled(ni, f, f_data, lb, ub, x, minf, 
213                                       &stop, 0.0, 
214                                       (algorithm != NLOPT_GN_DIRECT)
215                                       + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
216                                       + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
217               
218          case NLOPT_GN_ORIG_DIRECT:
219          case NLOPT_GN_ORIG_DIRECT_L: {
220               direct_return_code dret;
221               if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
222               opt->work = (double*) malloc(sizeof(double) *
223                                            nlopt_max_constraint_dim(opt->m,
224                                                                     opt->fc));
225               if (!opt->work) return NLOPT_OUT_OF_MEMORY;
226               dret = direct_optimize(f_direct, opt, ni, lb, ub, x, minf,
227                                      stop.maxeval, -1, 0.0, 0.0,
228                                      pow(stop.xtol_rel, (double) n), -1.0,
229                                      stop.minf_max, 0.0,
230                                      NULL, 
231                                      algorithm == NLOPT_GN_ORIG_DIRECT
232                                      ? DIRECT_ORIGINAL
233                                      : DIRECT_GABLONSKY);
234               free(opt->work); opt->work = NULL;
235               switch (dret) {
236                   case DIRECT_INVALID_BOUNDS:
237                   case DIRECT_MAXFEVAL_TOOBIG:
238                   case DIRECT_INVALID_ARGS:
239                        return NLOPT_INVALID_ARGS;
240                   case DIRECT_INIT_FAILED:
241                   case DIRECT_SAMPLEPOINTS_FAILED:
242                   case DIRECT_SAMPLE_FAILED:
243                        return NLOPT_FAILURE;
244                   case DIRECT_MAXFEVAL_EXCEEDED:
245                   case DIRECT_MAXITER_EXCEEDED:
246                        return NLOPT_MAXEVAL_REACHED;
247                   case DIRECT_GLOBAL_FOUND:
248                        return NLOPT_MINF_MAX_REACHED;
249                   case DIRECT_VOLTOL:
250                   case DIRECT_SIGMATOL:
251                        return NLOPT_XTOL_REACHED;
252                   case DIRECT_OUT_OF_MEMORY:
253                        return NLOPT_OUT_OF_MEMORY;
254               }
255               break;
256          }
257
258          case NLOPT_GD_STOGO:
259          case NLOPT_GD_STOGO_RAND:
260 #ifdef WITH_CXX
261               if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
262               if (!stogo_minimize(ni, f, f_data, x, minf, lb, ub, &stop,
263                                   algorithm == NLOPT_GD_STOGO
264                                   ? 0 : (int) POP(2*n)))
265                    return NLOPT_FAILURE;
266               break;
267 #else
268               return NLOPT_FAILURE;
269 #endif
270
271 #if 0
272               /* lacking a free/open-source license, we no longer use
273                  Rowan's code, and instead use by "sbplx" re-implementation */
274          case NLOPT_LN_SUBPLEX: {
275               int iret, freedx = 0;
276               if (!opt->dx) {
277                    freedx = 1;
278                    if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
279                         return NLOPT_OUT_OF_MEMORY;
280               }                
281               iret = nlopt_subplex(f_bound, minf, x, n, opt, &stop, opt->dx);
282               if (freedx) { free(opt->dx); opt->dx = NULL; }
283               switch (iret) {
284                   case -2: return NLOPT_INVALID_ARGS;
285                   case -20: return NLOPT_FORCED_STOP;
286                   case -10: return NLOPT_MAXTIME_REACHED;
287                   case -1: return NLOPT_MAXEVAL_REACHED;
288                   case 0: return NLOPT_XTOL_REACHED;
289                   case 1: return NLOPT_SUCCESS;
290                   case 2: return NLOPT_MINF_MAX_REACHED;
291                   case 20: return NLOPT_FTOL_REACHED;
292                   case -200: return NLOPT_OUT_OF_MEMORY;
293                   default: return NLOPT_FAILURE; /* unknown return code */
294               }
295               break;
296          }
297 #endif
298
299          case NLOPT_LN_PRAXIS: {
300               double step;
301               if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
302                    return NLOPT_OUT_OF_MEMORY;
303               return praxis_(0.0, DBL_EPSILON, 
304                              step, ni, x, f_bound, opt, &stop, minf);
305          }
306
307 #ifdef WITH_NOCEDAL
308          case NLOPT_LD_LBFGS_NOCEDAL: {
309               int iret, *nbd = (int *) malloc(sizeof(int) * n);
310               if (!nbd) return NLOPT_OUT_OF_MEMORY;
311               for (i = 0; i < n; ++i) {
312                    int linf = nlopt_isinf(lb[i]) && lb[i] < 0;
313                    int uinf = nlopt_isinf(ub[i]) && ub[i] > 0;
314                    nbd[i] = linf && uinf ? 0 : (uinf ? 1 : (linf ? 3 : 2));
315               }
316               iret = lbfgsb_minimize(ni, f, f_data, x, nbd, lb, ub,
317                                      ni < 5 ? ni : 5, 0.0, stop.ftol_rel, 
318                                      stop.xtol_abs[0] > 0 ? stop.xtol_abs[0]
319                                      : stop.xtol_rel,
320                                      stop.maxeval);
321               free(nbd);
322               if (iret <= 0) {
323                    switch (iret) {
324                        case -1: return NLOPT_INVALID_ARGS;
325                        case -2: default: return NLOPT_FAILURE;
326                    }
327               }
328               else {
329                    *minf = f(n, x, NULL, f_data);
330                    switch (iret) {
331                        case 5: return NLOPT_MAXEVAL_REACHED;
332                        case 2: return NLOPT_XTOL_REACHED;
333                        case 1: return NLOPT_FTOL_REACHED;
334                        default: return NLOPT_SUCCESS;
335                    }
336               }
337               break;
338          }
339 #endif
340
341          case NLOPT_LD_LBFGS: 
342               return luksan_plis(ni, f, f_data, lb, ub, x, minf, &stop);
343
344          case NLOPT_LD_VAR1: 
345          case NLOPT_LD_VAR2: 
346               return luksan_plip(ni, f, f_data, lb, ub, x, minf, &stop,
347                    algorithm == NLOPT_LD_VAR1 ? 1 : 2);
348
349          case NLOPT_LD_TNEWTON: 
350          case NLOPT_LD_TNEWTON_RESTART: 
351          case NLOPT_LD_TNEWTON_PRECOND: 
352          case NLOPT_LD_TNEWTON_PRECOND_RESTART: 
353               return luksan_pnet(ni, f, f_data, lb, ub, x, minf, &stop,
354                                  1 + (algorithm - NLOPT_LD_TNEWTON) % 2,
355                                  1 + (algorithm - NLOPT_LD_TNEWTON) / 2);
356
357          case NLOPT_GN_CRS2_LM:
358               if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
359               return crs_minimize(ni, f, f_data, lb, ub, x, minf, &stop, 
360                                   (int) POP(0), 0);
361
362          case NLOPT_G_MLSL:
363          case NLOPT_G_MLSL_LDS:
364          case NLOPT_GN_MLSL:
365          case NLOPT_GD_MLSL:
366          case NLOPT_GN_MLSL_LDS:
367          case NLOPT_GD_MLSL_LDS: {
368               if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
369               nlopt_opt local_opt = opt->local_opt;
370               nlopt_result ret;
371               if (!local_opt && (algorithm == NLOPT_G_MLSL 
372                                  || algorithm == NLOPT_G_MLSL_LDS))
373                    return NLOPT_INVALID_ARGS;
374               if (!local_opt) { /* default */
375                    nlopt_algorithm local_alg = (algorithm == NLOPT_GN_MLSL ||
376                                                 algorithm == NLOPT_GN_MLSL_LDS)
377                         ? nlopt_local_search_alg_nonderiv
378                         : nlopt_local_search_alg_deriv;
379                    /* don't call MLSL recursively! */
380                    if (local_alg >= NLOPT_GN_MLSL
381                        && local_alg <= NLOPT_GD_MLSL_LDS)
382                         local_alg = (algorithm == NLOPT_GN_MLSL ||
383                                      algorithm == NLOPT_GN_MLSL_LDS)
384                              ? NLOPT_LN_COBYLA : NLOPT_LD_MMA;
385                    local_opt = nlopt_create(local_alg, n);
386                    if (!local_opt) return NLOPT_FAILURE;
387                    nlopt_set_ftol_rel(local_opt, opt->ftol_rel);
388                    nlopt_set_ftol_abs(local_opt, opt->ftol_abs);
389                    nlopt_set_xtol_rel(local_opt, opt->xtol_rel);
390                    nlopt_set_xtol_abs(local_opt, opt->xtol_abs);
391                    nlopt_set_maxeval(local_opt, nlopt_local_search_maxeval);
392               }
393               if (opt->dx) nlopt_set_initial_step(local_opt, opt->dx);
394               for (i = 0; i < n && stop.xtol_abs[i] > 0; ++i) ;
395               if (local_opt->ftol_rel <= 0 && local_opt->ftol_abs <= 0 &&
396                   local_opt->xtol_rel <= 0 && i < n) {
397                    /* it is not sensible to call MLSL without *some*
398                       nonzero tolerance for the local search */
399                    nlopt_set_ftol_rel(local_opt, 1e-15);
400                    nlopt_set_xtol_rel(local_opt, 1e-7);
401               }
402               opt->force_stop_child = local_opt;
403               ret = mlsl_minimize(ni, f, f_data, lb, ub, x, minf, &stop,
404                                   local_opt, (int) POP(0),
405                                   algorithm >= NLOPT_GN_MLSL_LDS &&
406                                   algorithm != NLOPT_G_MLSL);
407               opt->force_stop_child = NULL;
408               if (!opt->local_opt) nlopt_destroy(local_opt);
409               return ret;
410          }
411
412          case NLOPT_LD_MMA: {
413               nlopt_opt dual_opt;
414               nlopt_result ret;
415 #define LO(param, def) (opt->local_opt ? opt->local_opt->param : (def))
416               dual_opt = nlopt_create(LO(algorithm,
417                                          nlopt_local_search_alg_deriv),
418                                       opt->m);
419               if (!dual_opt) return NLOPT_FAILURE;
420               nlopt_set_ftol_rel(dual_opt, LO(ftol_rel, 1e-12));
421               nlopt_set_ftol_abs(dual_opt, LO(ftol_abs, 0.0));
422               nlopt_set_maxeval(dual_opt, LO(maxeval, 100000));
423 #undef LO
424
425               ret = mma_minimize(ni, f, f_data, (int) (opt->m), opt->fc,
426                                  lb, ub, x, minf, &stop, dual_opt);
427               nlopt_destroy(dual_opt);
428               return ret;
429          }
430
431          case NLOPT_LN_COBYLA: {
432               double step;
433               if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
434                    return NLOPT_OUT_OF_MEMORY;
435               return cobyla_minimize(ni, f, f_data, 
436                                      opt->m, opt->fc,
437                                      opt->p, opt->h,
438                                      lb, ub, x, minf, &stop,
439                                      step);
440          }
441                                      
442          case NLOPT_LN_NEWUOA: {
443               double step;
444               if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
445                    return NLOPT_OUT_OF_MEMORY;
446               return newuoa(ni, 2*n+1, x, 0, 0, step,
447                             &stop, minf, f_noderiv, opt);
448          }
449                                      
450          case NLOPT_LN_NEWUOA_BOUND: {
451               double step;
452               if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
453                    return NLOPT_OUT_OF_MEMORY;
454               return newuoa(ni, 2*n+1, x, lb, ub, step,
455                             &stop, minf, f_noderiv, opt);
456          }
457
458          case NLOPT_LN_BOBYQA: {
459               double step;
460               if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
461                    return NLOPT_OUT_OF_MEMORY;
462               return bobyqa(ni, 2*n+1, x, lb, ub, step,
463                             &stop, minf, f_noderiv, opt);
464          }
465
466          case NLOPT_LN_NELDERMEAD: 
467          case NLOPT_LN_SBPLX: 
468          {
469               nlopt_result ret;
470               int freedx = 0;
471               if (!opt->dx) {
472                    freedx = 1;
473                    if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
474                         return NLOPT_OUT_OF_MEMORY;
475               }
476               if (algorithm == NLOPT_LN_NELDERMEAD)
477                    ret= nldrmd_minimize(ni,f,f_data,lb,ub,x,minf,opt->dx,&stop);
478               else
479                    ret= sbplx_minimize(ni,f,f_data,lb,ub,x,minf,opt->dx,&stop);
480               if (freedx) { free(opt->dx); opt->dx = NULL; }
481               return ret;
482          }
483
484          case NLOPT_AUGLAG:
485          case NLOPT_AUGLAG_EQ:
486          case NLOPT_LN_AUGLAG:
487          case NLOPT_LN_AUGLAG_EQ:
488          case NLOPT_LD_AUGLAG:
489          case NLOPT_LD_AUGLAG_EQ: {
490               nlopt_opt local_opt = opt->local_opt;
491               nlopt_result ret;
492               if ((algorithm == NLOPT_AUGLAG || algorithm == NLOPT_AUGLAG_EQ)
493                   && !local_opt)
494                    return NLOPT_INVALID_ARGS;
495               if (!local_opt) { /* default */
496                    local_opt = nlopt_create(
497                         algorithm == NLOPT_LN_AUGLAG || 
498                         algorithm == NLOPT_LN_AUGLAG_EQ
499                         ? nlopt_local_search_alg_nonderiv
500                         : nlopt_local_search_alg_deriv, n);
501                    if (!local_opt) return NLOPT_FAILURE;
502                    nlopt_set_ftol_rel(local_opt, opt->ftol_rel);
503                    nlopt_set_ftol_abs(local_opt, opt->ftol_abs);
504                    nlopt_set_xtol_rel(local_opt, opt->xtol_rel);
505                    nlopt_set_xtol_abs(local_opt, opt->xtol_abs);
506                    nlopt_set_maxeval(local_opt, nlopt_local_search_maxeval);
507               }
508               if (opt->dx) nlopt_set_initial_step(local_opt, opt->dx);
509               opt->force_stop_child = local_opt;
510               ret = auglag_minimize(ni, f, f_data, 
511                                     opt->m, opt->fc, 
512                                     opt->p, opt->h,
513                                     lb, ub, x, minf, &stop,
514                                     local_opt,
515                                     algorithm == NLOPT_AUGLAG_EQ
516                                     || algorithm == NLOPT_LN_AUGLAG_EQ
517                                     || algorithm == NLOPT_LD_AUGLAG_EQ);
518               opt->force_stop_child = NULL;
519               if (!opt->local_opt) nlopt_destroy(local_opt);
520               return ret;
521          }
522
523          case NLOPT_GN_ISRES:
524               if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
525               return isres_minimize(ni, f, f_data, 
526                                     (int) (opt->m), opt->fc,
527                                     (int) (opt->p), opt->h,
528                                     lb, ub, x, minf, &stop,
529                                     (int) POP(0));
530
531          default:
532               return NLOPT_INVALID_ARGS;
533      }
534
535      return NLOPT_SUCCESS; /* never reached */
536 }
537
538 /*********************************************************************/
539
540 typedef struct {
541      nlopt_func f;
542      void *f_data;
543 } f_max_data;
544
545 /* wrapper for maximizing: just flip the sign of f and grad */
546 static double f_max(unsigned n, const double *x, double *grad, void *data)
547 {
548      f_max_data *d = (f_max_data *) data;
549      double val = d->f(n, x, grad, d->f_data);
550      if (grad) {
551           unsigned i;
552           for (i = 0; i < n; ++i)
553                grad[i] = -grad[i];
554      }
555      return -val;
556 }
557
558 nlopt_result 
559 NLOPT_STDCALL nlopt_optimize(nlopt_opt opt, double *x, double *opt_f)
560 {
561      nlopt_func f; void *f_data;
562      f_max_data fmd;
563      int maximize;
564      nlopt_result ret;
565
566      if (!opt || !opt_f || !opt->f) return NLOPT_INVALID_ARGS;
567      f = opt->f; f_data = opt->f_data;
568
569      /* for maximizing, just minimize the f_max wrapper, which 
570         flips the sign of everything */
571      if ((maximize = opt->maximize)) {
572           fmd.f = f; fmd.f_data = f_data;
573           opt->f = f_max; opt->f_data = &fmd;
574           opt->stopval = -opt->stopval;
575           opt->maximize = 0;
576      }
577
578      ret = nlopt_optimize_(opt, x, opt_f);
579
580      if (maximize) { /* restore original signs */
581           opt->maximize = maximize;
582           opt->stopval = -opt->stopval;
583           opt->f = f; opt->f_data = f_data;
584           *opt_f = -*opt_f;
585      }
586
587      return ret;
588 }
589
590 /*********************************************************************/
591
592 nlopt_result nlopt_optimize_limited(nlopt_opt opt, double *x, double *minf,
593                                     int maxeval, double maxtime)
594 {
595      int save_maxeval;
596      double save_maxtime;
597      nlopt_result ret;
598
599      if (!opt) return NLOPT_INVALID_ARGS;
600
601      save_maxeval = nlopt_get_maxeval(opt);
602      save_maxtime = nlopt_get_maxtime(opt);
603      
604      /* override opt limits if maxeval and/or maxtime are more stringent */
605      if (save_maxeval <= 0 || (maxeval > 0 && maxeval < save_maxeval))
606           nlopt_set_maxeval(opt, maxeval);
607      if (save_maxtime <= 0 || (maxtime > 0 && maxtime < save_maxtime))
608           nlopt_set_maxtime(opt, maxtime);
609
610      ret = nlopt_optimize(opt, x, minf);
611
612      nlopt_set_maxeval(opt, save_maxeval);
613      nlopt_set_maxtime(opt, save_maxtime);
614
615      return ret;
616 }
617
618 /*********************************************************************/