Skip to content

Commit

Permalink
mnt: streamlined kwargs for minimizer in dual_annealing
Browse files Browse the repository at this point in the history
This fixes so that all global optimization strategies uses
minimizer_kwargs as the arguments for the minimizer.
This makes it compatible with shgo and bassin_hopping.

The dual_annealing previously only had local_search_options
as the kwarg, which is not compatible with the other ones.
  • Loading branch information
zerothi committed Jan 15, 2021
1 parent 5a42f06 commit 0780b60
Show file tree
Hide file tree
Showing 2 changed files with 67 additions and 10 deletions.
23 changes: 18 additions & 5 deletions scipy/optimize/_dual_annealing.py
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,8 @@
A Dual Annealing global optimization algorithm
"""

import warnings

import numpy as np
from scipy.optimize import OptimizeResult
from scipy.optimize import minimize
Expand Down Expand Up @@ -425,10 +427,10 @@ def local_search(self, x, e):


def dual_annealing(func, bounds, args=(), maxiter=1000,
local_search_options={}, initial_temp=5230.,
minimizer_kwargs={}, initial_temp=5230.,
restart_temp_ratio=2.e-5, visit=2.62, accept=-5.0,
maxfun=1e7, seed=None, no_local_search=False,
callback=None, x0=None):
callback=None, x0=None, local_search_options={}):
"""
Find the global minimum of a function using Dual Annealing.
Expand All @@ -447,7 +449,7 @@ def dual_annealing(func, bounds, args=(), maxiter=1000,
objective function.
maxiter : int, optional
The maximum number of global search iterations. Default value is 1000.
local_search_options : dict, optional
minimizer_kwargs : dict, optional
Extra keyword arguments to be passed to the local minimizer
(`minimize`). Some important options could be:
``method`` for the minimizer method to use and ``args`` for
Expand Down Expand Up @@ -503,6 +505,9 @@ def dual_annealing(func, bounds, args=(), maxiter=1000,
If the callback implementation returns True, the algorithm will stop.
x0 : ndarray, shape(n,), optional
Coordinates of a single N-D starting point.
local_search_options : dict, optional
Backwards compatible flag for `minimizer_kwargs`, only one of these
should be supplied.
Returns
-------
Expand Down Expand Up @@ -622,9 +627,17 @@ def dual_annealing(func, bounds, args=(), maxiter=1000,

# Wrapper for the objective function
func_wrapper = ObjectiveFunWrapper(func, maxfun, *args)
# Wrapper fot the minimizer
# Wrapper for the minimizer
if local_search_options and minimizer_kwargs:
raise ValueError("dual_annealing only allows either 'minimizer_kwargs' (preferred) or "
"'local_search_options' (deprecated); not both!")
if local_search_options:
warnings.warn("dual_annealing argument 'local_search_options' is "
"deprecated in favor of 'minimizer_kwargs'",
category=DeprecationWarning, stacklevel=2)
minimizer_kwargs = local_search_options
minimizer_wrapper = LocalSearchWrapper(
bounds, func_wrapper, **local_search_options)
bounds, func_wrapper, **minimizer_kwargs)
# Initialization of RandomState for reproducible runs if seed provided
rand_state = check_random_state(seed)
# Initialization of the energy state
Expand Down
54 changes: 49 additions & 5 deletions scipy/optimize/tests/test__dual_annealing.py
Original file line number Diff line number Diff line change
Expand Up @@ -174,22 +174,66 @@ def test_bound_validity(self):
assert_raises(ValueError, dual_annealing, self.func,
invalid_bounds)

def test_local_search_option_bounds(self):
def test_minimizer_kwargs_bounds(self):
func = lambda x: np.sum((x-5) * (x-1))
bounds = list(zip([-6, -5], [6, 5]))
# Test bounds can be passed (see gh-10831)
dual_annealing(
func,
bounds=bounds,
local_search_options={"method": "SLSQP", "bounds": bounds})
minimizer_kwargs={"method": "SLSQP", "bounds": bounds})

with np.testing.suppress_warnings() as sup:
sup.record(RuntimeWarning, "Method CG cannot handle ")

dual_annealing(
func,
bounds=bounds,
local_search_options={"method": "CG", "bounds": bounds})
minimizer_kwargs={"method": "CG", "bounds": bounds})

# Verify warning happened for Method cannot handle bounds.
assert sup.log


def test_deprecated_local_search_options_bounds(self):
func = lambda x: np.sum((x-5) * (x-1))
bounds = list(zip([-6, -5], [6, 5]))
# Test bounds can be passed (see gh-10831)
with np.testing.suppress_warnings() as sup:
sup.record(DeprecationWarning, "dual_annealing argument 'local_search_options'")
dual_annealing(
func,
bounds=bounds,
local_search_options={"method": "SLSQP", "bounds": bounds})

with np.testing.suppress_warnings() as sup:
sup.record(RuntimeWarning, "Method CG cannot handle ")

dual_annealing(
func,
bounds=bounds,
minimizer_kwargs={"method": "CG", "bounds": bounds})

# Verify warning happened for Method cannot handle bounds.
assert sup.log


def test_minimizer_kwargs_bounds(self):
func = lambda x: np.sum((x-5) * (x-1))
bounds = list(zip([-6, -5], [6, 5]))
# Test bounds can be passed (see gh-10831)
dual_annealing(
func,
bounds=bounds,
minimizer_kwargs={"method": "SLSQP", "bounds": bounds})

with np.testing.suppress_warnings() as sup:
sup.record(RuntimeWarning, "Method CG cannot handle ")

dual_annealing(
func,
bounds=bounds,
minimizer_kwargs={"method": "CG", "bounds": bounds})

# Verify warning happened for Method cannot handle bounds.
assert sup.log
Expand Down Expand Up @@ -249,7 +293,7 @@ def test_callback_stop(self):
])
def test_multi_ls_minimizer(self, method, atol):
ret = dual_annealing(self.func, self.ld_bounds,
local_search_options=dict(method=method),
minimizer_kwargs=dict(method=method),
seed=self.seed)
assert_allclose(ret.fun, 0., atol=atol)

Expand All @@ -264,7 +308,7 @@ def test_gradient_gnev(self):
'jac': self.rosen_der_wrapper,
}
ret = dual_annealing(rosen, self.ld_bounds,
local_search_options=minimizer_opts,
minimizer_kwargs=minimizer_opts,
seed=self.seed)
assert ret.njev == self.ngev

Expand Down

0 comments on commit 0780b60

Please sign in to comment.