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