chiark / gitweb /
add mconstraint to C API (not yet tested), thanks to Dmitrey of OpenOpt for the sugge...
[nlopt.git] / api / options.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 <stdio.h>
24 #include <stdlib.h>
25 #include <math.h>
26 #include <string.h>
27 #include <float.h>
28
29 #include "nlopt-internal.h"
30
31 /*************************************************************************/
32
33 void NLOPT_STDCALL nlopt_destroy(nlopt_opt opt)
34 {
35      if (opt) {
36           unsigned i;
37           if (opt->munge_on_destroy) {
38                nlopt_munge munge = opt->munge_on_destroy;
39                munge(opt->f_data);
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);
44           }
45           for (i = 0; i < opt->m; ++i)
46                free(opt->fc[i].tol);
47           for (i = 0; i < opt->p; ++i)
48                free(opt->h[i].tol);
49           free(opt->lb); free(opt->ub);
50           free(opt->xtol_abs);
51           free(opt->fc);
52           free(opt->h);
53           nlopt_destroy(opt->local_opt);
54           free(opt->dx);
55           free(opt->work);
56           free(opt);
57      }
58 }
59
60 nlopt_opt NLOPT_STDCALL nlopt_create(nlopt_algorithm algorithm, unsigned n)
61 {
62      nlopt_opt opt;
63
64      if (((int) algorithm) < 0 || algorithm >= NLOPT_NUM_ALGORITHMS)
65           return NULL;
66
67      opt = (nlopt_opt) malloc(sizeof(struct nlopt_opt_s));
68      if (opt) {
69           opt->algorithm = algorithm;
70           opt->n = n;
71           opt->f = NULL; opt->f_data = NULL;
72           opt->maximize = 0;
73           opt->munge_on_destroy = opt->munge_on_copy = NULL;
74
75           opt->lb = opt->ub = NULL;
76           opt->m = opt->m_alloc = 0;
77           opt->fc = NULL;
78           opt->p = opt->p_alloc = 0;
79           opt->h = NULL;
80
81           opt->stopval = -HUGE_VAL;
82           opt->ftol_rel = opt->ftol_abs = 0;
83           opt->xtol_rel = 0; opt->xtol_abs = NULL;
84           opt->maxeval = 0;
85           opt->maxtime = 0;
86           opt->force_stop = 0;
87           opt->force_stop_child = NULL;
88
89           opt->local_opt = NULL;
90           opt->stochastic_population = 0;
91           opt->dx = NULL;
92           opt->work = NULL;
93
94           if (n > 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);
104           }
105      }
106
107      return opt;
108
109 oom:
110      nlopt_destroy(opt);
111      return NULL;
112 }
113
114 nlopt_opt NLOPT_STDCALL nlopt_copy(const nlopt_opt opt)
115 {
116      nlopt_opt nopt = NULL;
117      unsigned i;
118      if (opt) {
119           nopt = (nlopt_opt) malloc(sizeof(struct nlopt_opt_s));
120           *nopt = *opt;
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;
125           nopt->dx = NULL;
126           nopt->work = NULL;
127           opt->force_stop_child = NULL;
128
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;
132
133           if (opt->n > 0) {
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;
140                
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));
144           }
145
146           if (opt->m) {
147                nopt->m_alloc = opt->m;
148                nopt->fc = (nlopt_constraint *) malloc(sizeof(nlopt_constraint)
149                                                       * (opt->m));
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;
153                if (munge)
154                     for (i = 0; i < opt->m; ++i)
155                          if (nopt->fc[i].f_data &&
156                              !(nopt->fc[i].f_data
157                                = munge(nopt->fc[i].f_data)))
158                               goto oom;
159                for (i = 0; i < opt->m; ++i)
160                     if (opt->fc[i].tol) {
161                          nopt->fc[i].tol = (double *) malloc(sizeof(double)
162                                                              * nopt->fc[i].m);
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);
166                     }
167           }
168
169           if (opt->p) {
170                nopt->p_alloc = opt->p;
171                nopt->h = (nlopt_constraint *) malloc(sizeof(nlopt_constraint)
172                                                      * (opt->p));
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;
176                if (munge)
177                     for (i = 0; i < opt->p; ++i)
178                          if (nopt->h[i].f_data &&
179                              !(nopt->h[i].f_data
180                                = munge(nopt->h[i].f_data)))
181                               goto oom;
182                for (i = 0; i < opt->p; ++i)
183                     if (opt->h[i].tol) {
184                          nopt->h[i].tol = (double *) malloc(sizeof(double)
185                                                              * nopt->h[i].m);
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);
189                     }
190           }
191
192           if (opt->local_opt) {
193                nopt->local_opt = nlopt_copy(opt->local_opt);
194                if (!nopt->local_opt) goto oom;
195           }
196
197           if (opt->dx) {
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));
201           }
202      }
203      return nopt;
204
205 oom:
206      nopt->munge_on_destroy = NULL; // better to leak mem than crash
207      nlopt_destroy(nopt);
208      return NULL;
209 }
210
211 /*************************************************************************/
212
213 nlopt_result NLOPT_STDCALL nlopt_set_min_objective(nlopt_opt opt,
214                                                    nlopt_func f, void *f_data)
215 {
216      if (opt) {
217           opt->f = f; opt->f_data = f_data;
218           opt->maximize = 0;
219           if (nlopt_isinf(opt->stopval) && opt->stopval > 0)
220                opt->stopval = -HUGE_VAL; /* switch default from max to min */
221           return NLOPT_SUCCESS;
222      }
223      return NLOPT_INVALID_ARGS;
224 }
225
226 nlopt_result NLOPT_STDCALL nlopt_set_max_objective(nlopt_opt opt, 
227                                                    nlopt_func f, void *f_data)
228 {
229      if (opt) {
230           opt->f = f; opt->f_data = f_data;
231           opt->maximize = 1;
232           if (nlopt_isinf(opt->stopval) && opt->stopval < 0)
233                opt->stopval = +HUGE_VAL; /* switch default from min to max */
234           return NLOPT_SUCCESS;
235      }
236      return NLOPT_INVALID_ARGS;
237 }
238
239 /*************************************************************************/
240
241 nlopt_result
242 NLOPT_STDCALL nlopt_set_lower_bounds(nlopt_opt opt, const double *lb)
243 {
244      if (opt && (opt->n == 0 || lb)) {
245           memcpy(opt->lb, lb, sizeof(double) * (opt->n));
246           return NLOPT_SUCCESS;
247      }
248      return NLOPT_INVALID_ARGS;
249 }
250
251 nlopt_result
252 NLOPT_STDCALL nlopt_set_lower_bounds1(nlopt_opt opt, double lb)
253 {
254      if (opt) {
255           unsigned i;
256           for (i = 0; i < opt->n; ++i)
257                opt->lb[i] = lb;
258           return NLOPT_SUCCESS;
259      }
260      return NLOPT_INVALID_ARGS;
261 }
262
263 nlopt_result
264 NLOPT_STDCALL nlopt_get_lower_bounds(nlopt_opt opt, double *lb)
265 {
266      if (opt && (opt->n == 0 || lb)) {
267           memcpy(lb, opt->lb, sizeof(double) * (opt->n));
268           return NLOPT_SUCCESS;
269      }
270      return NLOPT_INVALID_ARGS;
271 }
272
273 nlopt_result
274 NLOPT_STDCALL nlopt_set_upper_bounds(nlopt_opt opt, const double *ub)
275 {
276      if (opt && (opt->n == 0 || ub)) {
277           memcpy(opt->ub, ub, sizeof(double) * (opt->n));
278           return NLOPT_SUCCESS;
279      }
280      return NLOPT_INVALID_ARGS;
281 }
282
283 nlopt_result
284 NLOPT_STDCALL nlopt_set_upper_bounds1(nlopt_opt opt, double ub)
285 {
286      if (opt) {
287           unsigned i;
288           for (i = 0; i < opt->n; ++i)
289                opt->ub[i] = ub;
290           return NLOPT_SUCCESS;
291      }
292      return NLOPT_INVALID_ARGS;
293 }
294
295 nlopt_result
296 NLOPT_STDCALL nlopt_get_upper_bounds(nlopt_opt opt, double *ub)
297 {
298      if (opt && (opt->n == 0 || ub)) {
299           memcpy(ub, opt->ub, sizeof(double) * (opt->n));
300           return NLOPT_SUCCESS;
301      }
302      return NLOPT_INVALID_ARGS;
303 }
304
305 /*************************************************************************/
306
307 #define AUGLAG_ALG(a) ((a) == NLOPT_AUGLAG ||           \
308                        (a) == NLOPT_AUGLAG_EQ ||        \
309                        (a) == NLOPT_LN_AUGLAG ||        \
310                        (a) == NLOPT_LN_AUGLAG_EQ ||     \
311                        (a) == NLOPT_LD_AUGLAG ||        \
312                        (a) == NLOPT_LD_AUGLAG_EQ)
313
314 nlopt_result
315 NLOPT_STDCALL nlopt_remove_inequality_constraints(nlopt_opt opt)
316 {
317      unsigned i;
318      if (!opt) return NLOPT_INVALID_ARGS;
319      if (opt->munge_on_destroy) {
320           nlopt_munge munge = opt->munge_on_destroy;
321           for (i = 0; i < opt->m; ++i)
322                munge(opt->fc[i].f_data);
323      }
324      for (i = 0; i < opt->m; ++i)
325           free(opt->fc[i].tol);
326      free(opt->fc);
327      opt->fc = NULL;
328      opt->m = opt->m_alloc = 0;
329      return NLOPT_SUCCESS;
330 }
331
332 static nlopt_result add_constraint(unsigned *m, unsigned *m_alloc,
333                                    nlopt_constraint **c,
334                                    unsigned fm, nlopt_func fc, nlopt_mfunc mfc,
335                                    void *fc_data,
336                                    const double *tol)
337 {
338      double *tolcopy;
339
340      if ((fc && mfc) || (fc && fm != 1) || (!fc && !mfc) || tol < 0)
341           return NLOPT_INVALID_ARGS;
342      
343      tolcopy = (double *) malloc(sizeof(double) * fm);
344      if (fm && !tolcopy) return NLOPT_OUT_OF_MEMORY;
345      memcpy(tolcopy, tol, sizeof(double) * fm);
346
347      *m += 1;
348      if (*m > *m_alloc) {
349           /* allocate by repeated doubling so that 
350              we end up with O(log m) mallocs rather than O(m). */
351           *m_alloc = 2 * (*m);
352           *c = (nlopt_constraint *) realloc(*c,
353                                             sizeof(nlopt_constraint)
354                                             * (*m_alloc));
355           if (!*c) {
356                *m_alloc = *m = 0;
357                free(tolcopy);
358                return NLOPT_OUT_OF_MEMORY;
359           }
360      }
361      
362      (*c)[*m - 1].m = fm;
363      (*c)[*m - 1].f = fc;
364      (*c)[*m - 1].mf = mfc;
365      (*c)[*m - 1].f_data = fc_data;
366      (*c)[*m - 1].tol = tolcopy;
367      return NLOPT_SUCCESS;
368 }
369
370 static int inequality_ok(nlopt_algorithm algorithm) {
371      /* nonlinear constraints are only supported with some algorithms */
372      return (algorithm == NLOPT_LD_MMA 
373              || algorithm == NLOPT_LN_COBYLA
374              || AUGLAG_ALG(algorithm) 
375              || algorithm == NLOPT_GN_ISRES
376              || algorithm == NLOPT_GN_ORIG_DIRECT
377              || algorithm == NLOPT_GN_ORIG_DIRECT_L);
378 }
379
380 nlopt_result
381 NLOPT_STDCALL nlopt_add_inequality_mconstraint(nlopt_opt opt, unsigned m,
382                                               nlopt_mfunc fc, void *fc_data,
383                                                const double *tol)
384 {
385      if (!m) return NLOPT_SUCCESS; /* empty constraints are always ok */
386      if (!opt || !inequality_ok(opt->algorithm)) return NLOPT_INVALID_ARGS;
387      return add_constraint(&opt->m, &opt->m_alloc, &opt->fc,
388                            m, NULL, fc, fc_data, tol);
389 }
390
391 nlopt_result
392 NLOPT_STDCALL nlopt_add_inequality_constraint(nlopt_opt opt,
393                                               nlopt_func fc, void *fc_data,
394                                               double tol)
395 {
396      if (!opt || !inequality_ok(opt->algorithm)) return NLOPT_INVALID_ARGS;
397      return add_constraint(&opt->m, &opt->m_alloc, &opt->fc,
398                            1, fc, NULL, fc_data, &tol);
399 }
400
401 nlopt_result
402 NLOPT_STDCALL nlopt_remove_equality_constraints(nlopt_opt opt)
403 {
404      unsigned i;
405      if (!opt) return NLOPT_INVALID_ARGS;
406      if (opt->munge_on_destroy) {
407           nlopt_munge munge = opt->munge_on_destroy;
408           for (i = 0; i < opt->p; ++i)
409                munge(opt->h[i].f_data);
410      }
411      for (i = 0; i < opt->p; ++i)
412           free(opt->h[i].tol);
413      free(opt->h);
414      opt->h = NULL;
415      opt->p = opt->p_alloc = 0;
416      return NLOPT_SUCCESS;
417 }
418
419 static int equality_ok(nlopt_algorithm algorithm) {
420      /* equality constraints (h(x) = 0) only via some algorithms */
421      return (AUGLAG_ALG(algorithm) 
422              || algorithm == NLOPT_GN_ISRES
423              || algorithm == NLOPT_LN_COBYLA);
424 }
425
426 nlopt_result
427 NLOPT_STDCALL nlopt_add_equality_mconstraint(nlopt_opt opt, unsigned m,
428                                              nlopt_mfunc fc, void *fc_data,
429                                              const double *tol)
430 {
431      if (!m) return NLOPT_SUCCESS; /* empty constraints are always ok */
432      if (!opt || !equality_ok(opt->algorithm)) return NLOPT_INVALID_ARGS;
433      return add_constraint(&opt->p, &opt->p_alloc, &opt->h,
434                            m, NULL, fc, fc_data, tol);
435 }
436
437 nlopt_result
438 NLOPT_STDCALL nlopt_add_equality_constraint(nlopt_opt opt,
439                                             nlopt_func fc, void *fc_data,
440                                             double tol)
441 {
442      if (!opt || !equality_ok(opt->algorithm)) return NLOPT_INVALID_ARGS;
443      return add_constraint(&opt->p, &opt->p_alloc, &opt->h,
444                            1, fc, NULL, fc_data, &tol);
445 }
446
447 /*************************************************************************/
448
449 #define SET(param, T, arg)                                              \
450    nlopt_result NLOPT_STDCALL nlopt_set_##param(nlopt_opt opt, T arg)   \
451    {                                                                    \
452         if (opt) {                                                      \
453              opt->arg = arg;                                            \
454              return NLOPT_SUCCESS;                                      \
455         }                                                               \
456         return NLOPT_INVALID_ARGS;                                      \
457    }
458
459
460 #define GET(param, T, arg) T NLOPT_STDCALL      \
461    nlopt_get_##param(const nlopt_opt opt) {     \
462         return opt->arg;                        \
463    }
464
465 #define GETSET(param, T, arg) GET(param, T, arg) SET(param, T, arg)
466
467 GETSET(stopval, double, stopval)
468
469 GETSET(ftol_rel, double, ftol_rel)
470 GETSET(ftol_abs, double, ftol_abs)
471 GETSET(xtol_rel, double, xtol_rel)
472
473 nlopt_result
474 NLOPT_STDCALL nlopt_set_xtol_abs(nlopt_opt opt, const double *xtol_abs)
475 {
476      if (opt) {
477           memcpy(opt->xtol_abs, xtol_abs, opt->n & sizeof(double));
478           return NLOPT_SUCCESS;
479      }
480      return NLOPT_INVALID_ARGS;
481 }
482
483 nlopt_result
484 NLOPT_STDCALL nlopt_set_xtol_abs1(nlopt_opt opt, const double xtol_abs)
485 {
486      if (opt) {
487           unsigned i;
488           for (i = 0; i < opt->n; ++i)
489                opt->xtol_abs[i] = xtol_abs;
490           return NLOPT_SUCCESS;
491      }
492      return NLOPT_INVALID_ARGS;
493 }
494
495 nlopt_result
496 NLOPT_STDCALL nlopt_get_xtol_abs(const nlopt_opt opt, double *xtol_abs)
497 {
498      memcpy(xtol_abs, opt->xtol_abs, opt->n & sizeof(double));
499      return NLOPT_SUCCESS;
500 }
501
502 GETSET(maxeval, int, maxeval)
503
504 GETSET(maxtime, double, maxtime)
505
506 /*************************************************************************/
507
508 nlopt_result
509 NLOPT_STDCALL nlopt_set_force_stop(nlopt_opt opt, int force_stop)
510 {
511      if (opt) {
512           opt->force_stop = force_stop;
513           if (opt->force_stop_child)
514                return nlopt_set_force_stop(opt->force_stop_child, force_stop);
515           return NLOPT_SUCCESS;
516      }
517      return NLOPT_INVALID_ARGS;
518 }
519
520 GET(force_stop, int, force_stop)
521 nlopt_result NLOPT_STDCALL nlopt_force_stop(nlopt_opt opt) { 
522      return nlopt_set_force_stop(opt, 1); 
523 }
524
525 /*************************************************************************/
526
527 GET(algorithm, nlopt_algorithm, algorithm)
528 GET(dimension, unsigned, n)
529
530 /*************************************************************************/
531
532 nlopt_result
533 NLOPT_STDCALL nlopt_set_local_optimizer(nlopt_opt opt,
534                                         const nlopt_opt local_opt)
535 {
536      if (opt) {
537           if (local_opt && local_opt->n != opt->n) return NLOPT_INVALID_ARGS;
538           nlopt_destroy(opt->local_opt);
539           opt->local_opt = nlopt_copy(local_opt);
540           if (local_opt) {
541                if (!opt->local_opt) return NLOPT_OUT_OF_MEMORY;
542                nlopt_set_lower_bounds(opt->local_opt, opt->lb);
543                nlopt_set_upper_bounds(opt->local_opt, opt->ub);
544                nlopt_remove_inequality_constraints(opt->local_opt);
545                nlopt_remove_equality_constraints(opt->local_opt);
546                nlopt_set_min_objective(opt->local_opt, NULL, NULL);
547                opt->local_opt->force_stop = 0;
548           }
549           return NLOPT_SUCCESS;
550      }
551      return NLOPT_INVALID_ARGS;
552 }
553
554 /*************************************************************************/
555
556 GETSET(population, unsigned, stochastic_population)
557
558 /*************************************************************************/
559
560 nlopt_result NLOPT_STDCALL nlopt_set_initial_step1(nlopt_opt opt, double dx)
561 {
562      unsigned i;
563      if (!opt || dx == 0) return NLOPT_INVALID_ARGS;
564      if (!opt->dx && opt->n > 0) {
565           opt->dx = (double *) malloc(sizeof(double) * (opt->n));
566           if (!opt->dx) return NLOPT_OUT_OF_MEMORY;
567      }
568      for (i = 0; i < opt->n; ++i) opt->dx[i] = dx;
569      return NLOPT_SUCCESS;
570 }
571
572 nlopt_result
573 NLOPT_STDCALL nlopt_set_initial_step(nlopt_opt opt, const double *dx)
574 {
575      unsigned i;
576      if (!opt) return NLOPT_INVALID_ARGS;
577      if (!dx) {
578           free(opt->dx); opt->dx = NULL;
579           return NLOPT_SUCCESS;
580      }
581      for (i = 0; i < opt->n; ++i) if (dx[i] == 0) return NLOPT_INVALID_ARGS;
582      if (!opt->dx && nlopt_set_initial_step1(opt, 1) == NLOPT_OUT_OF_MEMORY)
583           return NLOPT_OUT_OF_MEMORY;
584      memcpy(opt->dx, dx, sizeof(double) * (opt->n));
585      return NLOPT_SUCCESS;
586 }
587
588 nlopt_result
589 NLOPT_STDCALL nlopt_get_initial_step(const nlopt_opt opt, const double *x, 
590                                      double *dx)
591 {
592      if (!opt) return NLOPT_INVALID_ARGS;
593      if (!opt->n) return NLOPT_SUCCESS;
594      if (!opt->dx) {
595           nlopt_opt o = (nlopt_opt) opt; /* discard const temporarily */
596           nlopt_result ret = nlopt_set_default_initial_step(o, x);
597           if (ret != NLOPT_SUCCESS) return ret;
598           memcpy(dx, o->dx, sizeof(double) * (opt->n));
599           free(o->dx); o->dx = NULL; /* don't save, since x-dependent */
600      }
601      else
602           memcpy(dx, opt->dx, sizeof(double) * (opt->n));
603      return NLOPT_SUCCESS;
604 }
605
606 nlopt_result
607 NLOPT_STDCALL nlopt_set_default_initial_step(nlopt_opt opt, const double *x)
608 {
609      const double *lb, *ub;
610      unsigned i;
611
612      if (!opt || !x) return NLOPT_INVALID_ARGS;
613      lb = opt->lb; ub = opt->ub;
614
615      if (!opt->dx && nlopt_set_initial_step1(opt, 1) == NLOPT_OUT_OF_MEMORY)
616           return NLOPT_OUT_OF_MEMORY;
617
618      /* crude heuristics for initial step size of nonderivative algorithms */
619      for (i = 0; i < opt->n; ++i) {
620           double step = HUGE_VAL;
621
622           if (!nlopt_isinf(ub[i]) && !nlopt_isinf(lb[i])
623               && (ub[i] - lb[i]) * 0.25 < step && ub[i] > lb[i])
624                step = (ub[i] - lb[i]) * 0.25;
625           if (!nlopt_isinf(ub[i]) 
626               && ub[i] - x[i] < step && ub[i] > x[i])
627                step = (ub[i] - x[i]) * 0.75;
628           if (!nlopt_isinf(lb[i]) 
629               && x[i] - lb[i] < step && x[i] > lb[i])
630                step = (x[i] - lb[i]) * 0.75;
631
632           if (nlopt_isinf(step)) {
633                if (!nlopt_isinf(ub[i]) 
634                    && fabs(ub[i] - x[i]) < fabs(step))
635                     step = (ub[i] - x[i]) * 1.1;
636                if (!nlopt_isinf(lb[i]) 
637                    && fabs(x[i] - lb[i]) < fabs(step))
638                     step = (x[i] - lb[i]) * 1.1;
639           }
640           if (nlopt_isinf(step) || step == 0) {
641                step = x[i];
642           }
643           if (nlopt_isinf(step) || step == 0)
644                step = 1;
645           
646           opt->dx[i] = step;
647      }
648      return NLOPT_SUCCESS;
649 }
650
651 /*************************************************************************/
652
653 void NLOPT_STDCALL nlopt_set_munge(nlopt_opt opt,
654                                    nlopt_munge munge_on_destroy,
655                                    nlopt_munge munge_on_copy) {
656      if (opt) {
657           opt->munge_on_destroy = munge_on_destroy;
658           opt->munge_on_copy = munge_on_copy;
659      }
660 }
661
662 /*************************************************************************/