5 In this tutorial, we illustrate the usage of NLopt in various languages via one or two trivial examples.
7 Example nonlinearly constrained problem
8 ---------------------------------------
11 ![right|thumb|400px|Feasible region for a simple example optimization problem with two nonlinear (cubic) constraints.](images/NLopt-example-constraints.png)
13 As a first example, we'll look at the following simple nonlinearly constrained minimization problem:
15 $$\min_{\mathbf{x}\in\mathbb{R}^2} \sqrt{x_2}$$
18 subject to $x_2 \geq 0$, $x_2 \geq (a_1 x_1 + b_1)^3$, and $x_2 \geq (a_2 x_1 + b_2)^3$
20 for parameters *a*<sub>1</sub>=2, *b*<sub>1</sub>=0, *a*<sub>2</sub>=-1, *b*<sub>2</sub>=1.
22 The feasible region defined by these constraints is plotted at right: *x*<sub>2</sub> is constrained to lie above the maximum of two cubics, and the optimum point is located at the intersection (1/3, 8/27) where the objective function takes on the value $\sqrt{8/27} \approx 0.5443310539518\ldots$.
24 (This problem is especially trivial, because by formulating it in terms of the cube root of *x*<sub>2</sub> you can turn it into a linear-programming problem, but we won't do that here.)
26 In principle, we don't need the bound constraint *x*<sub>2</sub>≥0, since the nonlinear constraints already imply a positive-*x*<sub>2</sub> feasible region. However, NLopt doesn't guarantee that, on the way to finding the optimum, it won't violate the nonlinear constraints at some intermediate steps, while it *does* guarantee that all intermediate steps will satisfy the bound constraints. So, we will explicitly impose *x*<sub>2</sub>≥0 in order to ensure that the √*x*<sub>2</sub> in our objective is real.
28 **Note:** The objective function here is not differentiable at *x*<sub>2</sub>=0. This doesn't cause problems in the examples below, but may cause problems with some other algorithms if they try to evaluate the gradient at *x*<sub>2</sub>=0 (e.g. I've seen it cause AUGLAG with a gradient-based solver to fail). To prevent this, you might want to use a small nonzero lower bound instead, e.g. *x*<sub>2</sub>≥10<sup>−6</sup>.
33 To implement the above example in C or C++, we would first do:
38 to include the NLopt header file as well as the standard math header file (needed for things like the `sqrt` function and the `HUGE_VAL` constant), then we would define our objective function as:
41 double myfunc(unsigned n, const double *x, double *grad, void *my_func_data)
45 grad[1] = 0.5 / sqrt(x[1]);
52 There are several things to notice here. First, since this is C, our indices are zero-based, so we have `x[0]` and `x[1]` instead of *x*<sub>1</sub> and *x*<sub>2</sub>. The return value of our function is the objective $\sqrt{x_2}$. Also, if the parameter `grad` is not `NULL`, then we set `grad[0]` and `grad[1]` to the partial derivatives of our objective with respect to `x[0]` and `x[1]`. The gradient is only needed for [gradient-based algorithms](NLopt_Introduction#Gradient-based_versus_derivative-free_algorithms.md); if you use a derivative-free optimization algorithm, `grad` will always be `NULL` and you need never compute any derivatives. Finally, we have an extra parameter `my_func_data` that can be used to pass additional data to `myfunc`, but no additional data is needed here so that parameter is unused.
54 For the constraints, on the other hand, we *will* have additional data. Each constraint is parameterized by two numbers *a* and *b*, so we will declare a data structure to hold this information:
63 Then, we implement our constraint function as follows.
66 double myconstraint(unsigned n, const double *x, double *grad, void *data)
68 my_constraint_data *d = (my_constraint_data *) data;
69 double a = d->a, b = d->b;
71 grad[0] = 3 * a * (a*x[0] + b) * (a*x[0] + b);
74 return ((a*x[0] + b) * (a*x[0] + b) * (a*x[0] + b) - x[1]);
79 The form of the constraint function is the same as that of the objective function. Here, the `data` parameter will actually be a pointer to `my_constraint_data` (because this is the type that we will pass to `nlopt_minimize_constrained` below), so we use a typecast to get the constraint data. NLopt always expects constraints to be of the form `myconstraint`(**x**) ≤ 0, so we implement the constraint *x*<sub>2</sub> ≥ (*a* *x*<sub>1</sub> + *b*)<sup>3</sup> as the function (*a* *x*<sub>1</sub> + *b*)<sup>3</sup> − *x*<sub>2</sub>. Again, we only compute the gradient if `grad` is non-`NULL`, which will never occur if we use a derivative-free optimization algorithm.
81 Now, to specify this optimization problem, we create an "object" of type `nlopt_opt` (an opaque pointer type) and set its various parameters:
84 double lb[2] = { -HUGE_VAL, 0 }; /* lower bounds */
90 opt = nlopt_create(NLOPT_LD_MMA, 2); /* algorithm and dimensionality */
91 nlopt_set_lower_bounds(opt, lb);
92 nlopt_set_min_objective(opt, myfunc, NULL);
96 Note that we do not need to set an upper bound (`nlopt_set_upper_bounds`), since we are happy with the default upper bounds (+∞). To add the two inequality constraints, we do:
99 my_constraint_data data[2] = { {2,0}, {-1,1} };
104 nlopt_add_inequality_constraint(opt, myconstraint, &data[0], 1e-8);
105 nlopt_add_inequality_constraint(opt, myconstraint, &data[1], 1e-8);
109 Here, the `1e-8` is an optional tolerance for the constraint: for purposes of convergence testing, a point will be considered feasible if the constraint is violated (is positive) by that tolerance (10<sup>−8</sup>). A nonzero tolerance is a good idea for many algorithms lest tiny errors prevent convergence. Speaking of convergence tests, we should also set one or more stopping criteria, e.g. a relative tolerance on the optimization parameters **x**:
112 nlopt_set_xtol_rel(opt, 1e-4);
116 There are many more possible parameters that you can set to control the optimization, which are described in detail by the [reference manual](NLopt_Reference.md), but these are enough for our example here (any unspecified parameters are set to innocuous defaults). At this point, we can call nlopt_optimize to actually perform the optimization, starting with some initial guess:
119 double x[2] = { 1.234, 5.678 }; /* `*`some` `initial` `guess`*` */
120 double minf; /* `*`the` `minimum` `objective` `value,` `upon` `return`*` */
121 if (nlopt_optimize(opt, x, &minf) < 0) {
122 printf("nlopt failed!\n");
125 printf("found minimum at f(%g,%g) = %0.10g\n", x[0], x[1], minf);
130 `nlopt_optimize` will return a negative result code on failure, but this usually only happens if you pass invalid parameters, it runs out of memory, or something like that. (Actually, most of the other NLopt functions also return an error code that you can check if you are paranoid.) (However, if it returns the failure code `NLOPT_ROUNDOFF_LIMITED`, indicating a breakdown due to roundoff errors, the minimum found may still be useful and you may want to still use it.) Otherwise, we print out the minimum function value and the corresponding parameters **x**.
132 Finally, we should call nlopt_destroy to dispose of the `nlopt_opt` object when we are done with it:
139 Assuming we save this in a file tutorial.c, we would compile and link (on Unix) with:
142 cc tutorial.c -o tutorial -lnlopt -lm
146 The result of running the program should then be something like:
149 found minimum at f(0.333334,0.296296) = 0.544330847
153 That is, it found the correct parameters to about 5 significant digits and the correct minimum function value to about 6 significant digits. (This is better than we specified; this often occurs because the local optimization routines usually try to be conservative in estimating the error.)
155 ### Number of evaluations
157 Let's modify our program to print out the number of function evaluations that were required to obtain this result. First, we'll change our objective function to:
161 double myfunc(int n, const double *x, double *grad, void *my_func_data)
166 grad[1] = 0.5 / sqrt(x[1]);
173 using a global variable `count` that is incremented for each function evaluation. (We could also pass a pointer to a counter variable as `my_func_data`, if we wanted to avoid global variables.) Then, adding a `printf`:
176 printf("found minimum after %d evaluations\n", count);
183 found minimum after 11 evaluations
184 found minimum at f(0.333334,0.296296) = 0.544330847
188 For such a simple problem, a gradient-based local optimization algorithm like MMA can converge very quickly!
190 ### Switching to a derivative-free algorithm
192 We can also try a derivative-free algorithm. Looking at the [NLopt Algorithms](NLopt_Algorithms.md) list, another algorithm in NLopt that handles nonlinear constraints is COBYLA, which is derivative-free. To use it, we just change `NLOPT_LD_MMA` ("LD" means local optimization, derivative/gradient-based) into `NLOPT_LN_COBYLA` ("LN" means local optimization, no derivatives), and obtain:
195 found minimum after 31 evaluations
196 found minimum at f(0.333329,0.2962) = 0.544242301
200 In such a low-dimensional problem, derivative-free algorithms usually work quite well—in this case, it only triples the number of function evaluations. However, the comparison is not perfect because, for the same relative **x** tolerance of 10<sup>−4</sup>, COBYLA is a bit less conservative and only finds the solution to 3 significant digits.
202 To do a fairer comparison of the two algorithms, we could set the **x** tolerance to zero and ask how many function evaluations each one requires to get the correct answer to three decimal places. We can specify this by using the `stopval` termination criterion, which allows us to halt the process as soon as a feasible point attains an objective function value less than `stopval`. In this case, we would set `stopval` to $\sqrt(8/27)+10^{-3}$, replacing `nlopt_set_xtol_rel` with the statement:
205 nlopt_set_stopval(opt, sqrt(8./27.)+1e-3);
209 corresponding to the last line of arguments to `nlopt_minimize_constrained` being `sqrt(8./27.)+1e-3,` `0.0,` `0.0,` `0.0,` `NULL,` `0,` `0.0`. If we do this, we find that COBYLA requires 25 evaluations while MMA requires 10.
211 The advantage of gradient-based algorithms over derivative-free algorithms typically grows for higher-dimensional problems. On the other hand, derivative-free algorithms are much easier to use because you don't need to worry about how to compute the gradient (which might be tricky if the function is very complicated).
216 Although it is perfectly possible to use the C interface from C++, many C++ programmers will find it more natural to use real C++ objects instead of opaque `nlopt_opt` pointers, `std::vector`<double> instead of arrays, and exceptions instead of error codes. NLopt provides a C++ header file `nlopt.hpp` that you can use for this purpose, which simply wraps a C++ object interface around the C interface above.
218 `#include `<nlopt.hpp>
220 The equivalent of the above example would then be:
223 nlopt::opt opt(nlopt::LD_MMA, 2);
224 std::vector`<double>` lb(2);
225 lb[0] = -HUGE_VAL; lb[1] = 0;
226 opt.set_lower_bounds(lb);
227 opt.set_min_objective(myfunc, NULL);
228 my_constraint_data data[2] = { {2,0}, {-1,1} };
229 opt.add_inequality_constraint(myconstraint, &data[0], 1e-8);
230 opt.add_inequality_constraint(myconstraint, &data[1], 1e-8);
231 opt.set_xtol_rel(1e-4);
232 std::vector`<double>` x(2);
233 x[0] = 1.234; x[1] = 5.678;
235 nlopt::result result = opt.optimize(x, minf);
239 There is no need to deallocate the `opt` object; its destructor will do that for you once it goes out of scope. Also, there is no longer any need to check for error codes; the NLopt C++ functions will throw exceptions if there is an error, which you can `catch` normally.
241 Here, we are using the same objective and constraint functions as in C, taking `double*` array arguments. Alternatively, you can define objective and constraint functions to take `std::vector`<double> arguments if you prefer. (Using `std::vector`<double> in the objective/constraint imposes a slight overhead because NLopt must copy the `double*` data to a `std::vector`<double>, but this overhead is unlikely to be significant in most real applications.) That is, you would do:
244 double myvfunc(const std::vector`<double>` &x, std::vector`<double>` &grad, void *my_func_data)
248 grad[1] = 0.5 / sqrt(x[1]);
256 double myvconstraint(const std::vector`<double>` &x, std::vector`<double>` &grad, void *data)
258 my_constraint_data *d = reinterpret_cast`<my_constraint_data*>`(data);
259 double a = d->a, b = d->b;
261 grad[0] = 3 * a * (a*x[0] + b) * (a*x[0] + b);
264 return ((a*x[0] + b) * (a*x[0] + b) * (a*x[0] + b) - x[1]);
269 Notice that, instead of checking whether `grad` is `NULL`, we check whether it is empty. (The vector arguments, if non-empty, are guaranteed to be of the same size as the dimension of the problem that you specified.) We then specify these in the same way as before:
272 opt.set_min_objective(myvfunc, NULL);
273 opt.add_inequality_constraint(myvconstraint, &data[0], 1e-8);
274 opt.add_inequality_constraint(myvconstraint, &data[1], 1e-8);
278 Note that the data pointers passed to these functions must remain valid (or rather, what they point to must remain valid) until you are done with `opt`. (It might have been nicer to use `shared_ptr`, but I don't like to rely on bleeding-edge language features.)
280 Instead of passing a separate data pointer, some users may wish to define a C++ [function object](https://en.wikipedia.org/wiki/Function_object) class that contains all of the data needed by their function, with an overloaded `operator()` method to implement the function call. You can easily do this with a two-line helper function. If your function class is MyFunction, then you could define a static member function:
283 static double wrap(const std::vector`<double>` &x, std::vector`<double>` &grad, void *data) {
284 return (*reinterpret_cast`<MyFunction*>`(data))(x, grad); }
288 which you would then use e.g. by `opt.set_min_objective(MyFunction::wrap,` `&some_MyFunction)`. Again, you have to make sure that `some_MyFunction` does not go out of scope before you are done calling `nlopt::opt::optimize`.
290 To link your program, just link to the C NLopt library (`-lnlopt` `-lm` on Unix).
292 Example in Matlab or GNU Octave
293 -------------------------------
295 To implement this objective function in Matlab (or GNU Octave), we would write a file myfunc.m that looks like:
298 function [val, gradient] = myfunc(x)
301 gradient = [0, 0.5 / val];
306 Notice that we check the Matlab builtin variable `nargout` (the number of output arguments) to decide whether to compute the gradient. If we use a derivative-free optimization algorithm below, then `nargout` will always be 1 and the gradient need never be computed.
308 Our constraint function looks similar, except that it is parameterized by the coefficients *a* and *b*. We can just add these on as extra parameters, in a file `myconstraint.m`:
311 function [val, gradient] = myconstraint(x,a,b)
312 val = (a*x(1) + b)^3 - x(2);
314 gradient = [3*a*(a*x(1) + b)^2, -1];
319 The equivalent of the `nlopt_opt` is just a [structure](http://www.mathworks.com/access/helpdesk/help/techdoc/matlab_prog/f2-88951.html), with fields corresponding to any parameters that we want to set. (Any structure fields that we don't include are equivalent to not setting those parameters, and using the defaults instead). You can get more information on the available parameters by typing `help` `nlopt_optimize` in Matlab. The equivalent of the C example above is to define an `opt` structure by:
322 opt.algorithm = NLOPT_LD_MMA
323 opt.lower_bounds = [-inf, 0]
324 opt.min_objective = @myfunc
325 opt.fc = { (@(x) myconstraint(x,2,0)), (@(x) myconstraint(x,-1,1)) }
326 opt.fc_tol = [1e-8, 1e-8];
331 We do not need to specify the dimension of the problem; this is implicitly specified by the size of the initial-guess vector passed to `nlopt_optimize` below (and must match the sizes of other vectors like `opt.lower_bounds`). The inequality constraints are specified as a [cell array](http://blogs.mathworks.com/loren/2006/06/21/cell-arrays-and-their-contents/) `opt.fc` of function handles (and the corresponding tolerances are in an array `opt.fc_tol`); notice how we use `@(x)` to define an anonymous/inline function in order to pass additional arguments to `myconstraint`.
333 Finally, we call `nlopt_optimize`:
336 [xopt, fmin, retcode] = nlopt_optimize(opt, [1.234 5.678])
340 `nlopt_optimize` returns three things: `xopt`, the optimal parameters found; `fmin`, the corresponding value of the objective function, and a return code `retcode` (positive on success and negative on failure).
342 The output of the above command is:
352 (The [return code](NLopt_Reference#Return_values.md) `4` corresponds to `NLOPT_XTOL_REACHED`, which means it converged to the specified *x* tolerance.) To switch to a derivative-free algorithm like COBYLA, we just change `opt.algorithm` parameter:
355 opt.algorithm = NLOPT_LN_COBYLA
359 ### Matlab verbose output
361 It is often useful to print out some status message to see what is happening, especially if your function evaluation is much slower or if a large number of evaluations are required (e.g. for global optimization). You can, of course, modify your function to print out whatever you want. As a shortcut, however, you can set a verbose option in NLopt's Matlab interface by:
368 If we do this, then running the MMA algorithm as above yields:
371 nlopt_minimize_constrained eval #1: 2.38286
372 nlopt_minimize_constrained eval #2: 2.35613
373 nlopt_minimize_constrained eval #3: 2.24586
374 nlopt_minimize_constrained eval #4: 2.0191
375 nlopt_minimize_constrained eval #5: 1.74093
376 nlopt_minimize_constrained eval #6: 1.40421
377 nlopt_minimize_constrained eval #7: 1.0223
378 nlopt_minimize_constrained eval #8: 0.685203
379 nlopt_minimize_constrained eval #9: 0.552985
380 nlopt_minimize_constrained eval #10: 0.544354
381 nlopt_minimize_constrained eval #11: 0.544331
385 This shows the objective function values at each intermediate step of the optimization. As in the C example above, it converges in 11 steps. The COBYLA algorithm requires a few more iterations, because it doesn't exploit the gradient information:
388 nlopt_optimize eval #1: 2.38286
389 nlopt_optimize eval #2: 2.38286
390 nlopt_optimize eval #3: 3.15222
391 nlopt_optimize eval #4: 1.20627
392 nlopt_optimize eval #5: 0.499441
393 nlopt_optimize eval #6: 0.709216
394 nlopt_optimize eval #7: 0.20341
395 nlopt_optimize eval #8: 0.745201
396 nlopt_optimize eval #9: 0.989693
397 nlopt_optimize eval #10: 0.324679
398 nlopt_optimize eval #11: 0.804318
399 nlopt_optimize eval #12: 0.541431
400 nlopt_optimize eval #13: 0.561137
401 nlopt_optimize eval #14: 0.531346
402 nlopt_optimize eval #15: 0.515784
403 nlopt_optimize eval #16: 0.541724
404 nlopt_optimize eval #17: 0.541422
405 nlopt_optimize eval #18: 0.541504
406 nlopt_optimize eval #19: 0.541499
407 nlopt_optimize eval #20: 0.541514
408 nlopt_optimize eval #21: 0.541352
409 nlopt_optimize eval #22: 0.542089
410 nlopt_optimize eval #23: 0.542575
411 nlopt_optimize eval #24: 0.543027
412 nlopt_optimize eval #25: 0.544911
413 nlopt_optimize eval #26: 0.541793
414 nlopt_optimize eval #27: 0.545225
415 nlopt_optimize eval #28: 0.544331
416 nlopt_optimize eval #29: 0.544256
417 nlopt_optimize eval #30: 0.544242
418 nlopt_optimize eval #31: 0.544116
422 Notice that some of the objective function values are below the minimum of 0.54433 — these are simply values of the objective function at infeasible points (violating the nonlinear constraints).
427 The same example in Python is:
435 grad[1] = 0.5 / sqrt(x[1])
437 def myconstraint(x, grad, a, b):
439 grad[0] = 3 * a * (a*x[0] + b)**2
441 return (a*x[0] + b)**3 - x[1]
442 opt = nlopt.opt(nlopt.LD_MMA, 2)
443 opt.set_lower_bounds([-float('inf'), 0])
444 opt.set_min_objective(myfunc)
445 opt.add_inequality_constraint(lambda x,grad: myconstraint(x,grad,2,0), 1e-8)
446 opt.add_inequality_constraint(lambda x,grad: myconstraint(x,grad,-1,1), 1e-8)
447 opt.set_xtol_rel(1e-4)
448 x = opt.optimize([1.234, 5.678])
449 minf = opt.last_optimum_value()
450 print "optimum at ", x[0],x[1]
451 print "minimum value = ", minf
452 print "result code = ", opt.last_optimize_result()
456 Notice that the `optimize` method returns only the location of the optimum (as a NumPy array), and that the value of the optimum and the result code are obtained by `last_optimum_value` and `last_optimize_result` values. Like in C++, the NLopt functions raise exceptions on errors, so we don't need to check return codes to look for errors.
458 The objective and constraint functions take NumPy arrays as arguments; if the `grad` argument is non-empty it must be modified *in-place* to the value of the gradient. Notice how we use Python's `lambda` construct to pass additional parameters to the constraints. Alternatively, we could define the objective/constraints as classes with a `__call__(self,` `x,` `grad)` method so that they can behave like functions.
460 The result of running the above code should be:
463 optimum at 0.333333331366 0.296296292697
464 minimum value = 0.544331050646
469 finding the same correct optimum as in the C interface (of course). (The [return code](NLopt_Reference#Return_values.md) `4` corresponds to `nlopt.XTOL_REACHED`, which means it converged to the specified *x* tolerance.)
471 ### Important: Modifying `grad` in-place
473 The grad argument of your objective/constraint functions must be modified *in-place*. If you use an operation like
480 however, Python allocates a new array to hold `2*x` and reassigns grad to point to it, rather than modifying the original contents of grad. **This will not work.** Instead, you should do:
487 which *overwrites* the old contents of grad with `2*x`. See also the [NLopt Python Reference](NLopt_Python_Reference#Assigning_results_in-place.md).
489 Example in GNU Guile (Scheme)
490 -----------------------------
492 In [GNU Guile](https://en.wikipedia.org/wiki/GNU_Guile), which is an implementation of the [Scheme programming language](https://en.wikipedia.org/wiki/Scheme_(programming_language)), the equivalent of the example above would be:
495 (use-modules (nlopt))
496 (define (myfunc x grad)
499 (vector-set! grad 0 0.0)
500 (vector-set! grad 1 (/ 0.5 (sqrt (vector-ref x 1))))))
501 (sqrt (vector-ref x 1)))
502 (define (myconstraint x grad a b)
503 (let ((x0 (vector-ref x 0)) (x1 (vector-ref x 1)))
506 (vector-set! grad 0 (* 3 a (expt (+ (* a x0) b) 2)))
507 (vector-set! grad 1 -1.0)))
508 (- (expt (+ (* a x0) b) 3) x1)))
509 (define opt (new-nlopt-opt NLOPT-LD-MMA 2))
510 (nlopt-opt-set-lower-bounds opt (vector (- (inf)) 0))
511 (nlopt-opt-set-min-objective opt myfunc)
512 (nlopt-opt-add-inequality-constraint opt (lambda (x grad)
513 (myconstraint x grad 2 0))
515 (nlopt-opt-add-inequality-constraint opt (lambda (x grad)
516 (myconstraint x grad -1 1))
518 (nlopt-opt-set-xtol-rel opt 1e-4)
519 (define x (nlopt-opt-optimize opt (vector 1.234 5.678)))
520 (define minf (nlopt-opt-last-optimum-value opt))
521 (define result (nlopt-opt-last-optimize-result opt))
525 Note that the objective/constraint functions take two arguments, `x` and `grad`, and return a number. `x` is a vector whose length is the dimension of the problem; grad is either false (`#f`) if it is not needed, or a `vector` that must be modified *in-place* to the gradient of the function.
527 On error conditions, the NLopt functions throw [exceptions](http://www.gnu.org/software/guile/manual/html_node/Exceptions.html) that can be caught by your Scheme code if you wish.
529 The heavy use of side-effects here is a bit unnatural in Scheme, but is used in order to closely map to the C++ interface. (Notice that `nlopt::` C++ functions map to `nlopt-` Guile functions, and `nlopt::opt::` methods map to `nlopt-opt-` functions that take the `opt` object as the first argument.) Of course, you are free to wrap your own Scheme-like functional interface around this if you wish.
534 In [Fortran](https://en.wikipedia.org/wiki/Fortran), the equivalent of the C example above would be as follows. First, we would write our functions as:
537 subroutine myfunc(val, n, x, grad, need_gradient, f_data)
538 double precision val, x(n), grad(n)
539 integer n, need_gradient
540 if (need_gradient.ne.0) then
542 grad(2) = 0.5 / dsqrt(x(2))
550 subroutine myconstraint(val, n, x, grad, need_gradient, d)
551 integer need_gradient
552 double precision val, x(n), grad(n), d(2), a, b
555 if (need_gradient.ne.0) then
556 grad(1) = 3. * a * (a*x(1) + b)**2
559 val = (a*x(1) + b)**3 - x(2)
564 Notice that that the "functions" are actually subroutines. This is because it turns out to be hard to call Fortran functions from C or vice versa in any remotely portable way. Therefore:
566 - In the NLopt Fortran interface, all C functions become subroutines in Fortran, with the return value becoming the first argument.
568 So, here the first argument `val` is used for the return value. Also, because there is no way in Fortran to pass `NULL` for the `grad` argument, we add an additional `need_gradient` argument which is nonzero if the gradient needs to be computed. Finally, the last argument is the equivalent of the `void*` argument in the C API, and can be used to pass a single argument of any type through to the objective/constraint functions: here, we use it in `myconstraint` to pass an array of two values for the constants *a* and *b*.
570 Then, to run the optimization, we can use the following Fortran program:
574 external myfunc, myconstraint
575 double precision lb(2)
577 double precision d1(2), d2(2)
578 double precision x(2), minf
582 call nlo_create(opt, NLOPT_LD_MMA, 2)
583 call nlo_get_lower_bounds(ires, opt, lb)
585 call nlo_set_lower_bounds(ires, opt, lb)
586 call nlo_set_min_objective(ires, opt, myfunc, 0)
590 call nlo_add_inequality_constraint(ires, opt,
591 $ myconstraint, d1, 1.D-8)
594 call nlo_add_inequality_constraint(ires, opt,
595 $ myconstraint, d2, 1.D-8)
597 call nlo_set_xtol_rel(ires, opt, 1.D-4)
601 call nlo_optimize(ires, opt, x, minf)
603 write(*,*) 'nlopt failed!'
605 write(*,*) 'found min at ', x(1), x(2)
606 write(*,*) 'min val = ', minf
609 call nlo_destroy(opt)
615 There are a few things to note here:
617 - All `nlopt_` functions are converted into `nlo_` subroutines, with return values converted into the first argument.
618 - The "`nlopt_opt`" variable `opt` is declared as `integer*8`. (Technically, we could use any type that is big enough to hold a pointer on all platforms; `integer*8` is big enough for pointers on both 32-bit and 64-bit machines.)
619 - The subroutines must be declared as `external`.
620 - We `include` `'nlopt.f'` in order to get the various constants like `NLOPT_LD_MMA`.
622 There is no standard Fortran 77 equivalent of C's `HUGE_VAL` constant, so instead we just call `nlo_get_lower_bounds` to get the default lower bounds (-∞) and then change one of them. In Fortran 90 (and supported as an extension in many Fortran 77 compilers), there is a `huge` intrinsic function that we could have used instead:
632 The same example in the [Julia programming language](https://en.wikipedia.org/wiki/Julia_(programming_language)) can be found at the [NLopt.jl](https://github.com/stevengj/NLopt.jl) web page.
634 [Category:NLopt](index.md)