chiark / gitweb /
work around broken solaris gcc
[nlopt.git] / api / nlopt.c
1 /* Copyright (c) 2007-2008 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 <stdlib.h>
24 #include <math.h>
25 #include <string.h>
26 #include <float.h>
27
28 #include "nlopt.h"
29 #include "nlopt-util.h"
30
31 #define MIN(a,b) ((a) < (b) ? (a) : (b))
32
33 /*************************************************************************/
34
35 int nlopt_isinf(double x) {
36      return fabs(x) >= HUGE_VAL * 0.99
37 #ifdef HAVE_ISINF
38           || isinf(x)
39 #endif
40           ;
41 }
42
43 #ifndef HAVE_ISNAN
44 static int my_isnan(double x) { return x != x; }
45 #  define isnan my_isnan
46 #endif
47
48 /*************************************************************************/
49
50 void nlopt_version(int *major, int *minor, int *bugfix)
51 {
52      *major = MAJOR_VERSION;
53      *minor = MINOR_VERSION;
54      *bugfix = BUGFIX_VERSION;
55 }
56
57 /*************************************************************************/
58
59 static const char nlopt_algorithm_names[NLOPT_NUM_ALGORITHMS][256] = {
60      "DIRECT (global, no-derivative)",
61      "DIRECT-L (global, no-derivative)",
62      "Randomized DIRECT-L (global, no-derivative)",
63      "Unscaled DIRECT (global, no-derivative)",
64      "Unscaled DIRECT-L (global, no-derivative)",
65      "Unscaled Randomized DIRECT-L (global, no-derivative)",
66      "Original DIRECT version (global, no-derivative)",
67      "Original DIRECT-L version (global, no-derivative)",
68 #ifdef WITH_CXX
69      "StoGO (global, derivative-based)",
70      "StoGO with randomized search (global, derivative-based)",
71 #else
72      "StoGO (NOT COMPILED)",
73      "StoGO randomized (NOT COMPILED)",
74 #endif
75 #ifdef WITH_NOCEDAL_LBFGS
76      "original NON-FREE L-BFGS code by Nocedal et al. (local, deriv.-based)",
77 #else
78      "original NON-FREE L-BFGS code by Nocedal et al. (NOT COMPILED)",
79 #endif
80      "Limited-memory BFGS (L-BFGS) (local, derivative-based)",
81      "Principal-axis, praxis (local, no-derivative)",
82      "Limited-memory variable-metric, rank 1 (local, derivative-based)",
83      "Limited-memory variable-metric, rank 2 (local, derivative-based)",
84      "Truncated Newton (local, derivative-based)",
85      "Truncated Newton with restarting (local, derivative-based)",
86      "Preconditioned truncated Newton (local, derivative-based)",
87      "Preconditioned truncated Newton with restarting (local, derivative-based)",
88      "Controlled random search (CRS2) with local mutation (global, no-derivative)",
89      "Multi-level single-linkage (MLSL), random (global, no-derivative)",
90      "Multi-level single-linkage (MLSL), random (global, derivative)",
91      "Multi-level single-linkage (MLSL), quasi-random (global, no-derivative)",
92      "Multi-level single-linkage (MLSL), quasi-random (global, derivative)",
93      "Method of Moving Asymptotes (MMA) (local, derivative)",
94      "COBYLA (Constrained Optimization BY Linear Approximations) (local, no-derivative)",
95      "NEWUOA unconstrained optimization via quadratic models (local, no-derivative)",
96      "Bound-constrained optimization via NEWUOA-based quadratic models (local, no-derivative)",
97      "Nelder-Mead simplex algorithm (local, no-derivative)",
98      "Sbplx variant of Nelder-Mead (re-implementation of Rowan's Subplex) (local, no-derivative)"
99 };
100
101 const char *nlopt_algorithm_name(nlopt_algorithm a)
102 {
103      if (a < 0 || a >= NLOPT_NUM_ALGORITHMS) return "UNKNOWN";
104      return nlopt_algorithm_names[a];
105 }
106
107 /*************************************************************************/
108
109 static int nlopt_srand_called = 0;
110 void nlopt_srand(unsigned long seed) {
111      nlopt_srand_called = 1;
112      nlopt_init_genrand(seed);
113 }
114
115 void nlopt_srand_time() {
116      nlopt_srand(nlopt_time_seed());
117 }
118
119 /*************************************************************************/
120
121 /* crude heuristics for initial step size of nonderivative algorithms */
122 static double initial_step(int n, 
123                            const double *lb, const double *ub, const double *x)
124 {
125      int i;
126      double step = HUGE_VAL;
127      for (i = 0; i < n; ++i) {
128           if (!nlopt_isinf(ub[i]) && !nlopt_isinf(lb[i])
129               && (ub[i] - lb[i]) * 0.25 < step && ub[i] > lb[i])
130                step = (ub[i] - lb[i]) * 0.25;
131           if (!nlopt_isinf(ub[i]) 
132               && ub[i] - x[i] < step && ub[i] > x[i])
133                step = ub[i] - x[i];
134           if (!nlopt_isinf(lb[i]) 
135               && x[i] - lb[i] < step && x[i] > lb[i])
136                step = x[i] - lb[i];
137      }
138      if (nlopt_isinf(step))
139           for (i = 0; i < n; ++i) {
140                if (!nlopt_isinf(ub[i]) 
141                    && ub[i] - x[i] < step && ub[i] > x[i] + 1e-4)
142                     step = ub[i] - x[i];
143                if (!nlopt_isinf(lb[i]) 
144                    && x[i] - lb[i] < step && x[i] > lb[i] + 1e-4)
145                     step = x[i] - lb[i];
146           }
147      if (nlopt_isinf(step)) {
148           step = 0;
149           for (i = 0; i < n; ++i)
150                if (fabs(x[i]) * 0.25 > step)
151                     step = fabs(x[i]) * 0.25;
152           if (step == 0)
153                step = 1;
154      }
155      return step;
156 }
157               
158 /*************************************************************************/
159
160 typedef struct {
161      nlopt_func f;
162      void *f_data;
163      const double *lb, *ub;
164 } nlopt_data;
165
166 #include "praxis.h"
167
168 static double f_bound(int n, const double *x, void *data_)
169 {
170      int i;
171      nlopt_data *data = (nlopt_data *) data_;
172      double f;
173
174      /* some methods do not support bound constraints, but support
175         discontinuous objectives so we can just return Inf for invalid x */
176      for (i = 0; i < n; ++i)
177           if (x[i] < data->lb[i] || x[i] > data->ub[i])
178                return HUGE_VAL;
179
180      f = data->f(n, x, NULL, data->f_data);
181      return (isnan(f) || nlopt_isinf(f) ? HUGE_VAL : f);
182 }
183
184 static double f_noderiv(int n, const double *x, void *data_)
185 {
186      nlopt_data *data = (nlopt_data *) data_;
187      return data->f(n, x, NULL, data->f_data);
188 }
189
190 #include "direct.h"
191
192 static double f_direct(int n, const double *x, int *undefined, void *data_)
193 {
194      nlopt_data *data = (nlopt_data *) data_;
195      double f;
196      f = data->f(n, x, NULL, data->f_data);
197      *undefined = isnan(f) || nlopt_isinf(f);
198      return f;
199 }
200
201 #ifdef WITH_CXX
202 #  include "stogo.h"
203 #endif
204
205 #include "cdirect.h"
206
207 #ifdef WITH_NOCEDAL
208 #  include "l-bfgs-b.h"
209 #endif
210
211 #include "luksan.h"
212
213 #include "crs.h"
214
215 #include "mlsl.h"
216 #include "mma.h"
217 #include "cobyla.h"
218 #include "newuoa.h"
219 #include "neldermead.h"
220
221 /*************************************************************************/
222
223 /* for "hybrid" algorithms that combine global search with some
224    local search algorithm, most of the time we anticipate that people
225    will stick with the default local search algorithm, so we
226    don't add this as a parameter to nlopt_minimize.  Instead, the user
227    can call nlopt_{set/get}_hybrid_local_algorithm to get/set the defaults. */
228
229 /* default local-search algorithms */
230 static nlopt_algorithm local_search_alg_deriv = NLOPT_LD_MMA;
231 static nlopt_algorithm local_search_alg_nonderiv = NLOPT_LN_COBYLA;
232
233 static int local_search_maxeval = -1; /* no maximum by default */
234
235 void nlopt_get_local_search_algorithm(nlopt_algorithm *deriv,
236                                       nlopt_algorithm *nonderiv,
237                                       int *maxeval)
238 {
239      *deriv = local_search_alg_deriv;
240      *nonderiv = local_search_alg_nonderiv;
241      *maxeval = local_search_maxeval;
242 }
243
244 void nlopt_set_local_search_algorithm(nlopt_algorithm deriv,
245                                       nlopt_algorithm nonderiv,
246                                       int maxeval)
247 {
248      local_search_alg_deriv = deriv;
249      local_search_alg_nonderiv = nonderiv;
250      local_search_maxeval = maxeval;
251 }
252
253 /*************************************************************************/
254
255 /* same as nlopt_minimize, but xtol_abs is required to be non-NULL */
256 static nlopt_result nlopt_minimize_(
257      nlopt_algorithm algorithm,
258      int n, nlopt_func f, void *f_data,
259      int m, nlopt_func fc, void *fc_data, ptrdiff_t fc_datum_size,
260      const double *lb, const double *ub, /* bounds */
261      double *x, /* in: initial guess, out: minimizer */
262      double *minf, /* out: minimum */
263      double minf_max, double ftol_rel, double ftol_abs,
264      double xtol_rel, const double *xtol_abs,
265      int maxeval, double maxtime)
266 {
267      int i;
268      nlopt_data d;
269      nlopt_stopping stop;
270
271      /* some basic argument checks */
272      if (!minf || !f || n < 0 || m < 0
273           || (m > 0 && !fc)) return NLOPT_INVALID_ARGS;
274      if (n == 0) { /* trivial case: no degrees of freedom */
275           *minf = f(n, x, NULL, f_data);
276           return NLOPT_SUCCESS;
277      }
278      else if (n < 0 || !lb || !ub || !x)
279           return NLOPT_INVALID_ARGS;
280
281      /* nonlinear constraints are only supported with MMA or COBYLA */
282      if (m != 0 && algorithm != NLOPT_LD_MMA && algorithm != NLOPT_LN_COBYLA) 
283           return NLOPT_INVALID_ARGS;
284
285      d.f = f;
286      d.f_data = f_data;
287      d.lb = lb;
288      d.ub = ub;
289
290      /* make sure rand generator is inited */
291      if (!nlopt_srand_called)
292           nlopt_srand_time(); /* default is non-deterministic */
293
294      /* check bound constraints */
295      for (i = 0; i < n; ++i)
296           if (lb[i] > ub[i] || x[i] < lb[i] || x[i] > ub[i])
297                return NLOPT_INVALID_ARGS;
298
299      stop.n = n;
300      stop.minf_max = (isnan(minf_max) 
301                       || (nlopt_isinf(minf_max) && minf_max < 0))
302           ? -HUGE_VAL : minf_max;
303      stop.ftol_rel = ftol_rel;
304      stop.ftol_abs = ftol_abs;
305      stop.xtol_rel = xtol_rel;
306      stop.xtol_abs = xtol_abs;
307      stop.nevals = 0;
308      stop.maxeval = maxeval;
309      stop.maxtime = maxtime;
310      stop.start = nlopt_seconds();
311
312      switch (algorithm) {
313          case NLOPT_GN_DIRECT:
314          case NLOPT_GN_DIRECT_L: 
315          case NLOPT_GN_DIRECT_L_RAND: 
316               return cdirect(n, f, f_data, lb, ub, x, minf, &stop, 0.0, 
317                              (algorithm != NLOPT_GN_DIRECT)
318                              + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
319                              + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
320               
321          case NLOPT_GN_DIRECT_NOSCAL:
322          case NLOPT_GN_DIRECT_L_NOSCAL: 
323          case NLOPT_GN_DIRECT_L_RAND_NOSCAL: 
324               return cdirect_unscaled(n, f, f_data, lb, ub, x, minf, 
325                                       &stop, 0.0, 
326                                       (algorithm != NLOPT_GN_DIRECT)
327                                       + 3 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 2 : (algorithm != NLOPT_GN_DIRECT))
328                                       + 9 * (algorithm == NLOPT_GN_DIRECT_L_RAND ? 1 : (algorithm != NLOPT_GN_DIRECT)));
329               
330          case NLOPT_GN_ORIG_DIRECT:
331          case NLOPT_GN_ORIG_DIRECT_L: 
332               switch (direct_optimize(f_direct, &d, n, lb, ub, x, minf,
333                                       maxeval, -1, 0.0, 0.0,
334                                       pow(xtol_rel, (double) n), -1.0,
335                                       stop.minf_max, 0.0,
336                                       NULL, 
337                                       algorithm == NLOPT_GN_ORIG_DIRECT
338                                       ? DIRECT_ORIGINAL
339                                       : DIRECT_GABLONSKY)) {
340                   case DIRECT_INVALID_BOUNDS:
341                   case DIRECT_MAXFEVAL_TOOBIG:
342                   case DIRECT_INVALID_ARGS:
343                        return NLOPT_INVALID_ARGS;
344                   case DIRECT_INIT_FAILED:
345                   case DIRECT_SAMPLEPOINTS_FAILED:
346                   case DIRECT_SAMPLE_FAILED:
347                        return NLOPT_FAILURE;
348                   case DIRECT_MAXFEVAL_EXCEEDED:
349                   case DIRECT_MAXITER_EXCEEDED:
350                        return NLOPT_MAXEVAL_REACHED;
351                   case DIRECT_GLOBAL_FOUND:
352                        return NLOPT_MINF_MAX_REACHED;
353                   case DIRECT_VOLTOL:
354                   case DIRECT_SIGMATOL:
355                        return NLOPT_XTOL_REACHED;
356                   case DIRECT_OUT_OF_MEMORY:
357                        return NLOPT_OUT_OF_MEMORY;
358               break;
359          }
360
361          case NLOPT_GD_STOGO:
362          case NLOPT_GD_STOGO_RAND:
363 #ifdef WITH_CXX
364               if (!stogo_minimize(n, f, f_data, x, minf, lb, ub, &stop,
365                                   algorithm == NLOPT_GD_STOGO
366                                   ? 0 : 2*n))
367                    return NLOPT_FAILURE;
368               break;
369 #else
370               return NLOPT_FAILURE;
371 #endif
372
373 #if 0
374               /* lacking a free/open-source license, we no longer use
375                  Rowan's code, and instead use by "sbplx" re-implementation */
376          case NLOPT_LN_SUBPLEX: {
377               int iret;
378               double *scale = (double *) malloc(sizeof(double) * n);
379               if (!scale) return NLOPT_OUT_OF_MEMORY;
380               for (i = 0; i < n; ++i)
381                    scale[i] = initial_step(1, lb+i, ub+i, x+i);
382               iret = nlopt_subplex(f_bound, minf, x, n, &d, &stop, scale);
383               free(scale);
384               switch (iret) {
385                   case -2: return NLOPT_INVALID_ARGS;
386                   case -10: return NLOPT_MAXTIME_REACHED;
387                   case -1: return NLOPT_MAXEVAL_REACHED;
388                   case 0: return NLOPT_XTOL_REACHED;
389                   case 1: return NLOPT_SUCCESS;
390                   case 2: return NLOPT_MINF_MAX_REACHED;
391                   case 20: return NLOPT_FTOL_REACHED;
392                   case -200: return NLOPT_OUT_OF_MEMORY;
393                   default: return NLOPT_FAILURE; /* unknown return code */
394               }
395               break;
396          }
397 #endif
398
399          case NLOPT_LN_PRAXIS:
400               return praxis_(0.0, DBL_EPSILON, 
401                              initial_step(n, lb, ub, x), n, x, f_bound, &d,
402                              &stop, minf);
403
404 #ifdef WITH_NOCEDAL
405          case NLOPT_LD_LBFGS_NOCEDAL: {
406               int iret, *nbd = (int *) malloc(sizeof(int) * n);
407               if (!nbd) return NLOPT_OUT_OF_MEMORY;
408               for (i = 0; i < n; ++i) {
409                    int linf = nlopt_isinf(lb[i]) && lb[i] < 0;
410                    int uinf = nlopt_isinf(ub[i]) && ub[i] > 0;
411                    nbd[i] = linf && uinf ? 0 : (uinf ? 1 : (linf ? 3 : 2));
412               }
413               iret = lbfgsb_minimize(n, f, f_data, x, nbd, lb, ub,
414                                      MIN(n, 5), 0.0, ftol_rel, 
415                                      xtol_abs ? *xtol_abs : xtol_rel,
416                                      maxeval);
417               free(nbd);
418               if (iret <= 0) {
419                    switch (iret) {
420                        case -1: return NLOPT_INVALID_ARGS;
421                        case -2: default: return NLOPT_FAILURE;
422                    }
423               }
424               else {
425                    *minf = f(n, x, NULL, f_data);
426                    switch (iret) {
427                        case 5: return NLOPT_MAXEVAL_REACHED;
428                        case 2: return NLOPT_XTOL_REACHED;
429                        case 1: return NLOPT_FTOL_REACHED;
430                        default: return NLOPT_SUCCESS;
431                    }
432               }
433               break;
434          }
435 #endif
436
437          case NLOPT_LD_LBFGS: 
438               return luksan_plis(n, f, f_data, lb, ub, x, minf, &stop);
439
440          case NLOPT_LD_VAR1: 
441          case NLOPT_LD_VAR2: 
442               return luksan_plip(n, f, f_data, lb, ub, x, minf, &stop,
443                    algorithm == NLOPT_LD_VAR1 ? 1 : 2);
444
445          case NLOPT_LD_TNEWTON: 
446          case NLOPT_LD_TNEWTON_RESTART: 
447          case NLOPT_LD_TNEWTON_PRECOND: 
448          case NLOPT_LD_TNEWTON_PRECOND_RESTART: 
449               return luksan_pnet(n, f, f_data, lb, ub, x, minf, &stop,
450                                  1 + (algorithm - NLOPT_LD_TNEWTON) % 2,
451                                  1 + (algorithm - NLOPT_LD_TNEWTON) / 2);
452
453          case NLOPT_GN_CRS2_LM:
454               return crs_minimize(n, f, f_data, lb, ub, x, minf, &stop, 0);
455
456          case NLOPT_GN_MLSL:
457          case NLOPT_GD_MLSL:
458          case NLOPT_GN_MLSL_LDS:
459          case NLOPT_GD_MLSL_LDS:
460               for (i = 0; i < n && stop.xtol_abs[i] <= 0; ++i) ;
461               if (stop.ftol_rel <= 0 && stop.ftol_abs <= 0 &&
462                   stop.xtol_rel <= 0 && i == n) {
463                    /* it is not sensible to call MLSL without *some*
464                       nonzero tolerance for the local search */
465                    stop.ftol_rel = 1e-15;
466                    stop.xtol_rel = 1e-7;
467               }
468               return mlsl_minimize(n, f, f_data, lb, ub, x, minf, &stop,
469                                    (algorithm == NLOPT_GN_MLSL ||
470                                     algorithm == NLOPT_GN_MLSL_LDS)
471                                    ? local_search_alg_nonderiv
472                                    : local_search_alg_deriv,
473                                    local_search_maxeval,
474                                    algorithm >= NLOPT_GN_MLSL_LDS);
475
476          case NLOPT_LD_MMA:
477               return mma_minimize(n, f, f_data, m, fc, fc_data, fc_datum_size,
478                                   lb, ub, x, minf, &stop,
479                                   local_search_alg_deriv, 1e-12, 100000);
480
481          case NLOPT_LN_COBYLA:
482               return cobyla_minimize(n, f, f_data, 
483                                      m, fc, fc_data, fc_datum_size,
484                                      lb, ub, x, minf, &stop,
485                                      initial_step(n, lb, ub, x));
486                                      
487          case NLOPT_LN_NEWUOA:
488               return newuoa(n, 2*n+1, x, 0, 0, initial_step(n, lb, ub, x),
489                             &stop, minf, f_noderiv, &d);
490                                      
491          case NLOPT_LN_NEWUOA_BOUND:
492               return newuoa(n, 2*n+1, x, lb, ub, initial_step(n, lb, ub, x),
493                             &stop, minf, f_noderiv, &d);
494
495          case NLOPT_LN_NELDERMEAD: 
496          case NLOPT_LN_SBPLX: 
497          {
498               nlopt_result ret;
499               double *scale = (double *) malloc(sizeof(double) * n);
500               if (!scale) return NLOPT_OUT_OF_MEMORY;
501               for (i = 0; i < n; ++i)
502                    scale[i] = initial_step(1, lb+i, ub+i, x+i);
503               if (algorithm == NLOPT_LN_NELDERMEAD)
504                    ret = nldrmd_minimize(n,f,f_data,lb,ub,x,minf,scale,&stop);
505               else
506                    ret = sbplx_minimize(n,f,f_data,lb,ub,x,minf,scale,&stop);
507               free(scale);
508               return ret;
509          }
510
511          default:
512               return NLOPT_INVALID_ARGS;
513      }
514
515      return NLOPT_SUCCESS;
516 }
517
518 nlopt_result nlopt_minimize_constrained(
519      nlopt_algorithm algorithm,
520      int n, nlopt_func f, void *f_data,
521      int m, nlopt_func fc, void *fc_data, ptrdiff_t fc_datum_size,
522      const double *lb, const double *ub, /* bounds */
523      double *x, /* in: initial guess, out: minimizer */
524      double *minf, /* out: minimum */
525      double minf_max, double ftol_rel, double ftol_abs,
526      double xtol_rel, const double *xtol_abs,
527      int maxeval, double maxtime)
528 {
529      nlopt_result ret;
530      if (xtol_abs)
531           ret = nlopt_minimize_(algorithm, n, f, f_data,
532                                 m, fc, fc_data, fc_datum_size, lb, ub,
533                                 x, minf, minf_max, ftol_rel, ftol_abs,
534                                 xtol_rel, xtol_abs, maxeval, maxtime);
535      else {
536           int i;
537           double *xtol = (double *) malloc(sizeof(double) * n);
538           if (!xtol) return NLOPT_OUT_OF_MEMORY;
539           for (i = 0; i < n; ++i) xtol[i] = -1;
540           ret = nlopt_minimize_(algorithm, n, f, f_data, 
541                                 m, fc, fc_data, fc_datum_size, lb, ub,
542                                 x, minf, minf_max, ftol_rel, ftol_abs,
543                                 xtol_rel, xtol, maxeval, maxtime);
544           free(xtol);
545      }
546      return ret;
547 }
548
549 nlopt_result nlopt_minimize(
550      nlopt_algorithm algorithm,
551      int n, nlopt_func f, void *f_data,
552      const double *lb, const double *ub, /* bounds */
553      double *x, /* in: initial guess, out: minimizer */
554      double *minf, /* out: minimum */
555      double minf_max, double ftol_rel, double ftol_abs,
556      double xtol_rel, const double *xtol_abs,
557      int maxeval, double maxtime)
558 {
559      return nlopt_minimize_constrained(
560           algorithm, n, f, f_data, 0, NULL, NULL, 0,
561           lb, ub, x, minf, minf_max, ftol_rel, ftol_abs,
562           xtol_rel, xtol_abs, maxeval, maxtime);
563 }