1 /* Copyright (c) 2007-2010 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.
29 #include "nlopt-internal.h"
31 /*************************************************************************/
33 void NLOPT_STDCALL nlopt_destroy(nlopt_opt opt)
37 if (opt->munge_on_destroy) {
38 nlopt_munge munge = opt->munge_on_destroy;
40 for (i = 0; i < opt->m; ++i)
41 munge(opt->fc[i].f_data);
42 for (i = 0; i < opt->p; ++i)
43 munge(opt->h[i].f_data);
45 for (i = 0; i < opt->m; ++i)
47 for (i = 0; i < opt->p; ++i)
49 free(opt->lb); free(opt->ub);
53 nlopt_destroy(opt->local_opt);
60 nlopt_opt NLOPT_STDCALL nlopt_create(nlopt_algorithm algorithm, unsigned n)
64 if (((int) algorithm) < 0 || algorithm >= NLOPT_NUM_ALGORITHMS)
67 opt = (nlopt_opt) malloc(sizeof(struct nlopt_opt_s));
69 opt->algorithm = algorithm;
71 opt->f = NULL; opt->f_data = NULL;
73 opt->munge_on_destroy = opt->munge_on_copy = NULL;
75 opt->lb = opt->ub = NULL;
76 opt->m = opt->m_alloc = 0;
78 opt->p = opt->p_alloc = 0;
81 opt->stopval = -HUGE_VAL;
82 opt->ftol_rel = opt->ftol_abs = 0;
83 opt->xtol_rel = 0; opt->xtol_abs = NULL;
87 opt->force_stop_child = NULL;
89 opt->local_opt = NULL;
90 opt->stochastic_population = 0;
95 opt->lb = (double *) malloc(sizeof(double) * (n));
96 if (!opt->lb) goto oom;
97 opt->ub = (double *) malloc(sizeof(double) * (n));
98 if (!opt->ub) goto oom;
99 opt->xtol_abs = (double *) malloc(sizeof(double) * (n));
100 if (!opt->xtol_abs) goto oom;
101 nlopt_set_lower_bounds1(opt, -HUGE_VAL);
102 nlopt_set_upper_bounds1(opt, +HUGE_VAL);
103 nlopt_set_xtol_abs1(opt, 0.0);
114 nlopt_opt NLOPT_STDCALL nlopt_copy(const nlopt_opt opt)
116 nlopt_opt nopt = NULL;
119 nopt = (nlopt_opt) malloc(sizeof(struct nlopt_opt_s));
121 nopt->lb = nopt->ub = nopt->xtol_abs = NULL;
122 nopt->fc = nopt->h = NULL;
123 nopt->m_alloc = nopt->p_alloc = 0;
124 nopt->local_opt = NULL;
127 opt->force_stop_child = NULL;
129 nlopt_munge munge = nopt->munge_on_copy;
130 if (munge && nopt->f_data)
131 if (!(nopt->f_data = munge(nopt->f_data))) goto oom;
134 nopt->lb = (double *) malloc(sizeof(double) * (opt->n));
135 if (!opt->lb) goto oom;
136 nopt->ub = (double *) malloc(sizeof(double) * (opt->n));
137 if (!opt->ub) goto oom;
138 nopt->xtol_abs = (double *) malloc(sizeof(double) * (opt->n));
139 if (!opt->xtol_abs) goto oom;
141 memcpy(nopt->lb, opt->lb, sizeof(double) * (opt->n));
142 memcpy(nopt->ub, opt->ub, sizeof(double) * (opt->n));
143 memcpy(nopt->xtol_abs, opt->xtol_abs, sizeof(double) * (opt->n));
147 nopt->m_alloc = opt->m;
148 nopt->fc = (nlopt_constraint *) malloc(sizeof(nlopt_constraint)
150 if (!nopt->fc) goto oom;
151 memcpy(nopt->fc, opt->fc, sizeof(nlopt_constraint) * (opt->m));
152 for (i = 0; i < opt->m; ++i) nopt->fc[i].tol = NULL;
154 for (i = 0; i < opt->m; ++i)
155 if (nopt->fc[i].f_data &&
157 = munge(nopt->fc[i].f_data)))
159 for (i = 0; i < opt->m; ++i)
160 if (opt->fc[i].tol) {
161 nopt->fc[i].tol = (double *) malloc(sizeof(double)
163 if (!nopt->fc[i].tol) goto oom;
164 memcpy(nopt->fc[i].tol, opt->fc[i].tol,
165 sizeof(double) * nopt->fc[i].m);
170 nopt->p_alloc = opt->p;
171 nopt->h = (nlopt_constraint *) malloc(sizeof(nlopt_constraint)
173 if (!nopt->h) goto oom;
174 memcpy(nopt->h, opt->h, sizeof(nlopt_constraint) * (opt->p));
175 for (i = 0; i < opt->p; ++i) nopt->h[i].tol = NULL;
177 for (i = 0; i < opt->p; ++i)
178 if (nopt->h[i].f_data &&
180 = munge(nopt->h[i].f_data)))
182 for (i = 0; i < opt->p; ++i)
184 nopt->h[i].tol = (double *) malloc(sizeof(double)
186 if (!nopt->h[i].tol) goto oom;
187 memcpy(nopt->h[i].tol, opt->h[i].tol,
188 sizeof(double) * nopt->h[i].m);
192 if (opt->local_opt) {
193 nopt->local_opt = nlopt_copy(opt->local_opt);
194 if (!nopt->local_opt) goto oom;
198 nopt->dx = (double *) malloc(sizeof(double) * (opt->n));
199 if (!nopt->dx) goto oom;
200 memcpy(nopt->dx, opt->dx, sizeof(double) * (opt->n));
206 nopt->munge_on_destroy = NULL; // better to leak mem than crash
211 /*************************************************************************/
213 nlopt_result NLOPT_STDCALL nlopt_set_min_objective(nlopt_opt opt,
214 nlopt_func f, void *f_data)
217 if (opt->munge_on_destroy) opt->munge_on_destroy(opt->f_data);
218 opt->f = f; opt->f_data = f_data;
220 if (nlopt_isinf(opt->stopval) && opt->stopval > 0)
221 opt->stopval = -HUGE_VAL; /* switch default from max to min */
222 return NLOPT_SUCCESS;
224 return NLOPT_INVALID_ARGS;
227 nlopt_result NLOPT_STDCALL nlopt_set_max_objective(nlopt_opt opt,
228 nlopt_func f, void *f_data)
231 if (opt->munge_on_destroy) opt->munge_on_destroy(opt->f_data);
232 opt->f = f; opt->f_data = f_data;
234 if (nlopt_isinf(opt->stopval) && opt->stopval < 0)
235 opt->stopval = +HUGE_VAL; /* switch default from min to max */
236 return NLOPT_SUCCESS;
238 return NLOPT_INVALID_ARGS;
241 /*************************************************************************/
244 NLOPT_STDCALL nlopt_set_lower_bounds(nlopt_opt opt, const double *lb)
246 if (opt && (opt->n == 0 || lb)) {
247 memcpy(opt->lb, lb, sizeof(double) * (opt->n));
248 return NLOPT_SUCCESS;
250 return NLOPT_INVALID_ARGS;
254 NLOPT_STDCALL nlopt_set_lower_bounds1(nlopt_opt opt, double lb)
258 for (i = 0; i < opt->n; ++i)
260 return NLOPT_SUCCESS;
262 return NLOPT_INVALID_ARGS;
266 NLOPT_STDCALL nlopt_get_lower_bounds(nlopt_opt opt, double *lb)
268 if (opt && (opt->n == 0 || lb)) {
269 memcpy(lb, opt->lb, sizeof(double) * (opt->n));
270 return NLOPT_SUCCESS;
272 return NLOPT_INVALID_ARGS;
276 NLOPT_STDCALL nlopt_set_upper_bounds(nlopt_opt opt, const double *ub)
278 if (opt && (opt->n == 0 || ub)) {
279 memcpy(opt->ub, ub, sizeof(double) * (opt->n));
280 return NLOPT_SUCCESS;
282 return NLOPT_INVALID_ARGS;
286 NLOPT_STDCALL nlopt_set_upper_bounds1(nlopt_opt opt, double ub)
290 for (i = 0; i < opt->n; ++i)
292 return NLOPT_SUCCESS;
294 return NLOPT_INVALID_ARGS;
298 NLOPT_STDCALL nlopt_get_upper_bounds(nlopt_opt opt, double *ub)
300 if (opt && (opt->n == 0 || ub)) {
301 memcpy(ub, opt->ub, sizeof(double) * (opt->n));
302 return NLOPT_SUCCESS;
304 return NLOPT_INVALID_ARGS;
307 /*************************************************************************/
309 #define AUGLAG_ALG(a) ((a) == NLOPT_AUGLAG || \
310 (a) == NLOPT_AUGLAG_EQ || \
311 (a) == NLOPT_LN_AUGLAG || \
312 (a) == NLOPT_LN_AUGLAG_EQ || \
313 (a) == NLOPT_LD_AUGLAG || \
314 (a) == NLOPT_LD_AUGLAG_EQ)
317 NLOPT_STDCALL nlopt_remove_inequality_constraints(nlopt_opt opt)
320 if (!opt) return NLOPT_INVALID_ARGS;
321 if (opt->munge_on_destroy) {
322 nlopt_munge munge = opt->munge_on_destroy;
323 for (i = 0; i < opt->m; ++i)
324 munge(opt->fc[i].f_data);
326 for (i = 0; i < opt->m; ++i)
327 free(opt->fc[i].tol);
330 opt->m = opt->m_alloc = 0;
331 return NLOPT_SUCCESS;
334 static nlopt_result add_constraint(unsigned *m, unsigned *m_alloc,
335 nlopt_constraint **c,
336 unsigned fm, nlopt_func fc, nlopt_mfunc mfc,
343 if ((fc && mfc) || (fc && fm != 1) || (!fc && !mfc) || !tol)
344 return NLOPT_INVALID_ARGS;
345 for (i = 0; i < fm; ++i) if (tol[i] < 0) return NLOPT_INVALID_ARGS;
347 tolcopy = (double *) malloc(sizeof(double) * fm);
348 if (fm && !tolcopy) return NLOPT_OUT_OF_MEMORY;
349 memcpy(tolcopy, tol, sizeof(double) * fm);
353 /* allocate by repeated doubling so that
354 we end up with O(log m) mallocs rather than O(m). */
356 *c = (nlopt_constraint *) realloc(*c,
357 sizeof(nlopt_constraint)
362 return NLOPT_OUT_OF_MEMORY;
368 (*c)[*m - 1].mf = mfc;
369 (*c)[*m - 1].f_data = fc_data;
370 (*c)[*m - 1].tol = tolcopy;
371 return NLOPT_SUCCESS;
374 static int inequality_ok(nlopt_algorithm algorithm) {
375 /* nonlinear constraints are only supported with some algorithms */
376 return (algorithm == NLOPT_LD_MMA
377 || algorithm == NLOPT_LD_SLSQP
378 || algorithm == NLOPT_LN_COBYLA
379 || AUGLAG_ALG(algorithm)
380 || algorithm == NLOPT_GN_ISRES
381 || algorithm == NLOPT_GN_ORIG_DIRECT
382 || algorithm == NLOPT_GN_ORIG_DIRECT_L);
386 NLOPT_STDCALL nlopt_add_inequality_mconstraint(nlopt_opt opt, unsigned m,
387 nlopt_mfunc fc, void *fc_data,
391 if (!m) { /* empty constraints are always ok */
392 if (opt && opt->munge_on_destroy) opt->munge_on_destroy(fc_data);
393 return NLOPT_SUCCESS;
395 if (!opt || !inequality_ok(opt->algorithm)) ret = NLOPT_INVALID_ARGS;
396 else ret = add_constraint(&opt->m, &opt->m_alloc, &opt->fc,
397 m, NULL, fc, fc_data, tol);
398 if (ret < 0 && opt && opt->munge_on_destroy)
399 opt->munge_on_destroy(fc_data);
404 NLOPT_STDCALL nlopt_add_inequality_constraint(nlopt_opt opt,
405 nlopt_func fc, void *fc_data,
409 if (!opt || !inequality_ok(opt->algorithm)) ret = NLOPT_INVALID_ARGS;
410 else ret = add_constraint(&opt->m, &opt->m_alloc, &opt->fc,
411 1, fc, NULL, fc_data, &tol);
412 if (ret < 0 && opt && opt->munge_on_destroy)
413 opt->munge_on_destroy(fc_data);
418 NLOPT_STDCALL nlopt_remove_equality_constraints(nlopt_opt opt)
421 if (!opt) return NLOPT_INVALID_ARGS;
422 if (opt->munge_on_destroy) {
423 nlopt_munge munge = opt->munge_on_destroy;
424 for (i = 0; i < opt->p; ++i)
425 munge(opt->h[i].f_data);
427 for (i = 0; i < opt->p; ++i)
431 opt->p = opt->p_alloc = 0;
432 return NLOPT_SUCCESS;
435 static int equality_ok(nlopt_algorithm algorithm) {
436 /* equality constraints (h(x) = 0) only via some algorithms */
437 return (AUGLAG_ALG(algorithm)
438 || algorithm == NLOPT_LD_SLSQP
439 || algorithm == NLOPT_GN_ISRES
440 || algorithm == NLOPT_LN_COBYLA);
444 NLOPT_STDCALL nlopt_add_equality_mconstraint(nlopt_opt opt, unsigned m,
445 nlopt_mfunc fc, void *fc_data,
449 if (!m) { /* empty constraints are always ok */
450 if (opt && opt->munge_on_destroy) opt->munge_on_destroy(fc_data);
451 return NLOPT_SUCCESS;
453 if (!opt || !equality_ok(opt->algorithm)
454 || nlopt_count_constraints(opt->p, opt->h) + m > opt->n)
455 ret = NLOPT_INVALID_ARGS;
456 else ret = add_constraint(&opt->p, &opt->p_alloc, &opt->h,
457 m, NULL, fc, fc_data, tol);
458 if (ret < 0 && opt && opt->munge_on_destroy)
459 opt->munge_on_destroy(fc_data);
464 NLOPT_STDCALL nlopt_add_equality_constraint(nlopt_opt opt,
465 nlopt_func fc, void *fc_data,
469 if (!opt || !equality_ok(opt->algorithm)
470 || nlopt_count_constraints(opt->p, opt->h) + 1 > opt->n)
471 ret = NLOPT_INVALID_ARGS;
472 else ret = add_constraint(&opt->p, &opt->p_alloc, &opt->h,
473 1, fc, NULL, fc_data, &tol);
474 if (ret < 0 && opt && opt->munge_on_destroy)
475 opt->munge_on_destroy(fc_data);
479 /*************************************************************************/
481 #define SET(param, T, arg) \
482 nlopt_result NLOPT_STDCALL nlopt_set_##param(nlopt_opt opt, T arg) \
486 return NLOPT_SUCCESS; \
488 return NLOPT_INVALID_ARGS; \
492 #define GET(param, T, arg) T NLOPT_STDCALL \
493 nlopt_get_##param(const nlopt_opt opt) { \
497 #define GETSET(param, T, arg) GET(param, T, arg) SET(param, T, arg)
499 GETSET(stopval, double, stopval)
501 GETSET(ftol_rel, double, ftol_rel)
502 GETSET(ftol_abs, double, ftol_abs)
503 GETSET(xtol_rel, double, xtol_rel)
506 NLOPT_STDCALL nlopt_set_xtol_abs(nlopt_opt opt, const double *xtol_abs)
509 memcpy(opt->xtol_abs, xtol_abs, opt->n & sizeof(double));
510 return NLOPT_SUCCESS;
512 return NLOPT_INVALID_ARGS;
516 NLOPT_STDCALL nlopt_set_xtol_abs1(nlopt_opt opt, const double xtol_abs)
520 for (i = 0; i < opt->n; ++i)
521 opt->xtol_abs[i] = xtol_abs;
522 return NLOPT_SUCCESS;
524 return NLOPT_INVALID_ARGS;
528 NLOPT_STDCALL nlopt_get_xtol_abs(const nlopt_opt opt, double *xtol_abs)
530 memcpy(xtol_abs, opt->xtol_abs, opt->n & sizeof(double));
531 return NLOPT_SUCCESS;
534 GETSET(maxeval, int, maxeval)
536 GETSET(maxtime, double, maxtime)
538 /*************************************************************************/
541 NLOPT_STDCALL nlopt_set_force_stop(nlopt_opt opt, int force_stop)
544 opt->force_stop = force_stop;
545 if (opt->force_stop_child)
546 return nlopt_set_force_stop(opt->force_stop_child, force_stop);
547 return NLOPT_SUCCESS;
549 return NLOPT_INVALID_ARGS;
552 GET(force_stop, int, force_stop)
553 nlopt_result NLOPT_STDCALL nlopt_force_stop(nlopt_opt opt) {
554 return nlopt_set_force_stop(opt, 1);
557 /*************************************************************************/
559 GET(algorithm, nlopt_algorithm, algorithm)
560 GET(dimension, unsigned, n)
562 /*************************************************************************/
565 NLOPT_STDCALL nlopt_set_local_optimizer(nlopt_opt opt,
566 const nlopt_opt local_opt)
569 if (local_opt && local_opt->n != opt->n) return NLOPT_INVALID_ARGS;
570 nlopt_destroy(opt->local_opt);
571 opt->local_opt = nlopt_copy(local_opt);
573 if (!opt->local_opt) return NLOPT_OUT_OF_MEMORY;
574 nlopt_set_lower_bounds(opt->local_opt, opt->lb);
575 nlopt_set_upper_bounds(opt->local_opt, opt->ub);
576 nlopt_remove_inequality_constraints(opt->local_opt);
577 nlopt_remove_equality_constraints(opt->local_opt);
578 nlopt_set_min_objective(opt->local_opt, NULL, NULL);
579 opt->local_opt->force_stop = 0;
581 return NLOPT_SUCCESS;
583 return NLOPT_INVALID_ARGS;
586 /*************************************************************************/
588 GETSET(population, unsigned, stochastic_population)
590 /*************************************************************************/
592 nlopt_result NLOPT_STDCALL nlopt_set_initial_step1(nlopt_opt opt, double dx)
595 if (!opt || dx == 0) return NLOPT_INVALID_ARGS;
596 if (!opt->dx && opt->n > 0) {
597 opt->dx = (double *) malloc(sizeof(double) * (opt->n));
598 if (!opt->dx) return NLOPT_OUT_OF_MEMORY;
600 for (i = 0; i < opt->n; ++i) opt->dx[i] = dx;
601 return NLOPT_SUCCESS;
605 NLOPT_STDCALL nlopt_set_initial_step(nlopt_opt opt, const double *dx)
608 if (!opt) return NLOPT_INVALID_ARGS;
610 free(opt->dx); opt->dx = NULL;
611 return NLOPT_SUCCESS;
613 for (i = 0; i < opt->n; ++i) if (dx[i] == 0) return NLOPT_INVALID_ARGS;
614 if (!opt->dx && nlopt_set_initial_step1(opt, 1) == NLOPT_OUT_OF_MEMORY)
615 return NLOPT_OUT_OF_MEMORY;
616 memcpy(opt->dx, dx, sizeof(double) * (opt->n));
617 return NLOPT_SUCCESS;
621 NLOPT_STDCALL nlopt_get_initial_step(const nlopt_opt opt, const double *x,
624 if (!opt) return NLOPT_INVALID_ARGS;
625 if (!opt->n) return NLOPT_SUCCESS;
627 nlopt_opt o = (nlopt_opt) opt; /* discard const temporarily */
628 nlopt_result ret = nlopt_set_default_initial_step(o, x);
629 if (ret != NLOPT_SUCCESS) return ret;
630 memcpy(dx, o->dx, sizeof(double) * (opt->n));
631 free(o->dx); o->dx = NULL; /* don't save, since x-dependent */
634 memcpy(dx, opt->dx, sizeof(double) * (opt->n));
635 return NLOPT_SUCCESS;
639 NLOPT_STDCALL nlopt_set_default_initial_step(nlopt_opt opt, const double *x)
641 const double *lb, *ub;
644 if (!opt || !x) return NLOPT_INVALID_ARGS;
645 lb = opt->lb; ub = opt->ub;
647 if (!opt->dx && nlopt_set_initial_step1(opt, 1) == NLOPT_OUT_OF_MEMORY)
648 return NLOPT_OUT_OF_MEMORY;
650 /* crude heuristics for initial step size of nonderivative algorithms */
651 for (i = 0; i < opt->n; ++i) {
652 double step = HUGE_VAL;
654 if (!nlopt_isinf(ub[i]) && !nlopt_isinf(lb[i])
655 && (ub[i] - lb[i]) * 0.25 < step && ub[i] > lb[i])
656 step = (ub[i] - lb[i]) * 0.25;
657 if (!nlopt_isinf(ub[i])
658 && ub[i] - x[i] < step && ub[i] > x[i])
659 step = (ub[i] - x[i]) * 0.75;
660 if (!nlopt_isinf(lb[i])
661 && x[i] - lb[i] < step && x[i] > lb[i])
662 step = (x[i] - lb[i]) * 0.75;
664 if (nlopt_isinf(step)) {
665 if (!nlopt_isinf(ub[i])
666 && fabs(ub[i] - x[i]) < fabs(step))
667 step = (ub[i] - x[i]) * 1.1;
668 if (!nlopt_isinf(lb[i])
669 && fabs(x[i] - lb[i]) < fabs(step))
670 step = (x[i] - lb[i]) * 1.1;
672 if (nlopt_isinf(step) || step == 0) {
675 if (nlopt_isinf(step) || step == 0)
680 return NLOPT_SUCCESS;
683 /*************************************************************************/
685 void NLOPT_STDCALL nlopt_set_munge(nlopt_opt opt,
686 nlopt_munge munge_on_destroy,
687 nlopt_munge munge_on_copy) {
689 opt->munge_on_destroy = munge_on_destroy;
690 opt->munge_on_copy = munge_on_copy;
694 /*************************************************************************/