1 /* Copyright (c) 2007-2011 Massachusetts Institute of Technology
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:
11 * The above copyright notice and this permission notice shall be
12 * included in all copies or substantial portions of the Software.
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.
27 #include "nlopt-internal.h"
29 /*********************************************************************/
32 static int my_isnan(double x) { return x != x; }
33 # define isnan my_isnan
36 /*********************************************************************/
48 # include "l-bfgs-b.h"
59 #include "neldermead.h"
65 /*********************************************************************/
67 static double f_bound(int n, const double *x, void *data_)
70 nlopt_opt data = (nlopt_opt) data_;
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])
79 f = data->f((unsigned) n, x, NULL, data->f_data);
80 return (isnan(f) || nlopt_isinf(f) ? HUGE_VAL : f);
83 static double f_noderiv(int n, const double *x, void *data_)
85 nlopt_opt data = (nlopt_opt) data_;
86 return data->f((unsigned) n, x, NULL, data->f_data);
89 static double f_direct(int n, const double *x, int *undefined, void *data_)
91 nlopt_opt data = (nlopt_opt) data_;
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)
107 /*********************************************************************/
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)
112 unsigned freedx = 0, i;
116 if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
117 return NLOPT_OUT_OF_MEMORY;
121 for (i = 0; i < opt->n; ++i)
122 if (*step > fabs(opt->dx[i]))
123 *step = fabs(opt->dx[i]);
125 if (freedx) { free(opt->dx); opt->dx = NULL; }
126 return NLOPT_SUCCESS;
129 /*********************************************************************/
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)
135 for (i = 0; i < n; ++i)
136 if (nlopt_isinf(ub[i] - lb[i])) return 0;
140 /*********************************************************************/
141 /* wrapper functions, only for derivative-free methods, that
142 eliminate dimensions with lb == ub. (The gradient-based methods
143 should handle this case directly, since they operate on much
144 larger vectors where I am loathe to make copies unnecessarily.) */
150 unsigned n; /* true dimension */
151 double *x; /* scratch vector of length n */
152 double *grad; /* optional scratch vector of length n */
153 const double *lb, *ub; /* bounds, of length n */
156 static void *elimdim_makedata(nlopt_func f, nlopt_mfunc mf, void *f_data,
157 unsigned n, double *x, const double *lb,
158 const double *ub, double *grad)
160 elimdim_data *d = (elimdim_data *) malloc(sizeof(elimdim_data));
162 d->f = f; d->mf = mf; d->f_data = f_data; d->n = n; d->x = x;
163 d->lb = lb; d->ub = ub;
168 static double elimdim_func(unsigned n0, const double *x0, double *grad, void *d_)
170 elimdim_data *d = (elimdim_data *) d_;
172 const double *lb = d->lb, *ub = d->ub;
174 unsigned n = d->n, i, j;
176 (void) n0; /* unused */
177 for (i = j = 0; i < n; ++i) {
180 else /* assert: j < n0 */
183 val = d->f(n, x, grad ? d->grad : NULL, d->f_data);
185 /* assert: d->grad != NULL */
186 for (i = j = 0; i < n; ++i)
188 grad[j++] = d->grad[i];
194 static void elimdim_mfunc(unsigned m, double *result,
195 unsigned n0, const double *x0, double *grad, void *d_)
197 elimdim_data *d = (elimdim_data *) d_;
199 const double *lb = d->lb, *ub = d->ub;
200 unsigned n = d->n, i, j;
202 (void) n0; /* unused */
203 (void) grad; /* assert: grad == NULL */
204 for (i = j = 0; i < n; ++i) {
207 else /* assert: j < n0 */
210 d->mf(m, result, n, x, NULL, d->f_data);
213 /* compute the eliminated dimension: number of dims with lb[i] != ub[i] */
214 static unsigned elimdim_dimension(unsigned n, const double *lb, const double *ub)
217 for (i = 0; i < n; ++i) n0 += lb[i] != ub[i] ? 1U : 0;
221 /* modify v to "shrunk" version, with dimensions for lb[i] == ub[i] elim'ed */
222 static void elimdim_shrink(unsigned n, double *v,
223 const double *lb, const double *ub)
227 for (i = j = 0; i < n; ++i)
232 /* inverse of elimdim_shrink */
233 static void elimdim_expand(unsigned n, double *v,
234 const double *lb, const double *ub)
238 j = elimdim_dimension(n, lb, ub) - 1;
239 for (i = n - 1; i > 0; --i) {
250 /* given opt, create a new opt with equal-constraint dimensions eliminated */
251 static nlopt_opt elimdim_create(nlopt_opt opt)
253 nlopt_opt opt0 = nlopt_copy(opt);
254 double *x, *grad = NULL;
257 if (!opt0) return NULL;
258 x = (double *) malloc(sizeof(double) * opt->n);
259 if (opt->n && !x) { nlopt_destroy(opt0); return NULL; }
261 if (opt->algorithm == NLOPT_GD_STOGO
262 || opt->algorithm == NLOPT_GD_STOGO_RAND) {
263 grad = (double *) malloc(sizeof(double) * opt->n);
264 if (opt->n && !grad) goto bad;
267 opt0->n = elimdim_dimension(opt->n, opt->lb, opt->ub);
268 elimdim_shrink(opt->n, opt0->lb, opt->lb, opt->ub);
269 elimdim_shrink(opt->n, opt0->ub, opt->lb, opt->ub);
270 elimdim_shrink(opt->n, opt0->xtol_abs, opt->lb, opt->ub);
271 elimdim_shrink(opt->n, opt0->dx, opt->lb, opt->ub);
273 opt0->munge_on_destroy = opt0->munge_on_copy = NULL;
275 opt0->f = elimdim_func;
276 opt0->f_data = elimdim_makedata(opt->f, NULL, opt->f_data,
277 opt->n, x, opt->lb, opt->ub, grad);
278 if (!opt0->f_data) goto bad;
280 for (i = 0; i < opt->m; ++i) {
281 opt0->fc[i].f = elimdim_func;
282 opt0->fc[i].mf = elimdim_mfunc;
283 opt0->fc[i].f_data = elimdim_makedata(opt->fc[i].f, opt->fc[i].mf,
285 opt->n, x, opt->lb, opt->ub,
287 if (!opt0->fc[i].f_data) goto bad;
290 for (i = 0; i < opt->p; ++i) {
291 opt0->h[i].f = elimdim_func;
292 opt0->h[i].mf = elimdim_mfunc;
293 opt0->h[i].f_data = elimdim_makedata(opt->h[i].f, opt->h[i].mf,
295 opt->n, x, opt->lb, opt->ub,
297 if (!opt0->h[i].f_data) goto bad;
308 /* like nlopt_destroy, but also frees elimdim_data */
309 static void elimdim_destroy(nlopt_opt opt)
314 free(((elimdim_data*) opt->f_data)->x);
315 free(((elimdim_data*) opt->f_data)->grad);
316 free(opt->f_data); opt->f_data = NULL;
318 for (i = 0; i < opt->m; ++i) {
319 free(opt->fc[i].f_data);
320 opt->fc[i].f_data = NULL;
322 for (i = 0; i < opt->p; ++i) {
323 free(opt->h[i].f_data);
324 opt->h[i].f_data = NULL;
330 /* return whether to use elimdim wrapping. */
331 static int elimdim_wrapcheck(nlopt_opt opt)
334 if (elimdim_dimension(opt->n, opt->lb, opt->ub) == opt->n) return 0;
335 switch (opt->algorithm) {
336 case NLOPT_GN_DIRECT:
337 case NLOPT_GN_DIRECT_L:
338 case NLOPT_GN_DIRECT_L_RAND:
339 case NLOPT_GN_DIRECT_NOSCAL:
340 case NLOPT_GN_DIRECT_L_NOSCAL:
341 case NLOPT_GN_DIRECT_L_RAND_NOSCAL:
342 case NLOPT_GN_ORIG_DIRECT:
343 case NLOPT_GN_ORIG_DIRECT_L:
344 case NLOPT_LN_PRAXIS:
345 case NLOPT_LN_COBYLA:
346 case NLOPT_LN_NEWUOA:
347 case NLOPT_LN_NEWUOA_BOUND:
348 case NLOPT_LN_BOBYQA:
349 case NLOPT_LN_NELDERMEAD:
353 case NLOPT_GD_STOGO_RAND:
360 /*********************************************************************/
362 #define POP(defaultpop) (opt->stochastic_population > 0 ? \
363 opt->stochastic_population : \
364 (nlopt_stochastic_population > 0 ? \
365 nlopt_stochastic_population : (defaultpop)))
367 /* unlike nlopt_optimize() below, only handles minimization case */
368 static nlopt_result nlopt_optimize_(nlopt_opt opt, double *x, double *minf)
370 const double *lb, *ub;
371 nlopt_algorithm algorithm;
372 nlopt_func f; void *f_data;
377 if (!opt || !x || !minf || !opt->f
378 || opt->maximize) return NLOPT_INVALID_ARGS;
380 /* reset stopping flag */
381 nlopt_set_force_stop(opt, 0);
382 opt->force_stop_child = NULL;
384 /* copy a few params to local vars for convenience */
386 ni = (int) n; /* most of the subroutines take "int" arg */
387 lb = opt->lb; ub = opt->ub;
388 algorithm = opt->algorithm;
389 f = opt->f; f_data = opt->f_data;
391 if (n == 0) { /* trivial case: no degrees of freedom */
392 *minf = opt->f(n, x, NULL, opt->f_data);
393 return NLOPT_SUCCESS;
398 /* make sure rand generator is inited */
399 nlopt_srand_time_default(); /* default is non-deterministic */
401 /* check bound constraints */
402 for (i = 0; i < n; ++i)
403 if (lb[i] > ub[i] || x[i] < lb[i] || x[i] > ub[i])
404 return NLOPT_INVALID_ARGS;
407 stop.minf_max = opt->stopval;
408 stop.ftol_rel = opt->ftol_rel;
409 stop.ftol_abs = opt->ftol_abs;
410 stop.xtol_rel = opt->xtol_rel;
411 stop.xtol_abs = opt->xtol_abs;
413 stop.maxeval = opt->maxeval;
414 stop.maxtime = opt->maxtime;
415 stop.start = nlopt_seconds();
416 stop.force_stop = &(opt->force_stop);
419 case NLOPT_GN_DIRECT:
420 case NLOPT_GN_DIRECT_L:
421 case NLOPT_GN_DIRECT_L_RAND:
422 if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
423 return cdirect(ni, f, f_data,
424 lb, ub, x, minf, &stop, 0.0,
425 (algorithm != NLOPT_GN_DIRECT)
426 + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND
427 ? 2 : (algorithm != NLOPT_GN_DIRECT))
428 + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND
429 ? 1 : (algorithm != NLOPT_GN_DIRECT)));
431 case NLOPT_GN_DIRECT_NOSCAL:
432 case NLOPT_GN_DIRECT_L_NOSCAL:
433 case NLOPT_GN_DIRECT_L_RAND_NOSCAL:
434 if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
435 return cdirect_unscaled(ni, f, f_data, lb, ub, x, minf,
437 (algorithm != NLOPT_GN_DIRECT)
438 + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
439 + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
441 case NLOPT_GN_ORIG_DIRECT:
442 case NLOPT_GN_ORIG_DIRECT_L: {
443 direct_return_code dret;
444 if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
445 opt->work = (double*) malloc(sizeof(double) *
446 nlopt_max_constraint_dim(opt->m,
448 if (!opt->work) return NLOPT_OUT_OF_MEMORY;
449 dret = direct_optimize(f_direct, opt, ni, lb, ub, x, minf,
451 stop.start, stop.maxtime,
453 pow(stop.xtol_rel, (double) n), -1.0,
457 algorithm == NLOPT_GN_ORIG_DIRECT
460 free(opt->work); opt->work = NULL;
462 case DIRECT_INVALID_BOUNDS:
463 case DIRECT_MAXFEVAL_TOOBIG:
464 case DIRECT_INVALID_ARGS:
465 return NLOPT_INVALID_ARGS;
466 case DIRECT_INIT_FAILED:
467 case DIRECT_SAMPLEPOINTS_FAILED:
468 case DIRECT_SAMPLE_FAILED:
469 return NLOPT_FAILURE;
470 case DIRECT_MAXFEVAL_EXCEEDED:
471 case DIRECT_MAXITER_EXCEEDED:
472 return NLOPT_MAXEVAL_REACHED;
473 case DIRECT_MAXTIME_EXCEEDED:
474 return NLOPT_MAXTIME_REACHED;
475 case DIRECT_GLOBAL_FOUND:
476 return NLOPT_MINF_MAX_REACHED;
478 case DIRECT_SIGMATOL:
479 return NLOPT_XTOL_REACHED;
480 case DIRECT_OUT_OF_MEMORY:
481 return NLOPT_OUT_OF_MEMORY;
482 case DIRECT_FORCED_STOP:
483 return NLOPT_FORCED_STOP;
489 case NLOPT_GD_STOGO_RAND:
491 if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
492 if (!stogo_minimize(ni, f, f_data, x, minf, lb, ub, &stop,
493 algorithm == NLOPT_GD_STOGO
494 ? 0 : (int) POP(2*n)))
495 return NLOPT_FAILURE;
498 return NLOPT_INVALID_ARGS;
502 /* lacking a free/open-source license, we no longer use
503 Rowan's code, and instead use by "sbplx" re-implementation */
504 case NLOPT_LN_SUBPLEX: {
505 int iret, freedx = 0;
508 if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
509 return NLOPT_OUT_OF_MEMORY;
511 iret = nlopt_subplex(f_bound, minf, x, n, opt, &stop, opt->dx);
512 if (freedx) { free(opt->dx); opt->dx = NULL; }
514 case -2: return NLOPT_INVALID_ARGS;
515 case -20: return NLOPT_FORCED_STOP;
516 case -10: return NLOPT_MAXTIME_REACHED;
517 case -1: return NLOPT_MAXEVAL_REACHED;
518 case 0: return NLOPT_XTOL_REACHED;
519 case 1: return NLOPT_SUCCESS;
520 case 2: return NLOPT_MINF_MAX_REACHED;
521 case 20: return NLOPT_FTOL_REACHED;
522 case -200: return NLOPT_OUT_OF_MEMORY;
523 default: return NLOPT_FAILURE; /* unknown return code */
529 case NLOPT_LN_PRAXIS: {
531 if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
532 return NLOPT_OUT_OF_MEMORY;
533 return praxis_(0.0, DBL_EPSILON,
534 step, ni, x, f_bound, opt, &stop, minf);
538 case NLOPT_LD_LBFGS_NOCEDAL: {
539 int iret, *nbd = (int *) malloc(sizeof(int) * n);
540 if (!nbd) return NLOPT_OUT_OF_MEMORY;
541 for (i = 0; i < n; ++i) {
542 int linf = nlopt_isinf(lb[i]) && lb[i] < 0;
543 int uinf = nlopt_isinf(ub[i]) && ub[i] > 0;
544 nbd[i] = linf && uinf ? 0 : (uinf ? 1 : (linf ? 3 : 2));
546 iret = lbfgsb_minimize(ni, f, f_data, x, nbd, lb, ub,
547 ni < 5 ? ni : 5, 0.0, stop.ftol_rel,
548 stop.xtol_abs[0] > 0 ? stop.xtol_abs[0]
554 case -1: return NLOPT_INVALID_ARGS;
555 case -2: default: return NLOPT_FAILURE;
559 *minf = f(n, x, NULL, f_data);
561 case 5: return NLOPT_MAXEVAL_REACHED;
562 case 2: return NLOPT_XTOL_REACHED;
563 case 1: return NLOPT_FTOL_REACHED;
564 default: return NLOPT_SUCCESS;
572 return luksan_plis(ni, f, f_data, lb, ub, x, minf,
573 &stop, opt->vector_storage);
577 return luksan_plip(ni, f, f_data, lb, ub, x, minf,
578 &stop, opt->vector_storage,
579 algorithm == NLOPT_LD_VAR1 ? 1 : 2);
581 case NLOPT_LD_TNEWTON:
582 case NLOPT_LD_TNEWTON_RESTART:
583 case NLOPT_LD_TNEWTON_PRECOND:
584 case NLOPT_LD_TNEWTON_PRECOND_RESTART:
585 return luksan_pnet(ni, f, f_data, lb, ub, x, minf,
586 &stop, opt->vector_storage,
587 1 + (algorithm - NLOPT_LD_TNEWTON) % 2,
588 1 + (algorithm - NLOPT_LD_TNEWTON) / 2);
590 case NLOPT_GN_CRS2_LM:
591 if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
592 return crs_minimize(ni, f, f_data, lb, ub, x, minf, &stop,
596 case NLOPT_G_MLSL_LDS:
599 case NLOPT_GN_MLSL_LDS:
600 case NLOPT_GD_MLSL_LDS: {
601 nlopt_opt local_opt = opt->local_opt;
603 if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
604 if (!local_opt && (algorithm == NLOPT_G_MLSL
605 || algorithm == NLOPT_G_MLSL_LDS))
606 return NLOPT_INVALID_ARGS;
607 if (!local_opt) { /* default */
608 nlopt_algorithm local_alg = (algorithm == NLOPT_GN_MLSL ||
609 algorithm == NLOPT_GN_MLSL_LDS)
610 ? nlopt_local_search_alg_nonderiv
611 : nlopt_local_search_alg_deriv;
612 /* don't call MLSL recursively! */
613 if (local_alg >= NLOPT_GN_MLSL
614 && local_alg <= NLOPT_GD_MLSL_LDS)
615 local_alg = (algorithm == NLOPT_GN_MLSL ||
616 algorithm == NLOPT_GN_MLSL_LDS)
617 ? NLOPT_LN_COBYLA : NLOPT_LD_MMA;
618 local_opt = nlopt_create(local_alg, n);
619 if (!local_opt) return NLOPT_FAILURE;
620 nlopt_set_ftol_rel(local_opt, opt->ftol_rel);
621 nlopt_set_ftol_abs(local_opt, opt->ftol_abs);
622 nlopt_set_xtol_rel(local_opt, opt->xtol_rel);
623 nlopt_set_xtol_abs(local_opt, opt->xtol_abs);
624 nlopt_set_maxeval(local_opt, nlopt_local_search_maxeval);
626 if (opt->dx) nlopt_set_initial_step(local_opt, opt->dx);
627 for (i = 0; i < n && stop.xtol_abs[i] > 0; ++i) ;
628 if (local_opt->ftol_rel <= 0 && local_opt->ftol_abs <= 0 &&
629 local_opt->xtol_rel <= 0 && i < n) {
630 /* it is not sensible to call MLSL without *some*
631 nonzero tolerance for the local search */
632 nlopt_set_ftol_rel(local_opt, 1e-15);
633 nlopt_set_xtol_rel(local_opt, 1e-7);
635 opt->force_stop_child = local_opt;
636 ret = mlsl_minimize(ni, f, f_data, lb, ub, x, minf, &stop,
637 local_opt, (int) POP(0),
638 algorithm >= NLOPT_GN_MLSL_LDS &&
639 algorithm != NLOPT_G_MLSL);
640 opt->force_stop_child = NULL;
641 if (!opt->local_opt) nlopt_destroy(local_opt);
645 case NLOPT_LD_MMA: case NLOPT_LD_CCSAQ: {
648 #define LO(param, def) (opt->local_opt ? opt->local_opt->param : (def))
649 dual_opt = nlopt_create(LO(algorithm,
650 nlopt_local_search_alg_deriv),
651 nlopt_count_constraints(opt->m,
653 if (!dual_opt) return NLOPT_FAILURE;
654 nlopt_set_ftol_rel(dual_opt, LO(ftol_rel, 1e-14));
655 nlopt_set_ftol_abs(dual_opt, LO(ftol_abs, 0.0));
656 nlopt_set_maxeval(dual_opt, LO(maxeval, 100000));
659 if (algorithm == NLOPT_LD_MMA)
660 ret = mma_minimize(n, f, f_data, opt->m, opt->fc,
661 lb, ub, x, minf, &stop, dual_opt);
663 ret = ccsa_quadratic_minimize(
664 n, f, f_data, opt->m, opt->fc, opt->pre,
665 lb, ub, x, minf, &stop, dual_opt);
666 nlopt_destroy(dual_opt);
670 case NLOPT_LN_COBYLA: {
675 if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
676 return NLOPT_OUT_OF_MEMORY;
678 return cobyla_minimize(n, f, f_data,
681 lb, ub, x, minf, &stop,
683 if (freedx) { free(opt->dx); opt->dx = NULL; }
687 case NLOPT_LN_NEWUOA: {
689 if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
690 return NLOPT_OUT_OF_MEMORY;
691 return newuoa(ni, 2*n+1, x, 0, 0, step,
692 &stop, minf, f_noderiv, opt);
695 case NLOPT_LN_NEWUOA_BOUND: {
697 if (initial_step(opt, x, &step) != NLOPT_SUCCESS)
698 return NLOPT_OUT_OF_MEMORY;
699 return newuoa(ni, 2*n+1, x, lb, ub, step,
700 &stop, minf, f_noderiv, opt);
703 case NLOPT_LN_BOBYQA: {
708 if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
709 return NLOPT_OUT_OF_MEMORY;
711 ret = bobyqa(ni, 2*n+1, x, lb, ub, opt->dx,
712 &stop, minf, opt->f, opt->f_data);
713 if (freedx) { free(opt->dx); opt->dx = NULL; }
717 case NLOPT_LN_NELDERMEAD:
724 if (nlopt_set_default_initial_step(opt, x) != NLOPT_SUCCESS)
725 return NLOPT_OUT_OF_MEMORY;
727 if (algorithm == NLOPT_LN_NELDERMEAD)
728 ret= nldrmd_minimize(ni,f,f_data,lb,ub,x,minf,opt->dx,&stop);
730 ret= sbplx_minimize(ni,f,f_data,lb,ub,x,minf,opt->dx,&stop);
731 if (freedx) { free(opt->dx); opt->dx = NULL; }
736 case NLOPT_AUGLAG_EQ:
737 case NLOPT_LN_AUGLAG:
738 case NLOPT_LN_AUGLAG_EQ:
739 case NLOPT_LD_AUGLAG:
740 case NLOPT_LD_AUGLAG_EQ: {
741 nlopt_opt local_opt = opt->local_opt;
743 if ((algorithm == NLOPT_AUGLAG || algorithm == NLOPT_AUGLAG_EQ)
745 return NLOPT_INVALID_ARGS;
746 if (!local_opt) { /* default */
747 local_opt = nlopt_create(
748 algorithm == NLOPT_LN_AUGLAG ||
749 algorithm == NLOPT_LN_AUGLAG_EQ
750 ? nlopt_local_search_alg_nonderiv
751 : nlopt_local_search_alg_deriv, n);
752 if (!local_opt) return NLOPT_FAILURE;
753 nlopt_set_ftol_rel(local_opt, opt->ftol_rel);
754 nlopt_set_ftol_abs(local_opt, opt->ftol_abs);
755 nlopt_set_xtol_rel(local_opt, opt->xtol_rel);
756 nlopt_set_xtol_abs(local_opt, opt->xtol_abs);
757 nlopt_set_maxeval(local_opt, nlopt_local_search_maxeval);
759 if (opt->dx) nlopt_set_initial_step(local_opt, opt->dx);
760 opt->force_stop_child = local_opt;
761 ret = auglag_minimize(ni, f, f_data,
764 lb, ub, x, minf, &stop,
766 algorithm == NLOPT_AUGLAG_EQ
767 || algorithm == NLOPT_LN_AUGLAG_EQ
768 || algorithm == NLOPT_LD_AUGLAG_EQ);
769 opt->force_stop_child = NULL;
770 if (!opt->local_opt) nlopt_destroy(local_opt);
775 if (!finite_domain(n, lb, ub)) return NLOPT_INVALID_ARGS;
776 return isres_minimize(ni, f, f_data,
777 (int) (opt->m), opt->fc,
778 (int) (opt->p), opt->h,
779 lb, ub, x, minf, &stop,
783 return nlopt_slsqp(n, f, f_data,
786 lb, ub, x, minf, &stop);
789 return NLOPT_INVALID_ARGS;
792 return NLOPT_SUCCESS; /* never reached */
795 /*********************************************************************/
803 /* wrapper for maximizing: just flip the sign of f and grad */
804 static double f_max(unsigned n, const double *x, double *grad, void *data)
806 f_max_data *d = (f_max_data *) data;
807 double val = d->f(n, x, grad, d->f_data);
810 for (i = 0; i < n; ++i)
816 static void pre_max(unsigned n, const double *x, const double *v,
817 double *vpre, void *data)
819 f_max_data *d = (f_max_data *) data;
821 d->pre(n, x, v, vpre, data);
822 for (i = 0; i < n; ++i) vpre[i] = -vpre[i];
826 NLOPT_STDCALL nlopt_optimize(nlopt_opt opt, double *x, double *opt_f)
828 nlopt_func f; void *f_data; nlopt_precond pre;
833 if (!opt || !opt_f || !opt->f) return NLOPT_INVALID_ARGS;
834 f = opt->f; f_data = opt->f_data; pre = opt->pre;
836 /* for maximizing, just minimize the f_max wrapper, which
837 flips the sign of everything */
838 if ((maximize = opt->maximize)) {
839 fmd.f = f; fmd.f_data = f_data; fmd.pre = pre;
840 opt->f = f_max; opt->f_data = &fmd;
841 if (opt->pre) opt->pre = pre_max;
842 opt->stopval = -opt->stopval;
846 { /* possibly eliminate lb == ub dimensions for some algorithms */
847 nlopt_opt elim_opt = opt;
848 if (elimdim_wrapcheck(opt)) {
849 elim_opt = elimdim_create(opt);
850 if (!elim_opt) { ret = NLOPT_OUT_OF_MEMORY; goto done; }
851 elimdim_shrink(opt->n, x, opt->lb, opt->ub);
854 ret = nlopt_optimize_(elim_opt, x, opt_f);
856 if (elim_opt != opt) {
857 elimdim_destroy(elim_opt);
858 elimdim_expand(opt->n, x, opt->lb, opt->ub);
863 if (maximize) { /* restore original signs */
864 opt->maximize = maximize;
865 opt->stopval = -opt->stopval;
866 opt->f = f; opt->f_data = f_data; opt->pre = pre;
873 /*********************************************************************/
875 nlopt_result nlopt_optimize_limited(nlopt_opt opt, double *x, double *minf,
876 int maxeval, double maxtime)
882 if (!opt) return NLOPT_INVALID_ARGS;
884 save_maxeval = nlopt_get_maxeval(opt);
885 save_maxtime = nlopt_get_maxtime(opt);
887 /* override opt limits if maxeval and/or maxtime are more stringent */
888 if (save_maxeval <= 0 || (maxeval > 0 && maxeval < save_maxeval))
889 nlopt_set_maxeval(opt, maxeval);
890 if (save_maxtime <= 0 || (maxtime > 0 && maxtime < save_maxtime))
891 nlopt_set_maxtime(opt, maxtime);
893 ret = nlopt_optimize(opt, x, minf);
895 nlopt_set_maxeval(opt, save_maxeval);
896 nlopt_set_maxtime(opt, save_maxtime);
901 /*********************************************************************/