scipy.optimize.direct#
- scipy.optimize.direct(func, bounds, *, args=(), eps=0.0001, maxfun=None, maxiter=1000, locally_biased=True, f_min=- inf, f_min_rtol=0.0001, vol_tol=1e-16, len_tol=1e-06, callback=None)[source]#
Finds the global minimum of a function using the DIRECT algorithm.
- Parameters
- funccallable
The objective function to be minimized.
func(x, *args) -> float
wherex
is an 1-D array with shape (n,) andargs
is a tuple of the fixed parameters needed to completely specify the function.- boundssequence or
Bounds
Bounds for variables. There are two ways to specify the bounds:
Instance of
Bounds
class.(min, max)
pairs for each element inx
.
- argstuple, optional
Any additional fixed parameters needed to completely specify the objective function.
- epsfloat, optional
Minimal required difference of the objective function values between the current best hyperrectangle and the next potentially optimal hyperrectangle to be divided. In consequence, eps serves as a tradeoff between local and global search: the smaller, the more local the search becomes. Default is 1e-4.
- maxfunint or None, optional
Approximate upper bound on objective function evaluations. If None, will be automatically set to
1000 * N
whereN
represents the number of dimensions. Will be capped if necessary to limit DIRECT’s RAM usage to app. 1GiB. This will only occur for very high dimensional problems and excessive max_fun. Default is None.- maxiterint, optional
Maximum number of iterations. Default is 1000.
- locally_biasedbool, optional
If True (default), use the locally biased variant of the algorithm known as DIRECT_L. If False, use the original unbiased DIRECT algorithm. For hard problems with many local minima, False is recommended.
- f_minfloat, optional
Function value of the global optimum. Set this value only if the global optimum is known. Default is
-np.inf
, so that this termination criterion is deactivated.- f_min_rtolfloat, optional
Terminate the optimization once the relative error between the current best minimum f and the supplied global minimum f_min is smaller than f_min_rtol. This parameter is only used if f_min is also set. Default is 1e-4.
- vol_tolfloat, optional
Terminate the optimization once the volume of the hyperrectangle containing the lowest function value is smaller than vol_tol of the complete search space. Must lie between 0 and 1. Default is 1e-16.
- len_tolfloat, optional
If locally_biased=True, terminate the optimization once half of the normalized maximal side length of the hyperrectangle containing the lowest function value is smaller than len_tol. If locally_biased=False, terminate the optimization once half of the normalized diagonal of the hyperrectangle containing the lowest function value is smaller than len_tol. Must lie between 0 and 1. Default is 1e-6.
- callbackcallable, optional
A callback function with signature
callback(xk)
wherexk
represents the best function value found so far.
- Returns
- resOptimizeResult
The optimization result represented as a
OptimizeResult
object. Important attributes are:x
the solution array,success
a Boolean flag indicating if the optimizer exited successfully andmessage
which describes the cause of the termination. SeeOptimizeResult
for a description of other attributes.
Notes
DIviding RECTangles (DIRECT) is a deterministic global optimization algorithm capable of minimizing a black box function with its variables subject to lower and upper bound constraints by sampling potential solutions in the search space [1]. The algorithm starts by normalising the search space to an n-dimensional unit hypercube. It samples the function at the center of this hypercube and at 2n (n is the number of variables) more points, 2 in each coordinate direction. Using these function values, DIRECT then divides the domain into hyperrectangles, each having exactly one of the sampling points as its center. In each iteration, DIRECT chooses, using the eps parameter which defaults to 1e-4, some of the existing hyperrectangles to be further divided. This division process continues until either the maximum number of iterations or maximum function evaluations allowed are exceeded, or the hyperrectangle containing the minimal value found so far becomes small enough. If f_min is specified, the optimization will stop once this function value is reached within a relative tolerance. The locally biased variant of DIRECT (originally called DIRECT_L) [2] is used by default. It makes the search more locally biased and more efficient for cases with only a few local minima.
A note about termination criteria: vol_tol refers to the volume of the hyperrectangle containing the lowest function value found so far. This volume decreases exponentially with increasing dimensionality of the problem. Therefore vol_tol should be decreased to avoid premature termination of the algorithm for higher dimensions. This does not hold for len_tol: it refers either to half of the maximal side length (for
locally_biased=True
) or half of the diagonal of the hyperrectangle (forlocally_biased=False
).This code is based on the DIRECT 2.0.4 Fortran code by Gablonsky et al. at https://ctk.math.ncsu.edu/SOFTWARE/DIRECTv204.tar.gz . This original version was initially converted via f2c and then cleaned up and reorganized by Steven G. Johnson, August 2007, for the NLopt project. The
direct
function wraps the C implementation.New in version 1.9.0.
References
- 1
Jones, D.R., Perttunen, C.D. & Stuckman, B.E. Lipschitzian optimization without the Lipschitz constant. J Optim Theory Appl 79, 157-181 (1993).
- 2
Gablonsky, J., Kelley, C. A Locally-Biased form of the DIRECT Algorithm. Journal of Global Optimization 21, 27-37 (2001).
Examples
The following example is a 2-D problem with four local minima: minimizing the Styblinski-Tang function (https://en.wikipedia.org/wiki/Test_functions_for_optimization).
>>> from scipy.optimize import direct, Bounds >>> def styblinski_tang(pos): ... x, y = pos ... return 0.5 * (x**4 - 16*x**2 + 5*x + y**4 - 16*y**2 + 5*y) >>> bounds = Bounds([-4., -4.], [4., 4.]) >>> result = direct(styblinski_tang, bounds) >>> result.x, result.fun, result.nfev array([-2.90321597, -2.90321597]), -78.3323279095383, 2011
The correct global minimum was found but with a huge number of function evaluations (2011). Loosening the termination tolerances vol_tol and len_tol can be used to stop DIRECT earlier.
>>> result = direct(styblinski_tang, bounds, len_tol=1e-3) >>> result.x, result.fun, result.nfev array([-2.9044353, -2.9044353]), -78.33230330754142, 207