chiark / gitweb /
Use trusty
[nlopt.git] / api / nlopt.3
1 .\" 
2 .\" Copyright (c) 2007 Massachusetts Institute of Technology
3 .\" 
4 .\" Copying and distribution of this file, with or without modification,
5 .\" are permitted in any medium without royalty provided the copyright
6 .\" notice and this notice are preserved.
7 .\"
8 .TH NLOPT 3  2007-08-23 "MIT" "NLopt programming manual"
9 .SH NAME
10 nlopt \- Nonlinear optimization library
11 .SH SYNOPSIS
12 .nf
13 .B #include <nlopt.h>
14 .sp
15 .BI "nlopt_opt " "opt" " = nlopt_create(" "algorithm" , " n" );
16 .BI "nlopt_set_min_objective(" "opt" , " f" , " f_data" );
17 .BI "nlopt_set_ftol_rel(" "opt" , " tol");
18 .BI "..."
19 .BI "nlopt_optimize(" "opt" , " x " , " &opt_f" );
20 .BI "nlopt_destroy(" "opt" );
21 .sp
22 The "..." indicates any number of calls to NLopt functions, below, to
23 set parameters of the optimization, constraints, and stopping
24 criteria.  Here, \fBnlopt_set_ftol_rel\fR is merely an example of a
25 possible stopping criterion.  You should link the resulting program
26 with the linker flags \-lnlopt \-lm on Unix.
27 .fi
28 .SH DESCRIPTION
29 NLopt is a library for nonlinear optimization.  It attempts to
30 minimize (or maximize) a given nonlinear objective function
31 .I f
32 of
33 .I n
34 design variables, using the specified
35 .IR algorithm ,
36 possibly subject to linear or nonlinear constraints.  The optimum
37 function value found is returned in \fIopt_f\fR (type double) with the
38 corresponding design variable values returned in the (double) array
39 .I x
40 of length
41 .IR n .
42 The input values in
43 .I x
44 should be a starting guess for the optimum.
45 .sp
46 The parameters of the optimization are controlled via the object
47 .I opt
48 of type
49 .BR nlopt_opt ,
50 which is created by the function
51 .B nlopt_create
52 and disposed of by
53 .BR nlopt_destroy .
54 By calling various functions in the NLopt library, one can specify
55 stopping criteria (e.g., a relative tolerance on the objective
56 function value is specified by
57 .BR nlopt_set_ftol_rel ),
58 upper and/or lower bounds on the design parameters 
59 .IR x ,
60 and even arbitrary nonlinear inequality and equality constraints.
61 .sp
62 By changing the parameter
63 .I algorithm
64 among several predefined constants described below, one can switch easily
65 between a variety of minimization algorithms.  Some of these algorithms
66 require the gradient (derivatives) of the function to be supplied via
67 .IR f ,
68 and other algorithms do not require derivatives.  Some of the
69 algorithms attempt to find a global optimum within the given bounds,
70 and others find only a local optimum.  Most of the algorithms only
71 handle the case where there are no nonlinear constraints.  The NLopt
72 library is a wrapper around several free/open-source minimization
73 packages, as well as some new implementations of published
74 optimization algorithms.  You could, of course, compile and call these
75 packages separately, and in some cases this will provide greater
76 flexibility than is available via NLopt.  However, depending upon the
77 specific function being optimized, the different algorithms will vary
78 in effectiveness.  The intent of NLopt is to allow you to quickly
79 switch between algorithms in order to experiment with them for your
80 problem, by providing a simple unified interface to these subroutines.
81 .SH OBJECTIVE FUNCTION
82 The objective function is specified by calling one of:
83 .sp
84 .BI "  nlopt_result nlopt_set_min_objective(nlopt_opt " "opt" , 
85 .br
86 .BI "                                       nlopt_func " "f" ,
87 .br
88 .BI "                                       void* " "f_data" );
89 .br
90 .BI "  nlopt_result nlopt_set_max_objective(nlopt_opt " "opt" , 
91 .br
92 .BI "                                       nlopt_func " "f" ,
93 .br
94 .BI "                                       void* " "f_data" );
95 .sp
96 depending on whether one wishes to minimize or maximize the objective
97 function
98 .IR f ,
99 respectively.  The function
100 .I f
101 should be of the form:
102 .sp
103 .BI "  double f(unsigned " "n" , 
104 .br
105 .BI "           const double* " "x" , 
106 .br
107 .BI "           double* " "grad" , 
108 .br
109 .BI "           void* " "f_data" );
110 .sp
111 The return value should be the value of the function at the point
112 .IR x ,
113 where
114 .I x
115 points to an array of length
116 .I n
117 of the design variables.  The dimension
118 .I n
119 is identical to the one passed to
120 .BR nlopt_create .
121 .sp
122 In addition, if the argument
123 .I grad
124 is not NULL, then
125 .I grad
126 points to an array of length
127 .I n
128 which should (upon return) be set to the gradient of the function with
129 respect to the design variables at
130 .IR x .
131 That is,
132 .IR grad[i]
133 should upon return contain the partial derivative df/dx[i],
134 for 0 <= i < n, if
135 .I grad
136 is non-NULL.
137 Not all of the optimization algorithms (below) use the gradient information:
138 for algorithms listed as "derivative-free," the 
139 .I grad
140 argument will always be NULL and need never be computed.  (For
141 algorithms that do use gradient information, however,
142 .I grad
143 may still be NULL for some calls.)
144 .sp
145 The 
146 .I f_data
147 argument is the same as the one passed to 
148 .B nlopt_set_min_objective
149 or
150 .BR nlopt_set_max_objective ,
151 and may be used to pass any additional data through to the function.
152 (That is, it may be a pointer to some caller-defined data
153 structure/type containing information your function needs, which you
154 convert from void* by a typecast.)
155 .SH BOUND CONSTRAINTS
156 Most of the algorithms in NLopt are designed for minimization of
157 functions with simple bound constraints on the inputs.  That is, the
158 input vectors x[i] are constrainted to lie in a hyperrectangle lb[i]
159 <= x[i] <= ub[i] for 0 <= i < n.  These bounds are specified by
160 passing arrays
161 .I lb
162 and
163 .I ub
164 of length
165 .I n
166 to one or both of the functions:
167 .sp
168 .BI "  nlopt_result nlopt_set_lower_bounds(nlopt_opt " "opt" , 
169 .br
170 .BI "                                      const double* " "lb" ); 
171 .br
172 .BI "  nlopt_result nlopt_set_upper_bounds(nlopt_opt " "opt" , 
173 .br
174 .BI "                                      const double* " "ub" ); 
175 .sp
176 If a lower/upper bound is not set, the default is no bound
177 (unconstrained, i.e. a bound of infinity); it is possible to have
178 lower bounds but not upper bounds or vice versa.  Alternatively, the
179 user can call one of the above functions and explicitly pass a lower
180 bound of \-HUGE_VAL and/or an upper bound of +HUGE_VAL for some design
181 variables to make them have no lower/upper bound, respectively.
182 (HUGE_VAL is the standard C constant for a floating-point infinity,
183 found in the math.h header file.)
184 .sp
185 Note, however, that some of the algorithms in NLopt, in particular
186 most of the global-optimization algorithms, do not support unconstrained
187 optimization and will return an error if you do not supply finite lower
188 and upper bounds.
189 .sp
190 For convenience, the following two functions are supplied in order to
191 set the lower/upper bounds for all design variables to a single
192 constant (so that you don't have to fill an array with a constant value):
193 .sp
194 .BI "  nlopt_result nlopt_set_lower_bounds1(nlopt_opt " "opt" , 
195 .br
196 .BI "                                       double " "lb" ); 
197 .br
198 .BI "  nlopt_result nlopt_set_upper_bounds1(nlopt_opt " "opt" , 
199 .br
200 .BI "                                       double " "ub" ); 
201 .sp
202 .SH NONLINEAR CONSTRAINTS
203 Several of the algorithms in NLopt (MMA and ORIG_DIRECT) also support
204 arbitrary nonlinear inequality constraints, and some also allow
205 nonlinear equality constraints (COBYLA, SLSQP, ISRES, and AUGLAG).
206 For these algorithms, you can specify as many nonlinear constraints as
207 you wish by calling the following functions multiple times.
208 .sp
209 In particular, a nonlinear inequality constraint of the form 
210 \fIfc\fR(\fIx\fR) <= 0, where the function
211 .I fc
212 is of the same form as the objective function described above,
213 can be specified by calling:
214 .sp
215 .BI "  nlopt_result nlopt_add_inequality_constraint(nlopt_opt " "opt" , 
216 .br
217 .BI "                                               nlopt_func " "fc" ,
218 .br
219 .BI "                                               void* " "fc_data" ,
220 .br
221 .BI "                                               double " "tol" );
222 .sp
223 Just as for the objective function, 
224 .I fc_data
225 is a pointer to arbitrary user data that will be passed through to the
226 .I fc
227 function whenever it is called.  The parameter
228 .I tol
229 is a tolerance that is used for the purpose of stopping criteria only:
230 a point
231 .I x
232 is considered feasible for judging whether to stop the optimization if
233 \fIfc\fR(\fIx\fR) <= \fItol\fR.  A tolerance of zero means that NLopt
234 will try not to consider any \fIx\fR to be converged unless
235 .I fc
236 is strictly non-positive; generally, at least a small positive tolerance is
237 advisable to reduce sensitivity to rounding errors.
238
239 A nonlinear equality constraint of the form 
240 \fIh\fR(\fIx\fR) = 0, where the function
241 .I h
242 is of the same form as the objective function described above,
243 can be specified by calling:
244 .sp
245 .BI "  nlopt_result nlopt_add_equality_constraint(nlopt_opt " "opt" , 
246 .br
247 .BI "                                             nlopt_func " "h" ,
248 .br
249 .BI "                                             void* " "h_data" ,
250 .br
251 .BI "                                             double " "tol" );
252 .sp
253 Just as for the objective function, 
254 .I h_data
255 is a pointer to arbitrary user data that will be passed through to the
256 .I h
257 function whenever it is called.  The parameter
258 .I tol
259 is a tolerance that is used for the purpose of stopping criteria only:
260 a point
261 .I x
262 is considered feasible for judging whether to stop the optimization if
263 |\fIh\fR(\fIx\fR)| <= \fItol\fR.  For equality constraints, a small
264 positive tolerance is strongly advised in order to allow NLopt to 
265 converge even if the equality constraint is slightly nonzero.
266 .sp
267 (For any algorithm listed as "derivative-free" below, the
268 .I grad
269 argument to \fIfc\fR or \fIh\fR will always be NULL and need never be
270 computed.)
271 .sp
272 To remove all of the inequality and/or equality constraints from 
273 a given problem \fIopt\fR, you can call the following functions:
274 .sp
275 .BI "  nlopt_result nlopt_remove_inequality_constraints(nlopt_opt " "opt" );
276 .br
277 .BI "  nlopt_result nlopt_remove_equality_constraints(nlopt_opt " "opt" );
278 .SH ALGORITHMS
279 The 
280 .I algorithm
281 parameter specifies the optimization algorithm (for more detail on
282 these, see the README files in the source-code subdirectories), and
283 can take on any of the following constant values.
284 .sp
285 Constants with
286 .B _G{N,D}_
287 in their names
288 refer to global optimization methods, whereas
289 .B _L{N,D}_
290 refers to local optimization methods (that try to find a local optimum
291 starting from the starting guess
292 .IR x ).
293 Constants with
294 .B _{G,L}N_
295 refer to non-gradient (derivative-free) algorithms that do not require the
296 objective function to supply a gradient, whereas
297 .B _{G,L}D_
298 refers to derivative-based algorithms that require the objective
299 function to supply a gradient.  (Especially for local optimization,
300 derivative-based algorithms are generally superior to derivative-free
301 ones: the gradient is good to have 
302 .I if 
303 you can compute it cheaply, e.g. via an adjoint method.)
304 .sp
305 The algorithm specified for a given problem
306 .I opt
307 is returned by the function:
308 .sp
309 .BI "  nlopt_algorithm nlopt_get_algorithm(nlopt_opt " "opt" );
310 .sp
311 The available algorithms are:
312 .TP 
313 .B NLOPT_GN_DIRECT_L
314 Perform a global (G) derivative-free (N) optimization using the
315 DIRECT-L search algorithm by Jones et al. as modified by Gablonsky et
316 al. to be more weighted towards local search.  Does not support
317 unconstrainted optimization.  There are also several other variants of
318 the DIRECT algorithm that are supported:
319 .BR NLOPT_GN_DIRECT ,
320 which is the original DIRECT algorithm;
321 .BR NLOPT_GN_DIRECT_L_RAND ,
322 a slightly randomized version of DIRECT-L that may be better in
323 high-dimensional search spaces;
324 .BR NLOPT_GN_DIRECT_NOSCAL ,
325 .BR NLOPT_GN_DIRECT_L_NOSCAL ,
326 and
327 .BR NLOPT_GN_DIRECT_L_RAND_NOSCAL ,
328 which are versions of DIRECT where the dimensions are not rescaled to
329 a unit hypercube (which means that dimensions with larger bounds are
330 given more weight).
331 .TP 
332 .B NLOPT_GN_ORIG_DIRECT_L
333 A global (G) derivative-free optimization using the DIRECT-L algorithm
334 as above, along with
335 .B NLOPT_GN_ORIG_DIRECT
336 which is the original DIRECT algorithm.  Unlike 
337 .B NLOPT_GN_DIRECT_L 
338 above, these two algorithms refer to code based on the original
339 Fortran code of Gablonsky et al., which has some hard-coded
340 limitations on the number of subdivisions etc. and does not support
341 all of the NLopt stopping criteria, but on the other hand it supports
342 arbitrary nonlinear inequality constraints.
343 .TP 
344 .B NLOPT_GD_STOGO
345 Global (G) optimization using the StoGO algorithm by Madsen et al.  StoGO
346 exploits gradient information (D) (which must be supplied by the
347 objective) for its local searches, and performs the global search by a
348 branch-and-bound technique.  Only bound-constrained optimization
349 is supported.  There is also another variant of this algorithm,
350 .BR NLOPT_GD_STOGO_RAND ,
351 which is a randomized version of the StoGO search scheme.  The StoGO
352 algorithms are only available if NLopt is compiled with C++ code
353 enabled, and should be linked via \-lnlopt_cxx instead of \-lnlopt (via
354 a C++ compiler, in order to link the C++ standard libraries).
355 .TP 
356 .B NLOPT_LN_NELDERMEAD
357 Perform a local (L) derivative-free (N) optimization, starting at
358 .IR x ,
359 using the Nelder-Mead simplex algorithm, modified to support bound
360 constraints.  Nelder-Mead, while popular, is known to occasionally
361 fail to converge for some objective functions, so it should be used
362 with caution.  Anecdotal evidence, on the other hand, suggests that it
363 works fairly well for some cases that are hard to handle otherwise,
364 e.g. noisy/discontinuous objectives.  See also
365 .B NLOPT_LN_SBPLX
366 below.
367 .TP 
368 .B NLOPT_LN_SBPLX
369 Perform a local (L) derivative-free (N) optimization, starting at
370 .IR x ,
371 using an algorithm based on the Subplex algorithm of Rowan et al.,
372 which is an improved variant of Nelder-Mead (above).  Our
373 implementation does not use Rowan's original code, and has some minor
374 modifications such as explicit support for bound constraints.  (Like
375 Nelder-Mead, Subplex often works well in practice, even for
376 noisy/discontinuous objectives, but there is no rigorous guarantee that it
377 will converge.)
378 .TP
379 .B NLOPT_LN_PRAXIS
380 Local (L) derivative-free (N) optimization using the principal-axis
381 method, based on code by Richard Brent.  Designed for unconstrained
382 optimization, although bound constraints are supported too (via the
383 inefficient method of returning +Inf when the constraints are violated).
384 .TP
385 .B NLOPT_LD_LBFGS
386 Local (L) gradient-based (D) optimization using the limited-memory BFGS
387 (L-BFGS) algorithm.  (The objective function must supply the
388 gradient.)  Unconstrained optimization is supported in addition to
389 simple bound constraints (see above).  Based on an implementation by
390 Luksan et al.
391 .TP
392 .B NLOPT_LD_VAR2
393 Local (L) gradient-based (D) optimization using a shifted limited-memory
394 variable-metric method based on code by Luksan et al., supporting both
395 unconstrained and bound-constrained optimization.  
396 .B NLOPT_LD_VAR2
397 uses a rank-2 method, while 
398 .B .B NLOPT_LD_VAR1
399 is another variant using a rank-1 method.
400 .TP
401 .B NLOPT_LD_TNEWTON_PRECOND_RESTART
402 Local (L) gradient-based (D) optimization using an
403 LBFGS-preconditioned truncated Newton method with steepest-descent
404 restarting, based on code by Luksan et al., supporting both
405 unconstrained and bound-constrained optimization.  There are several
406 other variants of this algorithm:
407 .B NLOPT_LD_TNEWTON_PRECOND 
408 (same without restarting), 
409 .B NLOPT_LD_TNEWTON_RESTART
410 (same without preconditioning), and
411 .B NLOPT_LD_TNEWTON
412 (same without restarting or preconditioning).
413 .TP
414 .B NLOPT_GN_CRS2_LM
415 Global (G) derivative-free (N) optimization using the controlled random
416 search (CRS2) algorithm of Price, with the "local mutation" (LM)
417 modification suggested by Kaelo and Ali.
418 .TP
419 .B NLOPT_GN_ISRES
420 Global (G) derivative-free (N) optimization using a genetic algorithm
421 (mutation and differential evolution), using a stochastic ranking to
422 handle nonlinear inequality and equality constraints as suggested by
423 Runarsson and Yao.
424 .TP
425 \fBNLOPT_G_MLSL_LDS\fR, \fBNLOPT_G_MLSL\fR 
426 Global (G) optimization using the multi-level single-linkage (MLSL)
427 algorithm with a low-discrepancy sequence (LDS) or pseudorandom
428 numbers, respectively.  This algorithm executes a low-discrepancy
429 or pseudorandom sequence of local searches, with a clustering
430 heuristic to avoid multiple local searches for the same local optimum.
431 The local search algorithm must be specified, along with termination
432 criteria/tolerances for the local searches, by
433 \fInlopt_set_local_optimizer\fR.  (This subsidiary algorithm can be
434 with or without derivatives, and determines whether the objective
435 function needs gradients.)
436 .TP
437 \fBNLOPT_LD_MMA\fR, \fBNLOPT_LD_CCSAQ\fR 
438 Local (L) gradient-based (D) optimization using the method of moving
439 asymptotes (MMA), or rather a refined version of the algorithm as
440 published by Svanberg (2002).  (NLopt uses an independent
441 free-software/open-source implementation of Svanberg's algorithm.) CCSAQ
442 is a related algorithm from Svanberg's paper which uses a local quadratic
443 approximation rather than the more-complicated MMA model; the two usually
444 have similar convergence rates.
445 The
446 .B NLOPT_LD_MMA
447 algorithm supports both bound-constrained and unconstrained
448 optimization, and also supports an arbitrary number (\fIm\fR) of
449 nonlinear inequality (not equality) constraints as described above.
450 .TP
451 .B NLOPT_LD_SLSQP
452 Local (L) gradient-based (D) optimization using sequential quadratic
453 programming and BFGS updates, supporting arbitrary nonlinear
454 inequality and equality constraints, based on the code by Dieter Kraft
455 (1988) adapted for use by the SciPy project.  Note that this algorithm
456 uses dense-matrix methods requiring O(\fIn\fR^2) storage and
457 O(\fIn\fR^3) time, making it less practical for problems involving
458 more than a few thousand parameters.
459 .TP
460 .B NLOPT_LN_COBYLA
461 Local (L) derivative-free (N) optimization using the COBYLA algorithm
462 of Powell (Constrained Optimization BY Linear Approximations).
463 The
464 .B NLOPT_LN_COBYLA
465 algorithm supports both bound-constrained and unconstrained
466 optimization, and also supports an arbitrary number (\fIm\fR) of
467 nonlinear inequality/equality constraints as described above.
468 .TP
469 .B NLOPT_LN_NEWUOA
470 Local (L) derivative-free (N) optimization using a variant of the
471 NEWUOA algorithm of Powell, based on successive quadratic
472 approximations of the objective function. We have modified the
473 algorithm to support bound constraints.  The original NEWUOA algorithm
474 is also available, as
475 .BR NLOPT_LN_NEWUOA ,
476 but this algorithm ignores the bound constraints
477 .I lb
478 and 
479 .IR ub ,
480 and so it should only be used for unconstrained problems.  Mostly
481 superseded by BOBYQA.
482 .TP
483 .B NLOPT_LN_BOBYQA
484 Local (L) derivative-free (N) optimization using the BOBYQA algorithm
485 of Powell, based on successive quadratic approximations of the
486 objective function, supporting bound constraints.
487 .TP
488 .B NLOPT_AUGLAG
489 Optimize an objective with nonlinear inequality/equality constraints
490 via an unconstrained (or bound-constrained) optimization algorithm,
491 using a gradually increasing "augmented Lagrangian" penalty for
492 violated constraints.  Requires you to specify another optimization
493 algorithm for optimizing the objective+penalty function, using
494 \fInlopt_set_local_optimizer\fR.  (This subsidiary algorithm can be
495 global or local and with or without derivatives, but you must specify
496 its own termination criteria.)  A variant, \fBNLOPT_AUGLAG_EQ\fR, only
497 uses the penalty approach for equality constraints, while inequality
498 constraints are handled directly by the subsidiary algorithm (restricting
499 the choice of subsidiary algorithms to those that can handle inequality
500 constraints).
501 .SH STOPPING CRITERIA
502 Multiple stopping criteria for the optimization are supported, as
503 specified by the functions to modify a given optimization problem
504 .BR opt .
505 The optimization halts whenever any one of these criteria is
506 satisfied.  In some cases, the precise interpretation of the stopping
507 criterion depends on the optimization algorithm above (although we
508 have tried to make them as consistent as reasonably possible), and
509 some algorithms do not support all of the stopping criteria.
510 .sp
511 Important: you do not need to use all of the stopping criteria!  In most
512 cases, you only need one or two, and can omit the remainder (all criteria
513 are disabled by default).
514 .TP
515 .BI "nlopt_result nlopt_set_stopval(nlopt_opt " "opt" ,
516 .br
517 .BI "                        double " stopval );
518 .sp
519 Stop when an objective value of at least
520 .I stopval
521 is found: stop minimizing when a value <= \fIstopval\fR is found, or
522 stop maximizing when a value >= \fIstopval\fR is found.  (Setting
523 \fIstopval\fR to \-HUGE_VAL for minimizing or +HUGE_VAL for maximizing
524 disables this stopping criterion.)
525 .TP
526 .BI "nlopt_result nlopt_set_ftol_rel(nlopt_opt " "opt" ,
527 .br
528 .BI "                         double " tol );
529 .sp
530 Set relative tolerance on function value: stop when an optimization step
531 (or an estimate of the optimum) changes the function value by less
532 than
533 .I tol
534 multiplied by the absolute value of the function value.  (If there is any chance that your optimum function value is close to zero, you might want to set an absolute tolerance with
535 .B nlopt_set_ftol_abs
536 as well.)  Criterion is disabled if \fItol\fR is non-positive.
537 .TP
538 .BI "nlopt_result nlopt_set_ftol_abs(nlopt_opt " "opt" ,
539 .br
540 .BI "                         double " tol );
541 .sp
542 Set absolute tolerance on function value: stop when an optimization step
543 (or an estimate of the optimum) changes the function value by less
544 than
545 .IR tol .
546 Criterion is disabled if \fItol\fR is non-positive.
547 .TP
548 .BI "nlopt_result nlopt_set_xtol_rel(nlopt_opt " "opt" ,
549 .br
550 .BI "                         double " tol );
551 .sp
552 Set relative tolerance on design variables: stop when an optimization step
553 (or an estimate of the optimum) changes every design variable by less
554 than
555 .I tol
556 multiplied by the absolute value of the design variable.  (If there is
557 any chance that an optimal design variable is close to zero, you
558 might want to set an absolute tolerance with
559 .B nlopt_set_xtol_abs
560 as well.)  Criterion is disabled if \fItol\fR is non-positive.
561 .TP
562 .BI "nlopt_result nlopt_set_xtol_abs(nlopt_opt " "opt" ,
563 .br
564 .BI "                         const double* " tol );
565 .sp
566 Set absolute tolerances on design variables.  \fItol\fR is a pointer
567 to an array of length
568 .I
569 n giving the tolerances: stop when an
570 optimization step (or an estimate of the optimum) changes every design
571 variable
572 .IR x [i]
573 by less than
574 .IR tol [i].
575 .sp
576 For convenience, the following function may be used to set the absolute tolerances in all \fIn\fR design variables to the same value:
577 .sp
578 .BI "  nlopt_result nlopt_set_xtol_abs1(nlopt_opt " "opt" ,
579 .br
580 .BI "                                   double " tol );
581 .sp
582 Criterion is disabled if \fItol\fR is non-positive.
583 .TP
584 .BI "nlopt_result nlopt_set_maxeval(nlopt_opt " "opt" ,
585 .br
586 .BI "                        int " maxeval );
587 .sp
588 Stop when the number of function evaluations exceeds
589 .IR maxeval .
590 (This is not a strict maximum: the number of function evaluations may
591 exceed
592 .I maxeval 
593 slightly, depending upon the algorithm.)  Criterion is disabled
594 if \fImaxeval\fR is non-positive.
595 .TP
596 .BI "nlopt_result nlopt_set_maxtime(nlopt_opt " "opt" ,
597 .br
598 .BI "                        double " maxtime );
599 .sp
600 Stop when the optimization time (in seconds) exceeds
601 .IR maxtime .
602 (This is not a strict maximum: the time may
603 exceed
604 .I maxtime
605 slightly, depending upon the algorithm and on how slow your function
606 evaluation is.)  Criterion is disabled if \fImaxtime\fR is non-positive.
607 .SH RETURN VALUE
608 Most of the NLopt functions return an enumerated constant
609 of type
610 .BR nlopt_result ,
611 which takes on one of the following values:
612 .SS Successful termination (positive return values):
613 .TP
614 .B NLOPT_SUCCESS
615 Generic success return value.
616 .TP
617 .B NLOPT_STOPVAL_REACHED
618 Optimization stopped because
619 .I stopval
620 (above) was reached.
621 .TP
622 .B NLOPT_FTOL_REACHED
623 Optimization stopped because
624 .I ftol_rel
625 or
626 .I ftol_abs
627 (above) was reached.
628 .TP
629 .B NLOPT_XTOL_REACHED
630 Optimization stopped because
631 .I xtol_rel
632 or
633 .I xtol_abs
634 (above) was reached.
635 .TP
636 .B NLOPT_MAXEVAL_REACHED
637 Optimization stopped because
638 .I maxeval
639 (above) was reached.
640 .TP
641 .B NLOPT_MAXTIME_REACHED
642 Optimization stopped because
643 .I maxtime
644 (above) was reached.
645 .SS Error codes (negative return values):
646 .TP
647 .B NLOPT_FAILURE
648 Generic failure code.
649 .TP
650 .B NLOPT_INVALID_ARGS
651 Invalid arguments (e.g. lower bounds are bigger than upper bounds, an
652 unknown algorithm was specified, etcetera).
653 .TP
654 .B NLOPT_OUT_OF_MEMORY
655 Ran out of memory.
656 .TP
657 .B NLOPT_ROUNDOFF_LIMITED
658 Halted because roundoff errors limited progress.
659 .TP
660 .B NLOPT_FORCED_STOP
661 Halted because the user called \fBnlopt_force_stop\fR(\fIopt\fR) on
662 the optimization's \fBnlopt_opt\fR object \fIopt\fR from the user's
663 objective function.
664 .SH LOCAL OPTIMIZER
665 Some of the algorithms, especially MLSL and AUGLAG, use a different
666 optimization algorithm as a subroutine, typically for local
667 optimization.  You can change the local search algorithm and its
668 tolerances by calling:
669 .sp
670 .BI "  nlopt_result nlopt_set_local_optimizer(nlopt_opt " "opt" , 
671 .br
672 .BI "                                         const nlopt_opt " "local_opt" );
673 .sp
674 Here, \fIlocal_opt\fR is another \fBnlopt_opt\fR object whose
675 parameters are used to determine the local search algorithm and
676 stopping criteria.  (The objective function, bounds, and
677 nonlinear-constraint parameters of \fIlocal_opt\fR are ignored.)  The
678 dimension \fIn\fR of \fIlocal_opt\fR must match that of \fIopt\fR.
679 .sp
680 This function makes a copy of the \fIlocal_opt\fR object, so you can
681 freely destroy your original \fIlocal_opt\fR afterwards.
682 .SH INITIAL STEP SIZE
683 For derivative-free local-optimization algorithms, the optimizer must
684 somehow decide on some initial step size to perturb \fIx\fR by when it
685 begins the optimization.  This step size should be big enough that the
686 value of the objective changes significantly, but not too big if you
687 want to find the local optimum nearest to \fIx\fR.  By default, NLopt
688 chooses this initial step size heuristically from the bounds,
689 tolerances, and other information, but this may not always be the best
690 choice.
691 .sp
692 You can modify the initial step size by calling:
693 .sp
694 .BI "  nlopt_result nlopt_set_initial_step(nlopt_opt " "opt" , 
695 .br
696 .BI "                                      const double* " "dx" );
697 .sp
698 Here, \fIdx\fR is an array of length \fIn\fR containing the (nonzero)
699 initial step size for each component of the design parameters \fIx\fR.
700 For convenience, if you want to set the step sizes in every direction
701 to be the same value, you can instead call:
702 .sp
703 .BI "  nlopt_result nlopt_set_initial_step1(nlopt_opt " "opt" , 
704 .br
705 .BI "                                       double " "dx" );
706 .SH STOCHASTIC POPULATION
707 Several of the stochastic search algorithms (e.g., CRS, MLSL, and
708 ISRES) start by generating some initial "population" of random points
709 \fIx\fR.  By default, this initial population size is chosen
710 heuristically in some algorithm-specific way, but the initial
711 population can by changed by calling:
712 .sp
713 .BI "  nlopt_result nlopt_set_population(nlopt_opt " "opt" , 
714 .br
715 .BI "                                    unsigned " "pop" );
716 .sp
717 (A \fIpop\fR of zero implies that the heuristic default will be used.)
718 .SH PSEUDORANDOM NUMBERS
719 For stochastic optimization algorithms, we use pseudorandom numbers generated
720 by the Mersenne Twister algorithm, based on code from Makoto Matsumoto.
721 By default, the seed for the random numbers is generated from the system
722 time, so that they will be different each time you run the program.  If
723 you want to use deterministic random numbers, you can set the seed by
724 calling:
725 .sp
726 .BI "            void nlopt_srand(unsigned long " "seed" );
727 .sp
728 Some of the algorithms also support using low-discrepancy sequences (LDS),
729 sometimes known as quasi-random numbers.  NLopt uses the Sobol LDS, which
730 is implemented for up to 1111 dimensions.
731 .SH AUTHORS
732 Written by Steven G. Johnson.
733 .PP
734 Copyright (c) 2007-2014 Massachusetts Institute of Technology.
735 .SH "SEE ALSO"
736 nlopt_minimize(3)