_minimize.py 36 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820
  1. """
  2. Unified interfaces to minimization algorithms.
  3. Functions
  4. ---------
  5. - minimize : minimization of a function of several variables.
  6. - minimize_scalar : minimization of a function of one variable.
  7. """
  8. from __future__ import division, print_function, absolute_import
  9. __all__ = ['minimize', 'minimize_scalar']
  10. from warnings import warn
  11. import numpy as np
  12. from scipy._lib.six import callable
  13. # unconstrained minimization
  14. from .optimize import (_minimize_neldermead, _minimize_powell, _minimize_cg,
  15. _minimize_bfgs, _minimize_newtoncg,
  16. _minimize_scalar_brent, _minimize_scalar_bounded,
  17. _minimize_scalar_golden, MemoizeJac)
  18. from ._trustregion_dogleg import _minimize_dogleg
  19. from ._trustregion_ncg import _minimize_trust_ncg
  20. from ._trustregion_krylov import _minimize_trust_krylov
  21. from ._trustregion_exact import _minimize_trustregion_exact
  22. from ._trustregion_constr import _minimize_trustregion_constr
  23. # constrained minimization
  24. from .lbfgsb import _minimize_lbfgsb
  25. from .tnc import _minimize_tnc
  26. from .cobyla import _minimize_cobyla
  27. from .slsqp import _minimize_slsqp
  28. from ._constraints import (old_bound_to_new, new_bounds_to_old,
  29. old_constraint_to_new, new_constraint_to_old,
  30. NonlinearConstraint, LinearConstraint, Bounds)
  31. def minimize(fun, x0, args=(), method=None, jac=None, hess=None,
  32. hessp=None, bounds=None, constraints=(), tol=None,
  33. callback=None, options=None):
  34. """Minimization of scalar function of one or more variables.
  35. Parameters
  36. ----------
  37. fun : callable
  38. The objective function to be minimized.
  39. ``fun(x, *args) -> float``
  40. where x is an 1-D array with shape (n,) and `args`
  41. is a tuple of the fixed parameters needed to completely
  42. specify the function.
  43. x0 : ndarray, shape (n,)
  44. Initial guess. Array of real elements of size (n,),
  45. where 'n' is the number of independent variables.
  46. args : tuple, optional
  47. Extra arguments passed to the objective function and its
  48. derivatives (`fun`, `jac` and `hess` functions).
  49. method : str or callable, optional
  50. Type of solver. Should be one of
  51. - 'Nelder-Mead' :ref:`(see here) <optimize.minimize-neldermead>`
  52. - 'Powell' :ref:`(see here) <optimize.minimize-powell>`
  53. - 'CG' :ref:`(see here) <optimize.minimize-cg>`
  54. - 'BFGS' :ref:`(see here) <optimize.minimize-bfgs>`
  55. - 'Newton-CG' :ref:`(see here) <optimize.minimize-newtoncg>`
  56. - 'L-BFGS-B' :ref:`(see here) <optimize.minimize-lbfgsb>`
  57. - 'TNC' :ref:`(see here) <optimize.minimize-tnc>`
  58. - 'COBYLA' :ref:`(see here) <optimize.minimize-cobyla>`
  59. - 'SLSQP' :ref:`(see here) <optimize.minimize-slsqp>`
  60. - 'trust-constr':ref:`(see here) <optimize.minimize-trustconstr>`
  61. - 'dogleg' :ref:`(see here) <optimize.minimize-dogleg>`
  62. - 'trust-ncg' :ref:`(see here) <optimize.minimize-trustncg>`
  63. - 'trust-exact' :ref:`(see here) <optimize.minimize-trustexact>`
  64. - 'trust-krylov' :ref:`(see here) <optimize.minimize-trustkrylov>`
  65. - custom - a callable object (added in version 0.14.0),
  66. see below for description.
  67. If not given, chosen to be one of ``BFGS``, ``L-BFGS-B``, ``SLSQP``,
  68. depending if the problem has constraints or bounds.
  69. jac : {callable, '2-point', '3-point', 'cs', bool}, optional
  70. Method for computing the gradient vector. Only for CG, BFGS,
  71. Newton-CG, L-BFGS-B, TNC, SLSQP, dogleg, trust-ncg, trust-krylov,
  72. trust-exact and trust-constr. If it is a callable, it should be a
  73. function that returns the gradient vector:
  74. ``jac(x, *args) -> array_like, shape (n,)``
  75. where x is an array with shape (n,) and `args` is a tuple with
  76. the fixed parameters. Alternatively, the keywords
  77. {'2-point', '3-point', 'cs'} select a finite
  78. difference scheme for numerical estimation of the gradient. Options
  79. '3-point' and 'cs' are available only to 'trust-constr'.
  80. If `jac` is a Boolean and is True, `fun` is assumed to return the
  81. gradient along with the objective function. If False, the gradient
  82. will be estimated using '2-point' finite difference estimation.
  83. hess : {callable, '2-point', '3-point', 'cs', HessianUpdateStrategy}, optional
  84. Method for computing the Hessian matrix. Only for Newton-CG, dogleg,
  85. trust-ncg, trust-krylov, trust-exact and trust-constr. If it is
  86. callable, it should return the Hessian matrix:
  87. ``hess(x, *args) -> {LinearOperator, spmatrix, array}, (n, n)``
  88. where x is a (n,) ndarray and `args` is a tuple with the fixed
  89. parameters. LinearOperator and sparse matrix returns are
  90. allowed only for 'trust-constr' method. Alternatively, the keywords
  91. {'2-point', '3-point', 'cs'} select a finite difference scheme
  92. for numerical estimation. Or, objects implementing
  93. `HessianUpdateStrategy` interface can be used to approximate
  94. the Hessian. Available quasi-Newton methods implementing
  95. this interface are:
  96. - `BFGS`;
  97. - `SR1`.
  98. Whenever the gradient is estimated via finite-differences,
  99. the Hessian cannot be estimated with options
  100. {'2-point', '3-point', 'cs'} and needs to be
  101. estimated using one of the quasi-Newton strategies.
  102. Finite-difference options {'2-point', '3-point', 'cs'} and
  103. `HessianUpdateStrategy` are available only for 'trust-constr' method.
  104. hessp : callable, optional
  105. Hessian of objective function times an arbitrary vector p. Only for
  106. Newton-CG, trust-ncg, trust-krylov, trust-constr.
  107. Only one of `hessp` or `hess` needs to be given. If `hess` is
  108. provided, then `hessp` will be ignored. `hessp` must compute the
  109. Hessian times an arbitrary vector:
  110. ``hessp(x, p, *args) -> ndarray shape (n,)``
  111. where x is a (n,) ndarray, p is an arbitrary vector with
  112. dimension (n,) and `args` is a tuple with the fixed
  113. parameters.
  114. bounds : sequence or `Bounds`, optional
  115. Bounds on variables for L-BFGS-B, TNC, SLSQP and
  116. trust-constr methods. There are two ways to specify the bounds:
  117. 1. Instance of `Bounds` class.
  118. 2. Sequence of ``(min, max)`` pairs for each element in `x`. None
  119. is used to specify no bound.
  120. constraints : {Constraint, dict} or List of {Constraint, dict}, optional
  121. Constraints definition (only for COBYLA, SLSQP and trust-constr).
  122. Constraints for 'trust-constr' are defined as a single object or a
  123. list of objects specifying constraints to the optimization problem.
  124. Available constraints are:
  125. - `LinearConstraint`
  126. - `NonlinearConstraint`
  127. Constraints for COBYLA, SLSQP are defined as a list of dictionaries.
  128. Each dictionary with fields:
  129. type : str
  130. Constraint type: 'eq' for equality, 'ineq' for inequality.
  131. fun : callable
  132. The function defining the constraint.
  133. jac : callable, optional
  134. The Jacobian of `fun` (only for SLSQP).
  135. args : sequence, optional
  136. Extra arguments to be passed to the function and Jacobian.
  137. Equality constraint means that the constraint function result is to
  138. be zero whereas inequality means that it is to be non-negative.
  139. Note that COBYLA only supports inequality constraints.
  140. tol : float, optional
  141. Tolerance for termination. For detailed control, use solver-specific
  142. options.
  143. options : dict, optional
  144. A dictionary of solver options. All methods accept the following
  145. generic options:
  146. maxiter : int
  147. Maximum number of iterations to perform.
  148. disp : bool
  149. Set to True to print convergence messages.
  150. For method-specific options, see :func:`show_options()`.
  151. callback : callable, optional
  152. Called after each iteration. For 'trust-constr' it is a callable with
  153. the signature:
  154. ``callback(xk, OptimizeResult state) -> bool``
  155. where ``xk`` is the current parameter vector. and ``state``
  156. is an `OptimizeResult` object, with the same fields
  157. as the ones from the return. If callback returns True
  158. the algorithm execution is terminated.
  159. For all the other methods, the signature is:
  160. ``callback(xk)``
  161. where ``xk`` is the current parameter vector.
  162. Returns
  163. -------
  164. res : OptimizeResult
  165. The optimization result represented as a ``OptimizeResult`` object.
  166. Important attributes are: ``x`` the solution array, ``success`` a
  167. Boolean flag indicating if the optimizer exited successfully and
  168. ``message`` which describes the cause of the termination. See
  169. `OptimizeResult` for a description of other attributes.
  170. See also
  171. --------
  172. minimize_scalar : Interface to minimization algorithms for scalar
  173. univariate functions
  174. show_options : Additional options accepted by the solvers
  175. Notes
  176. -----
  177. This section describes the available solvers that can be selected by the
  178. 'method' parameter. The default method is *BFGS*.
  179. **Unconstrained minimization**
  180. Method :ref:`Nelder-Mead <optimize.minimize-neldermead>` uses the
  181. Simplex algorithm [1]_, [2]_. This algorithm is robust in many
  182. applications. However, if numerical computation of derivative can be
  183. trusted, other algorithms using the first and/or second derivatives
  184. information might be preferred for their better performance in
  185. general.
  186. Method :ref:`Powell <optimize.minimize-powell>` is a modification
  187. of Powell's method [3]_, [4]_ which is a conjugate direction
  188. method. It performs sequential one-dimensional minimizations along
  189. each vector of the directions set (`direc` field in `options` and
  190. `info`), which is updated at each iteration of the main
  191. minimization loop. The function need not be differentiable, and no
  192. derivatives are taken.
  193. Method :ref:`CG <optimize.minimize-cg>` uses a nonlinear conjugate
  194. gradient algorithm by Polak and Ribiere, a variant of the
  195. Fletcher-Reeves method described in [5]_ pp. 120-122. Only the
  196. first derivatives are used.
  197. Method :ref:`BFGS <optimize.minimize-bfgs>` uses the quasi-Newton
  198. method of Broyden, Fletcher, Goldfarb, and Shanno (BFGS) [5]_
  199. pp. 136. It uses the first derivatives only. BFGS has proven good
  200. performance even for non-smooth optimizations. This method also
  201. returns an approximation of the Hessian inverse, stored as
  202. `hess_inv` in the OptimizeResult object.
  203. Method :ref:`Newton-CG <optimize.minimize-newtoncg>` uses a
  204. Newton-CG algorithm [5]_ pp. 168 (also known as the truncated
  205. Newton method). It uses a CG method to the compute the search
  206. direction. See also *TNC* method for a box-constrained
  207. minimization with a similar algorithm. Suitable for large-scale
  208. problems.
  209. Method :ref:`dogleg <optimize.minimize-dogleg>` uses the dog-leg
  210. trust-region algorithm [5]_ for unconstrained minimization. This
  211. algorithm requires the gradient and Hessian; furthermore the
  212. Hessian is required to be positive definite.
  213. Method :ref:`trust-ncg <optimize.minimize-trustncg>` uses the
  214. Newton conjugate gradient trust-region algorithm [5]_ for
  215. unconstrained minimization. This algorithm requires the gradient
  216. and either the Hessian or a function that computes the product of
  217. the Hessian with a given vector. Suitable for large-scale problems.
  218. Method :ref:`trust-krylov <optimize.minimize-trustkrylov>` uses
  219. the Newton GLTR trust-region algorithm [14]_, [15]_ for unconstrained
  220. minimization. This algorithm requires the gradient
  221. and either the Hessian or a function that computes the product of
  222. the Hessian with a given vector. Suitable for large-scale problems.
  223. On indefinite problems it requires usually less iterations than the
  224. `trust-ncg` method and is recommended for medium and large-scale problems.
  225. Method :ref:`trust-exact <optimize.minimize-trustexact>`
  226. is a trust-region method for unconstrained minimization in which
  227. quadratic subproblems are solved almost exactly [13]_. This
  228. algorithm requires the gradient and the Hessian (which is
  229. *not* required to be positive definite). It is, in many
  230. situations, the Newton method to converge in fewer iteraction
  231. and the most recommended for small and medium-size problems.
  232. **Bound-Constrained minimization**
  233. Method :ref:`L-BFGS-B <optimize.minimize-lbfgsb>` uses the L-BFGS-B
  234. algorithm [6]_, [7]_ for bound constrained minimization.
  235. Method :ref:`TNC <optimize.minimize-tnc>` uses a truncated Newton
  236. algorithm [5]_, [8]_ to minimize a function with variables subject
  237. to bounds. This algorithm uses gradient information; it is also
  238. called Newton Conjugate-Gradient. It differs from the *Newton-CG*
  239. method described above as it wraps a C implementation and allows
  240. each variable to be given upper and lower bounds.
  241. **Constrained Minimization**
  242. Method :ref:`COBYLA <optimize.minimize-cobyla>` uses the
  243. Constrained Optimization BY Linear Approximation (COBYLA) method
  244. [9]_, [10]_, [11]_. The algorithm is based on linear
  245. approximations to the objective function and each constraint. The
  246. method wraps a FORTRAN implementation of the algorithm. The
  247. constraints functions 'fun' may return either a single number
  248. or an array or list of numbers.
  249. Method :ref:`SLSQP <optimize.minimize-slsqp>` uses Sequential
  250. Least SQuares Programming to minimize a function of several
  251. variables with any combination of bounds, equality and inequality
  252. constraints. The method wraps the SLSQP Optimization subroutine
  253. originally implemented by Dieter Kraft [12]_. Note that the
  254. wrapper handles infinite values in bounds by converting them into
  255. large floating values.
  256. Method :ref:`trust-constr <optimize.minimize-trustconstr>` is a
  257. trust-region algorithm for constrained optimization. It swiches
  258. between two implementations depending on the problem definition.
  259. It is the most versatile constrained minimization algorithm
  260. implemented in SciPy and the most appropriate for large-scale problems.
  261. For equality constrained problems it is an implementation of Byrd-Omojokun
  262. Trust-Region SQP method described in [17]_ and in [5]_, p. 549. When
  263. inequality constraints are imposed as well, it swiches to the trust-region
  264. interior point method described in [16]_. This interior point algorithm,
  265. in turn, solves inequality constraints by introducing slack variables
  266. and solving a sequence of equality-constrained barrier problems
  267. for progressively smaller values of the barrier parameter.
  268. The previously described equality constrained SQP method is
  269. used to solve the subproblems with increasing levels of accuracy
  270. as the iterate gets closer to a solution.
  271. **Finite-Difference Options**
  272. For Method :ref:`trust-constr <optimize.minimize-trustconstr>`
  273. the gradient and the Hessian may be approximated using
  274. three finite-difference schemes: {'2-point', '3-point', 'cs'}.
  275. The scheme 'cs' is, potentially, the most accurate but it
  276. requires the function to correctly handles complex inputs and to
  277. be differentiable in the complex plane. The scheme '3-point' is more
  278. accurate than '2-point' but requires twice as much operations.
  279. **Custom minimizers**
  280. It may be useful to pass a custom minimization method, for example
  281. when using a frontend to this method such as `scipy.optimize.basinhopping`
  282. or a different library. You can simply pass a callable as the ``method``
  283. parameter.
  284. The callable is called as ``method(fun, x0, args, **kwargs, **options)``
  285. where ``kwargs`` corresponds to any other parameters passed to `minimize`
  286. (such as `callback`, `hess`, etc.), except the `options` dict, which has
  287. its contents also passed as `method` parameters pair by pair. Also, if
  288. `jac` has been passed as a bool type, `jac` and `fun` are mangled so that
  289. `fun` returns just the function values and `jac` is converted to a function
  290. returning the Jacobian. The method shall return an ``OptimizeResult``
  291. object.
  292. The provided `method` callable must be able to accept (and possibly ignore)
  293. arbitrary parameters; the set of parameters accepted by `minimize` may
  294. expand in future versions and then these parameters will be passed to
  295. the method. You can find an example in the scipy.optimize tutorial.
  296. .. versionadded:: 0.11.0
  297. References
  298. ----------
  299. .. [1] Nelder, J A, and R Mead. 1965. A Simplex Method for Function
  300. Minimization. The Computer Journal 7: 308-13.
  301. .. [2] Wright M H. 1996. Direct search methods: Once scorned, now
  302. respectable, in Numerical Analysis 1995: Proceedings of the 1995
  303. Dundee Biennial Conference in Numerical Analysis (Eds. D F
  304. Griffiths and G A Watson). Addison Wesley Longman, Harlow, UK.
  305. 191-208.
  306. .. [3] Powell, M J D. 1964. An efficient method for finding the minimum of
  307. a function of several variables without calculating derivatives. The
  308. Computer Journal 7: 155-162.
  309. .. [4] Press W, S A Teukolsky, W T Vetterling and B P Flannery.
  310. Numerical Recipes (any edition), Cambridge University Press.
  311. .. [5] Nocedal, J, and S J Wright. 2006. Numerical Optimization.
  312. Springer New York.
  313. .. [6] Byrd, R H and P Lu and J. Nocedal. 1995. A Limited Memory
  314. Algorithm for Bound Constrained Optimization. SIAM Journal on
  315. Scientific and Statistical Computing 16 (5): 1190-1208.
  316. .. [7] Zhu, C and R H Byrd and J Nocedal. 1997. L-BFGS-B: Algorithm
  317. 778: L-BFGS-B, FORTRAN routines for large scale bound constrained
  318. optimization. ACM Transactions on Mathematical Software 23 (4):
  319. 550-560.
  320. .. [8] Nash, S G. Newton-Type Minimization Via the Lanczos Method.
  321. 1984. SIAM Journal of Numerical Analysis 21: 770-778.
  322. .. [9] Powell, M J D. A direct search optimization method that models
  323. the objective and constraint functions by linear interpolation.
  324. 1994. Advances in Optimization and Numerical Analysis, eds. S. Gomez
  325. and J-P Hennart, Kluwer Academic (Dordrecht), 51-67.
  326. .. [10] Powell M J D. Direct search algorithms for optimization
  327. calculations. 1998. Acta Numerica 7: 287-336.
  328. .. [11] Powell M J D. A view of algorithms for optimization without
  329. derivatives. 2007.Cambridge University Technical Report DAMTP
  330. 2007/NA03
  331. .. [12] Kraft, D. A software package for sequential quadratic
  332. programming. 1988. Tech. Rep. DFVLR-FB 88-28, DLR German Aerospace
  333. Center -- Institute for Flight Mechanics, Koln, Germany.
  334. .. [13] Conn, A. R., Gould, N. I., and Toint, P. L.
  335. Trust region methods. 2000. Siam. pp. 169-200.
  336. .. [14] F. Lenders, C. Kirches, A. Potschka: "trlib: A vector-free
  337. implementation of the GLTR method for iterative solution of
  338. the trust region problem", https://arxiv.org/abs/1611.04718
  339. .. [15] N. Gould, S. Lucidi, M. Roma, P. Toint: "Solving the
  340. Trust-Region Subproblem using the Lanczos Method",
  341. SIAM J. Optim., 9(2), 504--525, (1999).
  342. .. [16] Byrd, Richard H., Mary E. Hribar, and Jorge Nocedal. 1999.
  343. An interior point algorithm for large-scale nonlinear programming.
  344. SIAM Journal on Optimization 9.4: 877-900.
  345. .. [17] Lalee, Marucha, Jorge Nocedal, and Todd Plantega. 1998. On the
  346. implementation of an algorithm for large-scale equality constrained
  347. optimization. SIAM Journal on Optimization 8.3: 682-706.
  348. Examples
  349. --------
  350. Let us consider the problem of minimizing the Rosenbrock function. This
  351. function (and its respective derivatives) is implemented in `rosen`
  352. (resp. `rosen_der`, `rosen_hess`) in the `scipy.optimize`.
  353. >>> from scipy.optimize import minimize, rosen, rosen_der
  354. A simple application of the *Nelder-Mead* method is:
  355. >>> x0 = [1.3, 0.7, 0.8, 1.9, 1.2]
  356. >>> res = minimize(rosen, x0, method='Nelder-Mead', tol=1e-6)
  357. >>> res.x
  358. array([ 1., 1., 1., 1., 1.])
  359. Now using the *BFGS* algorithm, using the first derivative and a few
  360. options:
  361. >>> res = minimize(rosen, x0, method='BFGS', jac=rosen_der,
  362. ... options={'gtol': 1e-6, 'disp': True})
  363. Optimization terminated successfully.
  364. Current function value: 0.000000
  365. Iterations: 26
  366. Function evaluations: 31
  367. Gradient evaluations: 31
  368. >>> res.x
  369. array([ 1., 1., 1., 1., 1.])
  370. >>> print(res.message)
  371. Optimization terminated successfully.
  372. >>> res.hess_inv
  373. array([[ 0.00749589, 0.01255155, 0.02396251, 0.04750988, 0.09495377], # may vary
  374. [ 0.01255155, 0.02510441, 0.04794055, 0.09502834, 0.18996269],
  375. [ 0.02396251, 0.04794055, 0.09631614, 0.19092151, 0.38165151],
  376. [ 0.04750988, 0.09502834, 0.19092151, 0.38341252, 0.7664427 ],
  377. [ 0.09495377, 0.18996269, 0.38165151, 0.7664427, 1.53713523]])
  378. Next, consider a minimization problem with several constraints (namely
  379. Example 16.4 from [5]_). The objective function is:
  380. >>> fun = lambda x: (x[0] - 1)**2 + (x[1] - 2.5)**2
  381. There are three constraints defined as:
  382. >>> cons = ({'type': 'ineq', 'fun': lambda x: x[0] - 2 * x[1] + 2},
  383. ... {'type': 'ineq', 'fun': lambda x: -x[0] - 2 * x[1] + 6},
  384. ... {'type': 'ineq', 'fun': lambda x: -x[0] + 2 * x[1] + 2})
  385. And variables must be positive, hence the following bounds:
  386. >>> bnds = ((0, None), (0, None))
  387. The optimization problem is solved using the SLSQP method as:
  388. >>> res = minimize(fun, (2, 0), method='SLSQP', bounds=bnds,
  389. ... constraints=cons)
  390. It should converge to the theoretical solution (1.4 ,1.7).
  391. """
  392. x0 = np.asarray(x0)
  393. if x0.dtype.kind in np.typecodes["AllInteger"]:
  394. x0 = np.asarray(x0, dtype=float)
  395. if not isinstance(args, tuple):
  396. args = (args,)
  397. if method is None:
  398. # Select automatically
  399. if constraints:
  400. method = 'SLSQP'
  401. elif bounds is not None:
  402. method = 'L-BFGS-B'
  403. else:
  404. method = 'BFGS'
  405. if callable(method):
  406. meth = "_custom"
  407. else:
  408. meth = method.lower()
  409. if options is None:
  410. options = {}
  411. # check if optional parameters are supported by the selected method
  412. # - jac
  413. if meth in ('nelder-mead', 'powell', 'cobyla') and bool(jac):
  414. warn('Method %s does not use gradient information (jac).' % method,
  415. RuntimeWarning)
  416. # - hess
  417. if meth not in ('newton-cg', 'dogleg', 'trust-ncg', 'trust-constr',
  418. 'trust-krylov', 'trust-exact', '_custom') and hess is not None:
  419. warn('Method %s does not use Hessian information (hess).' % method,
  420. RuntimeWarning)
  421. # - hessp
  422. if meth not in ('newton-cg', 'dogleg', 'trust-ncg', 'trust-constr',
  423. 'trust-krylov', '_custom') \
  424. and hessp is not None:
  425. warn('Method %s does not use Hessian-vector product '
  426. 'information (hessp).' % method, RuntimeWarning)
  427. # - constraints or bounds
  428. if (meth in ('nelder-mead', 'powell', 'cg', 'bfgs', 'newton-cg', 'dogleg',
  429. 'trust-ncg') and (bounds is not None or np.any(constraints))):
  430. warn('Method %s cannot handle constraints nor bounds.' % method,
  431. RuntimeWarning)
  432. if meth in ('l-bfgs-b', 'tnc') and np.any(constraints):
  433. warn('Method %s cannot handle constraints.' % method,
  434. RuntimeWarning)
  435. if meth == 'cobyla' and bounds is not None:
  436. warn('Method %s cannot handle bounds.' % method,
  437. RuntimeWarning)
  438. # - callback
  439. if (meth in ('cobyla',) and callback is not None):
  440. warn('Method %s does not support callback.' % method, RuntimeWarning)
  441. # - return_all
  442. if (meth in ('l-bfgs-b', 'tnc', 'cobyla', 'slsqp') and
  443. options.get('return_all', False)):
  444. warn('Method %s does not support the return_all option.' % method,
  445. RuntimeWarning)
  446. # check gradient vector
  447. if meth == 'trust-constr':
  448. if type(jac) is bool:
  449. if jac:
  450. fun = MemoizeJac(fun)
  451. jac = fun.derivative
  452. else:
  453. jac = '2-point'
  454. elif jac is None:
  455. jac = '2-point'
  456. elif not callable(jac) and jac not in ('2-point', '3-point', 'cs'):
  457. raise ValueError("Unsupported jac definition.")
  458. else:
  459. if jac in ('2-point', '3-point', 'cs'):
  460. if jac in ('3-point', 'cs'):
  461. warn("Only 'trust-constr' method accept %s "
  462. "options for 'jac'. Using '2-point' instead." % jac)
  463. jac = None
  464. elif not callable(jac):
  465. if bool(jac):
  466. fun = MemoizeJac(fun)
  467. jac = fun.derivative
  468. else:
  469. jac = None
  470. # set default tolerances
  471. if tol is not None:
  472. options = dict(options)
  473. if meth == 'nelder-mead':
  474. options.setdefault('xatol', tol)
  475. options.setdefault('fatol', tol)
  476. if meth in ('newton-cg', 'powell', 'tnc'):
  477. options.setdefault('xtol', tol)
  478. if meth in ('powell', 'l-bfgs-b', 'tnc', 'slsqp'):
  479. options.setdefault('ftol', tol)
  480. if meth in ('bfgs', 'cg', 'l-bfgs-b', 'tnc', 'dogleg',
  481. 'trust-ncg', 'trust-exact', 'trust-krylov'):
  482. options.setdefault('gtol', tol)
  483. if meth in ('cobyla', '_custom'):
  484. options.setdefault('tol', tol)
  485. if meth == 'trust-constr':
  486. options.setdefault('xtol', tol)
  487. options.setdefault('gtol', tol)
  488. options.setdefault('barrier_tol', tol)
  489. if bounds is not None:
  490. bounds = standardize_bounds(bounds, x0, meth)
  491. if constraints is not None:
  492. constraints = standardize_constraints(constraints, x0, meth)
  493. if meth == '_custom':
  494. return method(fun, x0, args=args, jac=jac, hess=hess, hessp=hessp,
  495. bounds=bounds, constraints=constraints,
  496. callback=callback, **options)
  497. elif meth == 'nelder-mead':
  498. return _minimize_neldermead(fun, x0, args, callback, **options)
  499. elif meth == 'powell':
  500. return _minimize_powell(fun, x0, args, callback, **options)
  501. elif meth == 'cg':
  502. return _minimize_cg(fun, x0, args, jac, callback, **options)
  503. elif meth == 'bfgs':
  504. return _minimize_bfgs(fun, x0, args, jac, callback, **options)
  505. elif meth == 'newton-cg':
  506. return _minimize_newtoncg(fun, x0, args, jac, hess, hessp, callback,
  507. **options)
  508. elif meth == 'l-bfgs-b':
  509. return _minimize_lbfgsb(fun, x0, args, jac, bounds,
  510. callback=callback, **options)
  511. elif meth == 'tnc':
  512. return _minimize_tnc(fun, x0, args, jac, bounds, callback=callback,
  513. **options)
  514. elif meth == 'cobyla':
  515. return _minimize_cobyla(fun, x0, args, constraints, **options)
  516. elif meth == 'slsqp':
  517. return _minimize_slsqp(fun, x0, args, jac, bounds,
  518. constraints, callback=callback, **options)
  519. elif meth == 'trust-constr':
  520. return _minimize_trustregion_constr(fun, x0, args, jac, hess, hessp,
  521. bounds, constraints,
  522. callback=callback, **options)
  523. elif meth == 'dogleg':
  524. return _minimize_dogleg(fun, x0, args, jac, hess,
  525. callback=callback, **options)
  526. elif meth == 'trust-ncg':
  527. return _minimize_trust_ncg(fun, x0, args, jac, hess, hessp,
  528. callback=callback, **options)
  529. elif meth == 'trust-krylov':
  530. return _minimize_trust_krylov(fun, x0, args, jac, hess, hessp,
  531. callback=callback, **options)
  532. elif meth == 'trust-exact':
  533. return _minimize_trustregion_exact(fun, x0, args, jac, hess,
  534. callback=callback, **options)
  535. else:
  536. raise ValueError('Unknown solver %s' % method)
  537. def minimize_scalar(fun, bracket=None, bounds=None, args=(),
  538. method='brent', tol=None, options=None):
  539. """Minimization of scalar function of one variable.
  540. Parameters
  541. ----------
  542. fun : callable
  543. Objective function.
  544. Scalar function, must return a scalar.
  545. bracket : sequence, optional
  546. For methods 'brent' and 'golden', `bracket` defines the bracketing
  547. interval and can either have three items ``(a, b, c)`` so that
  548. ``a < b < c`` and ``fun(b) < fun(a), fun(c)`` or two items ``a`` and
  549. ``c`` which are assumed to be a starting interval for a downhill
  550. bracket search (see `bracket`); it doesn't always mean that the
  551. obtained solution will satisfy ``a <= x <= c``.
  552. bounds : sequence, optional
  553. For method 'bounded', `bounds` is mandatory and must have two items
  554. corresponding to the optimization bounds.
  555. args : tuple, optional
  556. Extra arguments passed to the objective function.
  557. method : str or callable, optional
  558. Type of solver. Should be one of:
  559. - 'Brent' :ref:`(see here) <optimize.minimize_scalar-brent>`
  560. - 'Bounded' :ref:`(see here) <optimize.minimize_scalar-bounded>`
  561. - 'Golden' :ref:`(see here) <optimize.minimize_scalar-golden>`
  562. - custom - a callable object (added in version 0.14.0), see below
  563. tol : float, optional
  564. Tolerance for termination. For detailed control, use solver-specific
  565. options.
  566. options : dict, optional
  567. A dictionary of solver options.
  568. maxiter : int
  569. Maximum number of iterations to perform.
  570. disp : bool
  571. Set to True to print convergence messages.
  572. See :func:`show_options()` for solver-specific options.
  573. Returns
  574. -------
  575. res : OptimizeResult
  576. The optimization result represented as a ``OptimizeResult`` object.
  577. Important attributes are: ``x`` the solution array, ``success`` a
  578. Boolean flag indicating if the optimizer exited successfully and
  579. ``message`` which describes the cause of the termination. See
  580. `OptimizeResult` for a description of other attributes.
  581. See also
  582. --------
  583. minimize : Interface to minimization algorithms for scalar multivariate
  584. functions
  585. show_options : Additional options accepted by the solvers
  586. Notes
  587. -----
  588. This section describes the available solvers that can be selected by the
  589. 'method' parameter. The default method is *Brent*.
  590. Method :ref:`Brent <optimize.minimize_scalar-brent>` uses Brent's
  591. algorithm to find a local minimum. The algorithm uses inverse
  592. parabolic interpolation when possible to speed up convergence of
  593. the golden section method.
  594. Method :ref:`Golden <optimize.minimize_scalar-golden>` uses the
  595. golden section search technique. It uses analog of the bisection
  596. method to decrease the bracketed interval. It is usually
  597. preferable to use the *Brent* method.
  598. Method :ref:`Bounded <optimize.minimize_scalar-bounded>` can
  599. perform bounded minimization. It uses the Brent method to find a
  600. local minimum in the interval x1 < xopt < x2.
  601. **Custom minimizers**
  602. It may be useful to pass a custom minimization method, for example
  603. when using some library frontend to minimize_scalar. You can simply
  604. pass a callable as the ``method`` parameter.
  605. The callable is called as ``method(fun, args, **kwargs, **options)``
  606. where ``kwargs`` corresponds to any other parameters passed to `minimize`
  607. (such as `bracket`, `tol`, etc.), except the `options` dict, which has
  608. its contents also passed as `method` parameters pair by pair. The method
  609. shall return an ``OptimizeResult`` object.
  610. The provided `method` callable must be able to accept (and possibly ignore)
  611. arbitrary parameters; the set of parameters accepted by `minimize` may
  612. expand in future versions and then these parameters will be passed to
  613. the method. You can find an example in the scipy.optimize tutorial.
  614. .. versionadded:: 0.11.0
  615. Examples
  616. --------
  617. Consider the problem of minimizing the following function.
  618. >>> def f(x):
  619. ... return (x - 2) * x * (x + 2)**2
  620. Using the *Brent* method, we find the local minimum as:
  621. >>> from scipy.optimize import minimize_scalar
  622. >>> res = minimize_scalar(f)
  623. >>> res.x
  624. 1.28077640403
  625. Using the *Bounded* method, we find a local minimum with specified
  626. bounds as:
  627. >>> res = minimize_scalar(f, bounds=(-3, -1), method='bounded')
  628. >>> res.x
  629. -2.0000002026
  630. """
  631. if not isinstance(args, tuple):
  632. args = (args,)
  633. if callable(method):
  634. meth = "_custom"
  635. else:
  636. meth = method.lower()
  637. if options is None:
  638. options = {}
  639. if tol is not None:
  640. options = dict(options)
  641. if meth == 'bounded' and 'xatol' not in options:
  642. warn("Method 'bounded' does not support relative tolerance in x; "
  643. "defaulting to absolute tolerance.", RuntimeWarning)
  644. options['xatol'] = tol
  645. elif meth == '_custom':
  646. options.setdefault('tol', tol)
  647. else:
  648. options.setdefault('xtol', tol)
  649. if meth == '_custom':
  650. return method(fun, args=args, bracket=bracket, bounds=bounds, **options)
  651. elif meth == 'brent':
  652. return _minimize_scalar_brent(fun, bracket, args, **options)
  653. elif meth == 'bounded':
  654. if bounds is None:
  655. raise ValueError('The `bounds` parameter is mandatory for '
  656. 'method `bounded`.')
  657. # replace boolean "disp" option, if specified, by an integer value, as
  658. # expected by _minimize_scalar_bounded()
  659. disp = options.get('disp')
  660. if isinstance(disp, bool):
  661. options['disp'] = 2 * int(disp)
  662. return _minimize_scalar_bounded(fun, bounds, args, **options)
  663. elif meth == 'golden':
  664. return _minimize_scalar_golden(fun, bracket, args, **options)
  665. else:
  666. raise ValueError('Unknown solver %s' % method)
  667. def standardize_bounds(bounds, x0, meth):
  668. """Converts bounds to the form required by the solver."""
  669. if meth == 'trust-constr':
  670. if not isinstance(bounds, Bounds):
  671. lb, ub = old_bound_to_new(bounds)
  672. bounds = Bounds(lb, ub)
  673. elif meth in ('l-bfgs-b', 'tnc', 'slsqp'):
  674. if isinstance(bounds, Bounds):
  675. bounds = new_bounds_to_old(bounds.lb, bounds.ub, x0.shape[0])
  676. return bounds
  677. def standardize_constraints(constraints, x0, meth):
  678. """Converts constraints to the form required by the solver."""
  679. all_constraint_types = (NonlinearConstraint, LinearConstraint, dict)
  680. new_constraint_types = all_constraint_types[:-1]
  681. if isinstance(constraints, all_constraint_types):
  682. constraints = [constraints]
  683. constraints = list(constraints) # ensure it's a mutable sequence
  684. if meth == 'trust-constr':
  685. for i, con in enumerate(constraints):
  686. if not isinstance(con, new_constraint_types):
  687. constraints[i] = old_constraint_to_new(i, con)
  688. else:
  689. # iterate over copy, changing original
  690. for i, con in enumerate(list(constraints)):
  691. if isinstance(con, new_constraint_types):
  692. old_constraints = new_constraint_to_old(con, x0)
  693. constraints[i] = old_constraints[0]
  694. constraints.extend(old_constraints[1:]) # appends 1 if present
  695. return constraints