Find the global minimum of a function using the basin-hopping algorithm
New in version 0.12.0.
Parameters : | func : callable f(x, *args)
x0 : ndarray
niter : integer, optional
T : float, optional
stepsize : float, optional
minimizer_kwargs : dict, optional
take_step : callable take_step(x), optional
accept_test : callable, accept_test(f_new=f_new, x_new=x_new, f_old=fold, x_old=x_old), optional
callback : callable, callback(x, f, accept), optional
interval : integer, optional
disp : bool, optional
niter_success : integer, optional
|
---|---|
Returns : | res : Result
|
See also
Notes
Basin-hopping is a stochastic algorithm which attempts to find the global minimum of a smooth scalar function of one or more variables [R72] [R73] [R74] [R75]. The algorithm in its current form was described by David Wales and Jonathan Doye [R73] http://www-wales.ch.cam.ac.uk/.
The algorithm is iterative with each cycle composed of the following features
The acceptance test used here is the Metropolis criterion of standard Monte Carlo algorithms, although there are many other possibilities [R74].
This global minimization method has been shown to be extremely efficient for a wide variety of problems in physics and chemistry. It is particularly useful when the function has many minima separated by large barriers. See the Cambridge Cluster Database http://www-wales.ch.cam.ac.uk/CCD.html for databases of molecular systems that have been optimized primarily using basin-hopping. This database includes minimization problems exceeding 300 degrees of freedom.
See the free software program GMIN (http://www-wales.ch.cam.ac.uk/GMIN) for a Fortran implementation of basin-hopping. This implementation has many different variations of the procedure described above, including more advanced step taking algorithms and alternate acceptance criterion.
For stochastic global optimization there is no way to determine if the true global minimum has actually been found. Instead, as a consistency check, the algorithm can be run from a number of different random starting points to ensure the lowest minimum found in each example has converged to the global minimum. For this reason basinhopping will by default simply run for the number of iterations niter and return the lowest minimum found. It is left to the user to ensure that this is in fact the global minimum.
Choosing stepsize: This is a crucial parameter in basinhopping and depends on the problem being solved. Ideally it should be comparable to the typical separation between local minima of the function being optimized. basinhopping will, by default, adjust stepsize to find an optimal value, but this may take many iterations. You will get quicker results if you set a sensible value for stepsize.
Choosing T: The parameter T is the temperature used in the metropolis criterion. Basinhopping steps are accepted with probability 1 if func(xnew) < func(xold), or otherwise with probability:
exp( -(func(xnew) - func(xold)) / T )
So, for best results, T should to be comparable to the typical difference in function value between between local minima
References
[R72] | (1, 2) Wales, David J. 2003, Energy Landscapes, Cambridge University Press, Cambridge, UK. |
[R73] | (1, 2, 3) Wales, D J, and Doye J P K, Global Optimization by Basin-Hopping and the Lowest Energy Structures of Lennard-Jones Clusters Containing up to 110 Atoms. Journal of Physical Chemistry A, 1997, 101, 5111. |
[R74] | (1, 2, 3) Li, Z. and Scheraga, H. A., Monte Carlo-minimization approach to the multiple-minima problem in protein folding, Proc. Natl. Acad. Sci. USA, 1987, 84, 6611. |
[R75] | (1, 2) Wales, D. J. and Scheraga, H. A., Global optimization of clusters, crystals, and biomolecules, Science, 1999, 285, 1368. |
Examples
The following example is a one-dimensional minimization problem, with many local minima superimposed on a parabola.
>>> func = lambda x: cos(14.5 * x - 0.3) + (x + 0.2) * x
>>> x0=[1.]
Basinhopping, internally, uses a local minimization algorithm. We will use the parameter minimizer_kwargs to tell basinhopping which algorithm to use and how to set up that minimizer. This parameter will be passed to scipy.optimize.minimize().
>>> minimizer_kwargs = {"method": "BFGS"}
>>> ret = basinhopping(func, x0, minimizer_kwargs=minimizer_kwargs,
... niter=200)
>>> print("global minimum: x = %.4f, f(x0) = %.4f" % (ret.x, ret.fun))
global minimum: x = -0.1951, f(x0) = -1.0009
Next consider a two-dimensional minimization problem. Also, this time we will use gradient information to significantly speed up the search.
>>> def func2d(x):
... f = cos(14.5 * x[0] - 0.3) + (x[1] + 0.2) * x[1] + (x[0] +
... 0.2) * x[0]
... df = np.zeros(2)
... df[0] = -14.5 * sin(14.5 * x[0] - 0.3) + 2. * x[0] + 0.2
... df[1] = 2. * x[1] + 0.2
... return f, df
We’ll also use a different local minimization algorithm. Also we must tell the minimizer that our function returns both energy and gradient (jacobian)
>>> minimizer_kwargs = {"method":"L-BFGS-B", "jac":True}
>>> x0 = [1.0, 1.0]
>>> ret = basinhopping(func2d, x0, minimizer_kwargs=minimizer_kwargs,
... niter=200)
>>> print("global minimum: x = [%.4f, %.4f], f(x0) = %.4f" % (ret.x[0],
... ret.x[1],
... ret.fun))
global minimum: x = [-0.1951, -0.1000], f(x0) = -1.0109
Here is an example using a custom step taking routine. Imagine you want the first coordinate to take larger steps then the rest of the coordinates. This can be implemented like so:
>>> class MyTakeStep(object):
... def __init__(self, stepsize=0.5):
... self.stepsize = stepsize
... def __call__(self, x):
... s = self.stepsize
... x[0] += np.random.uniform(-2.*s, 2.*s)
... x[1:] += np.random.uniform(-s, s, x[1:].shape)
... return x
Since MyTakeStep.stepsize exists basinhopping will adjust the magnitude of stepsize to optimize the search. We’ll use the same 2-D function as before
>>> mytakestep = MyTakeStep()
>>> ret = basinhopping(func2d, x0, minimizer_kwargs=minimizer_kwargs,
... niter=200, take_step=mytakestep)
>>> print("global minimum: x = [%.4f, %.4f], f(x0) = %.4f" % (ret.x[0],
... ret.x[1],
... ret.fun))
global minimum: x = [-0.1951, -0.1000], f(x0) = -1.0109
Now let’s do an example using a custom callback function which prints the value of every minimum found
>>> def print_fun(x, f, accepted):
... print("at minima %.4f accepted %d" % (f, int(accepted)))
We’ll run it for only 10 basinhopping steps this time.
>>> np.random.seed(1)
>>> ret = basinhopping(func2d, x0, minimizer_kwargs=minimizer_kwargs,
... niter=10, callback=print_fun)
at minima 0.4159 accepted 1
at minima -0.9073 accepted 1
at minima -0.1021 accepted 1
at minima -0.1021 accepted 1
at minima 0.9102 accepted 1
at minima 0.9102 accepted 1
at minima 2.2945 accepted 0
at minima -0.1021 accepted 1
at minima -1.0109 accepted 1
at minima -1.0109 accepted 1
The minima at -1.0109 is actually the global minimum, found already on the 8th iteration.
Now let’s implement bounds on the problem using a custom accept_test:
>>> class MyBounds(object):
... def __init__(self, xmax=[1.1,1.1], xmin=[-1.1,-1.1] ):
... self.xmax = np.array(xmax)
... self.xmin = np.array(xmin)
... def __call__(self, **kwargs):
... x = kwargs["x_new"]
... tmax = bool(np.all(x <= self.xmax))
... tmin = bool(np.all(x >= self.xmin))
... return tmax and tmin
>>> mybounds = MyBounds()
>>> ret = basinhopping(func2d, x0, minimizer_kwargs=minimizer_kwargs,
... niter=10, accept_test=mybounds)