Skip to content

Commit

Permalink
ENH: define and use Result to store results of optimization wrappers
Browse files Browse the repository at this point in the history
`Result` is a subclass of `dict` with attribute accessor.
It is now the only return from optimization wrappers. The solution is
available through the `x` attribute (previously named `solution` in
`InfoDict`).

All tests updated.
  • Loading branch information
dlax committed May 15, 2012
1 parent 1f92156 commit 9cc3d37
Show file tree
Hide file tree
Showing 10 changed files with 279 additions and 285 deletions.
40 changes: 18 additions & 22 deletions scipy/optimize/anneal.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,8 +4,8 @@

import numpy
from numpy import asarray, tan, exp, ones, squeeze, sign, \
all, log, sqrt, pi, shape, array, minimum, where
from numpy import random
all, log, sqrt, pi, shape, array, minimum, where, random
from optimize import Result

__all__ = ['anneal']

Expand Down Expand Up @@ -303,13 +303,13 @@ def anneal(func, x0, args=(), schedule='fast', full_output=0,
'dwell' : dwell,
'disp' : disp}

x, info = _minimize_anneal(func, x0, args, opts)
res = _minimize_anneal(func, x0, args, opts)

if full_output:
return x, info['fun'], info['T'], info['nfev'], info['nit'], \
info['accept'], info['status']
return res['x'], res['fun'], res['T'], res['nfev'], res['nit'], \
res['accept'], res['status']
else:
return x, info['status']
return res['x'], res['status']

def _minimize_anneal(func, x0, args=(), options=None):
"""
Expand Down Expand Up @@ -449,22 +449,18 @@ def _minimize_anneal(func, x0, args=(), options=None):
retval = 4
break

info = {'solution': best_state.x,
'fun' : best_state.cost,
'T' : schedule.T,
'nfev' : schedule.feval,
'nit' : iters,
'accept' : schedule.accepted,
'status' : retval,
'success' : retval <= 1}
info['message'] = {0: 'Points no longer changing',
1: 'Cooled to final temperature',
2: 'Maximum function evaluations',
3: 'Maximum cooling iterations reached',
4: 'Maximum accepted query locations reached',
5: 'Final point not the minimum amongst '
'encountered points'}[retval]
return best_state.x, info
result = Result(x=best_state.x, fun=best_state.cost,
T=schedule.T, nfev=schedule.feval, nit=iters,
accept=schedule.accepted, status=retval,
success=(retval <= 1),
message={0: 'Points no longer changing',
1: 'Cooled to final temperature',
2: 'Maximum function evaluations',
3: 'Maximum cooling iterations reached',
4: 'Maximum accepted query locations reached',
5: 'Final point not the minimum amongst '
'encountered points'}[retval])
return result



Expand Down
5 changes: 3 additions & 2 deletions scipy/optimize/cobyla.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,6 +9,7 @@
"""

from scipy.optimize import _cobyla
from optimize import Result
from numpy import copy
from warnings import warn

Expand Down Expand Up @@ -160,7 +161,7 @@ def fmin_cobyla(func, x0, cons, args=(), consargs=None, rhobeg=1.0, rhoend=1e-4,
'maxfev': maxfun}

return _minimize_cobyla(func, x0, args, constraints=con,
options=opts)[0]
options=opts)['x']

def _minimize_cobyla(fun, x0, args=(), constraints=(), options=None):
"""
Expand Down Expand Up @@ -233,7 +234,7 @@ def calcfc(x, con):
xopt = _cobyla.minimize(calcfc, m=m, x=copy(x0), rhobeg=rhobeg,
rhoend=rhoend, iprint=iprint, maxfun=maxfun)

return xopt, dict()
return Result(x=xopt)

if __name__ == '__main__':

Expand Down
27 changes: 11 additions & 16 deletions scipy/optimize/lbfgsb.py
Original file line number Diff line number Diff line change
Expand Up @@ -26,7 +26,7 @@

from numpy import array, asarray, float64, int32, zeros
import _lbfgsb
from optimize import approx_fprime, MemoizeJac
from optimize import approx_fprime, MemoizeJac, Result
from numpy.compat import asbytes

__all__ = ['fmin_l_bfgs_b']
Expand Down Expand Up @@ -160,13 +160,14 @@ def fmin_l_bfgs_b(func, x0, fprime=None, args=(),
'eps' : epsilon,
'maxfev': maxfun}

x, info = _minimize_lbfgsb(fun, x0, args=args, jac=jac, bounds=bounds,
options=opts)
d = {'grad': info['jac'],
'task': info['message'],
'funcalls': info['nfev'],
'warnflag': info['status']}
f = info['fun']
res = _minimize_lbfgsb(fun, x0, args=args, jac=jac, bounds=bounds,
options=opts)
d = {'grad': res['jac'],
'task': res['message'],
'funcalls': res['nfev'],
'warnflag': res['status']}
f = res['fun']
x = res['x']

return x, f, d

Expand Down Expand Up @@ -304,14 +305,8 @@ def func_and_grad(x):
'warnflag' : warnflag
}

info = {'fun': f,
'jac': g,
'nfev': n_function_evals,
'status': warnflag,
'message': task_str,
'solution': x,
'success': warnflag==0}
return x, info
return Result(fun=f, jac=g, nfev=n_function_evals, status=warnflag,
message=task_str, x=x, success=(warnflag==0))


if __name__ == '__main__':
Expand Down
Loading

0 comments on commit 9cc3d37

Please sign in to comment.