_minimize.py 52 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199
  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. __all__ = ['minimize', 'minimize_scalar']
  9. from warnings import warn
  10. import numpy as np
  11. from scipy._lib._util import wrapped_inspect_signature
  12. # unconstrained minimization
  13. from ._optimize import (_minimize_neldermead, _minimize_powell, _minimize_cg,
  14. _minimize_bfgs, _minimize_newtoncg,
  15. _minimize_scalar_brent, _minimize_scalar_bounded,
  16. _minimize_scalar_golden, MemoizeJac, OptimizeResult,
  17. _wrap_callback, _recover_from_bracket_error)
  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_py import _minimize_lbfgsb
  25. from ._tnc import _minimize_tnc
  26. from ._cobyla_py import _minimize_cobyla
  27. from ._cobyqa_py import _minimize_cobyqa
  28. from ._slsqp_py import _minimize_slsqp
  29. from ._constraints import (old_bound_to_new, new_bounds_to_old,
  30. old_constraint_to_new, new_constraint_to_old,
  31. NonlinearConstraint, LinearConstraint, Bounds,
  32. PreparedConstraint)
  33. from ._differentiable_functions import FD_METHODS
  34. MINIMIZE_METHODS = ['nelder-mead', 'powell', 'cg', 'bfgs', 'newton-cg',
  35. 'l-bfgs-b', 'tnc', 'cobyla', 'cobyqa', 'slsqp',
  36. 'trust-constr', 'dogleg', 'trust-ncg', 'trust-exact',
  37. 'trust-krylov']
  38. # These methods support the new callback interface (passed an OptimizeResult)
  39. MINIMIZE_METHODS_NEW_CB = ['nelder-mead', 'powell', 'cg', 'bfgs', 'newton-cg',
  40. 'l-bfgs-b', 'trust-constr', 'dogleg', 'trust-ncg',
  41. 'trust-exact', 'trust-krylov', 'cobyqa', 'cobyla', 'slsqp']
  42. MINIMIZE_SCALAR_METHODS = ['brent', 'bounded', 'golden']
  43. def minimize(fun, x0, args=(), method=None, jac=None, hess=None,
  44. hessp=None, bounds=None, constraints=(), tol=None,
  45. callback=None, options=None):
  46. """Minimization of scalar function of one or more variables.
  47. Parameters
  48. ----------
  49. fun : callable
  50. The objective function to be minimized::
  51. fun(x, *args) -> float
  52. where ``x`` is a 1-D array with shape (n,) and ``args``
  53. is a tuple of the fixed parameters needed to completely
  54. specify the function.
  55. Suppose the callable has signature ``f0(x, *my_args, **my_kwargs)``, where
  56. ``my_args`` and ``my_kwargs`` are required positional and keyword arguments.
  57. Rather than passing ``f0`` as the callable, wrap it to accept
  58. only ``x``; e.g., pass ``fun=lambda x: f0(x, *my_args, **my_kwargs)`` as the
  59. callable, where ``my_args`` (tuple) and ``my_kwargs`` (dict) have been
  60. gathered before invoking this function.
  61. x0 : ndarray, shape (n,)
  62. Initial guess. Array of real elements of size (n,),
  63. where ``n`` is the number of independent variables.
  64. args : tuple, optional
  65. Extra arguments passed to the objective function and its
  66. derivatives (`fun`, `jac` and `hess` functions).
  67. method : str or callable, optional
  68. Type of solver. Should be one of
  69. - 'Nelder-Mead' :ref:`(see here) <optimize.minimize-neldermead>`
  70. - 'Powell' :ref:`(see here) <optimize.minimize-powell>`
  71. - 'CG' :ref:`(see here) <optimize.minimize-cg>`
  72. - 'BFGS' :ref:`(see here) <optimize.minimize-bfgs>`
  73. - 'Newton-CG' :ref:`(see here) <optimize.minimize-newtoncg>`
  74. - 'L-BFGS-B' :ref:`(see here) <optimize.minimize-lbfgsb>`
  75. - 'TNC' :ref:`(see here) <optimize.minimize-tnc>`
  76. - 'COBYLA' :ref:`(see here) <optimize.minimize-cobyla>`
  77. - 'COBYQA' :ref:`(see here) <optimize.minimize-cobyqa>`
  78. - 'SLSQP' :ref:`(see here) <optimize.minimize-slsqp>`
  79. - 'trust-constr':ref:`(see here) <optimize.minimize-trustconstr>`
  80. - 'dogleg' :ref:`(see here) <optimize.minimize-dogleg>`
  81. - 'trust-ncg' :ref:`(see here) <optimize.minimize-trustncg>`
  82. - 'trust-exact' :ref:`(see here) <optimize.minimize-trustexact>`
  83. - 'trust-krylov' :ref:`(see here) <optimize.minimize-trustkrylov>`
  84. - custom - a callable object, see below for description.
  85. If not given, chosen to be one of ``BFGS``, ``L-BFGS-B``, ``SLSQP``,
  86. depending on whether or not the problem has constraints or bounds.
  87. jac : {callable, '2-point', '3-point', 'cs', bool}, optional
  88. Method for computing the gradient vector. Only for CG, BFGS,
  89. Newton-CG, L-BFGS-B, TNC, SLSQP, dogleg, trust-ncg, trust-krylov,
  90. trust-exact and trust-constr.
  91. If it is a callable, it should be a function that returns the gradient
  92. vector::
  93. jac(x, *args) -> array_like, shape (n,)
  94. where ``x`` is an array with shape (n,) and ``args`` is a tuple with
  95. the fixed parameters. If `jac` is a Boolean and is True, `fun` is
  96. assumed to return a tuple ``(f, g)`` containing the objective
  97. function and the gradient.
  98. Methods 'Newton-CG', 'trust-ncg', 'dogleg', 'trust-exact', and
  99. 'trust-krylov' require that either a callable be supplied, or that
  100. `fun` return the objective and gradient.
  101. If None or False, the gradient will be estimated using 2-point finite
  102. difference estimation with an absolute step size.
  103. Alternatively, the keywords {'2-point', '3-point', 'cs'} can be used
  104. to select a finite difference scheme for numerical estimation of the
  105. gradient with a relative step size. These finite difference schemes
  106. obey any specified `bounds`.
  107. hess : {callable, '2-point', '3-point', 'cs', HessianUpdateStrategy}, optional
  108. Method for computing the Hessian matrix. Only for Newton-CG, dogleg,
  109. trust-ncg, trust-krylov, trust-exact and trust-constr.
  110. If it is callable, it should return the Hessian matrix::
  111. hess(x, *args) -> {LinearOperator, spmatrix, array}, (n, n)
  112. where ``x`` is a (n,) ndarray and ``args`` is a tuple with the fixed
  113. parameters.
  114. The keywords {'2-point', '3-point', 'cs'} can also be used to select
  115. a finite difference scheme for numerical estimation of the hessian.
  116. Alternatively, objects implementing the `HessianUpdateStrategy`
  117. interface can be used to approximate the Hessian. Available
  118. quasi-Newton methods implementing this interface are:
  119. - `BFGS`
  120. - `SR1`
  121. Not all of the options are available for each of the methods; for
  122. availability refer to the notes.
  123. hessp : callable, optional
  124. Hessian of objective function times an arbitrary vector p. Only for
  125. Newton-CG, trust-ncg, trust-krylov, trust-constr.
  126. Only one of `hessp` or `hess` needs to be given. If `hess` is
  127. provided, then `hessp` will be ignored. `hessp` must compute the
  128. Hessian times an arbitrary vector::
  129. hessp(x, p, *args) -> ndarray shape (n,)
  130. where ``x`` is a (n,) ndarray, ``p`` is an arbitrary vector with
  131. dimension (n,) and ``args`` is a tuple with the fixed
  132. parameters.
  133. bounds : sequence or `Bounds`, optional
  134. Bounds on variables for Nelder-Mead, L-BFGS-B, TNC, SLSQP, Powell,
  135. trust-constr, COBYLA, and COBYQA methods. There are two ways to specify
  136. the bounds:
  137. 1. Instance of `Bounds` class.
  138. 2. Sequence of ``(min, max)`` pairs for each element in `x`. None
  139. is used to specify no bound.
  140. constraints : {Constraint, dict} or List of {Constraint, dict}, optional
  141. Constraints definition. Only for COBYLA, COBYQA, SLSQP and trust-constr.
  142. Constraints for 'trust-constr', 'cobyqa', and 'cobyla' are defined as a single
  143. object or a list of objects specifying constraints to the optimization problem.
  144. Available constraints are:
  145. - `LinearConstraint`
  146. - `NonlinearConstraint`
  147. Constraints for COBYLA, SLSQP are defined as a list of dictionaries.
  148. Each dictionary with fields:
  149. type : str
  150. Constraint type: 'eq' for equality, 'ineq' for inequality.
  151. fun : callable
  152. The function defining the constraint.
  153. jac : callable, optional
  154. The Jacobian of `fun` (only for SLSQP).
  155. args : sequence, optional
  156. Extra arguments to be passed to the function and Jacobian.
  157. Equality constraint means that the constraint function result is to
  158. be zero whereas inequality means that it is to be non-negative.
  159. tol : float, optional
  160. Tolerance for termination. When `tol` is specified, the selected
  161. minimization algorithm sets some relevant solver-specific tolerance(s)
  162. equal to `tol`. For detailed control, use solver-specific
  163. options.
  164. options : dict, optional
  165. A dictionary of solver options. All methods except `TNC` accept the
  166. following generic options:
  167. maxiter : int
  168. Maximum number of iterations to perform. Depending on the
  169. method each iteration may use several function evaluations.
  170. For `TNC` use `maxfun` instead of `maxiter`.
  171. disp : bool
  172. Set to True to print convergence messages.
  173. For method-specific options, see :func:`show_options()`.
  174. callback : callable, optional
  175. A callable called after each iteration.
  176. All methods except TNC support a callable with
  177. the signature::
  178. callback(intermediate_result: OptimizeResult)
  179. where ``intermediate_result`` is a keyword parameter containing an
  180. `OptimizeResult` with attributes ``x`` and ``fun``, the present values
  181. of the parameter vector and objective function. Not all attributes of
  182. `OptimizeResult` may be present. The name of the parameter must be
  183. ``intermediate_result`` for the callback to be passed an `OptimizeResult`.
  184. These methods will also terminate if the callback raises ``StopIteration``.
  185. All methods except trust-constr (also) support a signature like::
  186. callback(xk)
  187. where ``xk`` is the current parameter vector.
  188. Introspection is used to determine which of the signatures above to
  189. invoke.
  190. Returns
  191. -------
  192. res : OptimizeResult
  193. The optimization result represented as a ``OptimizeResult`` object.
  194. Important attributes are: ``x`` the solution array, ``success`` a
  195. Boolean flag indicating if the optimizer exited successfully and
  196. ``message`` which describes the cause of the termination. See
  197. `OptimizeResult` for a description of other attributes.
  198. See also
  199. --------
  200. minimize_scalar : Interface to minimization algorithms for scalar
  201. univariate functions
  202. show_options : Additional options accepted by the solvers
  203. Notes
  204. -----
  205. This section describes the available solvers that can be selected by the
  206. 'method' parameter. The default method is *BFGS*.
  207. **Unconstrained minimization**
  208. Method :ref:`CG <optimize.minimize-cg>` uses a nonlinear conjugate
  209. gradient algorithm by Polak and Ribiere, a variant of the
  210. Fletcher-Reeves method described in [5]_ pp.120-122. Only the
  211. first derivatives are used.
  212. Method :ref:`BFGS <optimize.minimize-bfgs>` uses the quasi-Newton
  213. method of Broyden, Fletcher, Goldfarb, and Shanno (BFGS) [5]_
  214. pp. 136. It uses the first derivatives only. BFGS has proven good
  215. performance even for non-smooth optimizations. This method also
  216. returns an approximation of the Hessian inverse, stored as
  217. `hess_inv` in the OptimizeResult object.
  218. Method :ref:`Newton-CG <optimize.minimize-newtoncg>` uses a
  219. Newton-CG algorithm [5]_ pp. 168 (also known as the truncated
  220. Newton method). It uses a CG method to the compute the search
  221. direction. See also *TNC* method for a box-constrained
  222. minimization with a similar algorithm. Suitable for large-scale
  223. problems.
  224. Method :ref:`dogleg <optimize.minimize-dogleg>` uses the dog-leg
  225. trust-region algorithm [5]_ for unconstrained minimization. This
  226. algorithm requires the gradient and Hessian; furthermore the
  227. Hessian is required to be positive definite.
  228. Method :ref:`trust-ncg <optimize.minimize-trustncg>` uses the
  229. Newton conjugate gradient trust-region algorithm [5]_ for
  230. unconstrained minimization. This algorithm requires the gradient
  231. and either the Hessian or a function that computes the product of
  232. the Hessian with a given vector. Suitable for large-scale problems.
  233. Method :ref:`trust-krylov <optimize.minimize-trustkrylov>` uses
  234. the Newton GLTR trust-region algorithm [14]_, [15]_ for unconstrained
  235. minimization. This algorithm requires the gradient
  236. and either the Hessian or a function that computes the product of
  237. the Hessian with a given vector. Suitable for large-scale problems.
  238. On indefinite problems it requires usually less iterations than the
  239. `trust-ncg` method and is recommended for medium and large-scale problems.
  240. Method :ref:`trust-exact <optimize.minimize-trustexact>`
  241. is a trust-region method for unconstrained minimization in which
  242. quadratic subproblems are solved almost exactly [13]_. This
  243. algorithm requires the gradient and the Hessian (which is
  244. *not* required to be positive definite). It is, in many
  245. situations, the Newton method to converge in fewer iterations
  246. and the most recommended for small and medium-size problems.
  247. **Bound-Constrained minimization**
  248. Method :ref:`Nelder-Mead <optimize.minimize-neldermead>` uses the
  249. Simplex algorithm [1]_, [2]_. This algorithm is robust in many
  250. applications. However, if numerical computation of derivative can be
  251. trusted, other algorithms using the first and/or second derivatives
  252. information might be preferred for their better performance in
  253. general.
  254. Method :ref:`L-BFGS-B <optimize.minimize-lbfgsb>` uses the L-BFGS-B
  255. algorithm [6]_, [7]_ for bound constrained minimization.
  256. Method :ref:`Powell <optimize.minimize-powell>` is a modification
  257. of Powell's method [3]_, [4]_ which is a conjugate direction
  258. method. It performs sequential one-dimensional minimizations along
  259. each vector of the directions set (`direc` field in `options` and
  260. `info`), which is updated at each iteration of the main
  261. minimization loop. The function need not be differentiable, and no
  262. derivatives are taken. If bounds are not provided, then an
  263. unbounded line search will be used. If bounds are provided and
  264. the initial guess is within the bounds, then every function
  265. evaluation throughout the minimization procedure will be within
  266. the bounds. If bounds are provided, the initial guess is outside
  267. the bounds, and `direc` is full rank (default has full rank), then
  268. some function evaluations during the first iteration may be
  269. outside the bounds, but every function evaluation after the first
  270. iteration will be within the bounds. If `direc` is not full rank,
  271. then some parameters may not be optimized and the solution is not
  272. guaranteed to be within the bounds.
  273. Method :ref:`TNC <optimize.minimize-tnc>` uses a truncated Newton
  274. algorithm [5]_, [8]_ to minimize a function with variables subject
  275. to bounds. This algorithm uses gradient information; it is also
  276. called Newton Conjugate-Gradient. It differs from the *Newton-CG*
  277. method described above as it wraps a C implementation and allows
  278. each variable to be given upper and lower bounds.
  279. **Constrained Minimization**
  280. Method :ref:`COBYLA <optimize.minimize-cobyla>` uses the PRIMA
  281. implementation [19]_ of the
  282. Constrained Optimization BY Linear Approximation (COBYLA) method
  283. [9]_, [10]_, [11]_. The algorithm is based on linear
  284. approximations to the objective function and each constraint.
  285. Method :ref:`COBYQA <optimize.minimize-cobyqa>` uses the Constrained
  286. Optimization BY Quadratic Approximations (COBYQA) method [18]_. The
  287. algorithm is a derivative-free trust-region SQP method based on quadratic
  288. approximations to the objective function and each nonlinear constraint. The
  289. bounds are treated as unrelaxable constraints, in the sense that the
  290. algorithm always respects them throughout the optimization process.
  291. Method :ref:`SLSQP <optimize.minimize-slsqp>` uses Sequential
  292. Least SQuares Programming to minimize a function of several
  293. variables with any combination of bounds, equality and inequality
  294. constraints. The method wraps the SLSQP Optimization subroutine
  295. originally implemented by Dieter Kraft [12]_. Note that the
  296. wrapper handles infinite values in bounds by converting them into
  297. large floating values.
  298. Method :ref:`trust-constr <optimize.minimize-trustconstr>` is a
  299. trust-region algorithm for constrained optimization. It switches
  300. between two implementations depending on the problem definition.
  301. It is the most versatile constrained minimization algorithm
  302. implemented in SciPy and the most appropriate for large-scale problems.
  303. For equality constrained problems it is an implementation of Byrd-Omojokun
  304. Trust-Region SQP method described in [17]_ and in [5]_, p. 549. When
  305. inequality constraints are imposed as well, it switches to the trust-region
  306. interior point method described in [16]_. This interior point algorithm,
  307. in turn, solves inequality constraints by introducing slack variables
  308. and solving a sequence of equality-constrained barrier problems
  309. for progressively smaller values of the barrier parameter.
  310. The previously described equality constrained SQP method is
  311. used to solve the subproblems with increasing levels of accuracy
  312. as the iterate gets closer to a solution.
  313. **Finite-Difference Options**
  314. For Method :ref:`trust-constr <optimize.minimize-trustconstr>`
  315. the gradient and the Hessian may be approximated using
  316. three finite-difference schemes: {'2-point', '3-point', 'cs'}.
  317. The scheme 'cs' is, potentially, the most accurate but it
  318. requires the function to correctly handle complex inputs and to
  319. be differentiable in the complex plane. The scheme '3-point' is more
  320. accurate than '2-point' but requires twice as many operations. If the
  321. gradient is estimated via finite-differences the Hessian must be
  322. estimated using one of the quasi-Newton strategies.
  323. **Method specific options for the** `hess` **keyword**
  324. +--------------+------+----------+-------------------------+-----+
  325. | method/Hess | None | callable | '2-point/'3-point'/'cs' | HUS |
  326. +==============+======+==========+=========================+=====+
  327. | Newton-CG | x | (n, n) | x | x |
  328. | | | LO | | |
  329. +--------------+------+----------+-------------------------+-----+
  330. | dogleg | | (n, n) | | |
  331. +--------------+------+----------+-------------------------+-----+
  332. | trust-ncg | | (n, n) | x | x |
  333. +--------------+------+----------+-------------------------+-----+
  334. | trust-krylov | | (n, n) | x | x |
  335. +--------------+------+----------+-------------------------+-----+
  336. | trust-exact | | (n, n) | | |
  337. +--------------+------+----------+-------------------------+-----+
  338. | trust-constr | x | (n, n) | x | x |
  339. | | | LO | | |
  340. | | | sp | | |
  341. +--------------+------+----------+-------------------------+-----+
  342. where LO=LinearOperator, sp=Sparse matrix, HUS=HessianUpdateStrategy
  343. **Custom minimizers**
  344. It may be useful to pass a custom minimization method, for example
  345. when using a frontend to this method such as `scipy.optimize.basinhopping`
  346. or a different library. You can simply pass a callable as the ``method``
  347. parameter.
  348. The callable is called as ``method(fun, x0, args, **kwargs, **options)``
  349. where ``kwargs`` corresponds to any other parameters passed to `minimize`
  350. (such as `callback`, `hess`, etc.), except the `options` dict, which has
  351. its contents also passed as `method` parameters pair by pair. Also, if
  352. `jac` has been passed as a bool type, `jac` and `fun` are mangled so that
  353. `fun` returns just the function values and `jac` is converted to a function
  354. returning the Jacobian. The method shall return an `OptimizeResult`
  355. object.
  356. The provided `method` callable must be able to accept (and possibly ignore)
  357. arbitrary parameters; the set of parameters accepted by `minimize` may
  358. expand in future versions and then these parameters will be passed to
  359. the method. You can find an example in the scipy.optimize tutorial.
  360. References
  361. ----------
  362. .. [1] Nelder, J A, and R Mead. 1965. A Simplex Method for Function
  363. Minimization. The Computer Journal 7: 308-13.
  364. .. [2] Wright M H. 1996. Direct search methods: Once scorned, now
  365. respectable, in Numerical Analysis 1995: Proceedings of the 1995
  366. Dundee Biennial Conference in Numerical Analysis (Eds. D F
  367. Griffiths and G A Watson). Addison Wesley Longman, Harlow, UK.
  368. 191-208.
  369. .. [3] Powell, M J D. 1964. An efficient method for finding the minimum of
  370. a function of several variables without calculating derivatives. The
  371. Computer Journal 7: 155-162.
  372. .. [4] Press W, S A Teukolsky, W T Vetterling and B P Flannery.
  373. Numerical Recipes (any edition), Cambridge University Press.
  374. .. [5] Nocedal, J, and S J Wright. 2006. Numerical Optimization.
  375. Springer New York.
  376. .. [6] Byrd, R H and P Lu and J. Nocedal. 1995. A Limited Memory
  377. Algorithm for Bound Constrained Optimization. SIAM Journal on
  378. Scientific and Statistical Computing 16 (5): 1190-1208.
  379. .. [7] Zhu, C and R H Byrd and J Nocedal. 1997. L-BFGS-B: Algorithm
  380. 778: L-BFGS-B, FORTRAN routines for large scale bound constrained
  381. optimization. ACM Transactions on Mathematical Software 23 (4):
  382. 550-560.
  383. .. [8] Nash, S G. Newton-Type Minimization Via the Lanczos Method.
  384. 1984. SIAM Journal of Numerical Analysis 21: 770-778.
  385. .. [9] Powell, M J D. A direct search optimization method that models
  386. the objective and constraint functions by linear interpolation.
  387. 1994. Advances in Optimization and Numerical Analysis, eds. S. Gomez
  388. and J-P Hennart, Kluwer Academic (Dordrecht), 51-67.
  389. .. [10] Powell M J D. Direct search algorithms for optimization
  390. calculations. 1998. Acta Numerica 7: 287-336.
  391. .. [11] Powell M J D. A view of algorithms for optimization without
  392. derivatives. 2007.Cambridge University Technical Report DAMTP
  393. 2007/NA03
  394. .. [12] Kraft, D. A software package for sequential quadratic
  395. programming. 1988. Tech. Rep. DFVLR-FB 88-28, DLR German Aerospace
  396. Center -- Institute for Flight Mechanics, Koln, Germany.
  397. .. [13] Conn, A. R., Gould, N. I., and Toint, P. L.
  398. Trust region methods. 2000. Siam. pp. 169-200.
  399. .. [14] F. Lenders, C. Kirches, A. Potschka: "trlib: A vector-free
  400. implementation of the GLTR method for iterative solution of
  401. the trust region problem", :arxiv:`1611.04718`
  402. .. [15] N. Gould, S. Lucidi, M. Roma, P. Toint: "Solving the
  403. Trust-Region Subproblem using the Lanczos Method",
  404. SIAM J. Optim., 9(2), 504--525, (1999).
  405. .. [16] Byrd, Richard H., Mary E. Hribar, and Jorge Nocedal. 1999.
  406. An interior point algorithm for large-scale nonlinear programming.
  407. SIAM Journal on Optimization 9.4: 877-900.
  408. .. [17] Lalee, Marucha, Jorge Nocedal, and Todd Plantenga. 1998. On the
  409. implementation of an algorithm for large-scale equality constrained
  410. optimization. SIAM Journal on Optimization 8.3: 682-706.
  411. .. [18] Ragonneau, T. M. *Model-Based Derivative-Free Optimization Methods
  412. and Software*. PhD thesis, Department of Applied Mathematics, The Hong
  413. Kong Polytechnic University, Hong Kong, China, 2022. URL:
  414. https://theses.lib.polyu.edu.hk/handle/200/12294.
  415. .. [19] Zhang, Z. "PRIMA: Reference Implementation for Powell's Methods with
  416. Modernization and Amelioration", https://www.libprima.net,
  417. :doi:`10.5281/zenodo.8052654`
  418. .. [20] Karush-Kuhn-Tucker conditions,
  419. https://en.wikipedia.org/wiki/Karush%E2%80%93Kuhn%E2%80%93Tucker_conditions
  420. Examples
  421. --------
  422. Let us consider the problem of minimizing the Rosenbrock function. This
  423. function (and its respective derivatives) is implemented in `rosen`
  424. (resp. `rosen_der`, `rosen_hess`) in the `scipy.optimize`.
  425. >>> from scipy.optimize import minimize, rosen, rosen_der
  426. A simple application of the *Nelder-Mead* method is:
  427. >>> x0 = [1.3, 0.7, 0.8, 1.9, 1.2]
  428. >>> res = minimize(rosen, x0, method='Nelder-Mead', tol=1e-6)
  429. >>> res.x
  430. array([ 1., 1., 1., 1., 1.])
  431. Now using the *BFGS* algorithm, using the first derivative and a few
  432. options:
  433. >>> res = minimize(rosen, x0, method='BFGS', jac=rosen_der,
  434. ... options={'gtol': 1e-6, 'disp': True})
  435. Optimization terminated successfully.
  436. Current function value: 0.000000
  437. Iterations: 26
  438. Function evaluations: 31
  439. Gradient evaluations: 31
  440. >>> res.x
  441. array([ 1., 1., 1., 1., 1.])
  442. >>> print(res.message)
  443. Optimization terminated successfully.
  444. >>> res.hess_inv
  445. array([
  446. [ 0.00749589, 0.01255155, 0.02396251, 0.04750988, 0.09495377], # may vary
  447. [ 0.01255155, 0.02510441, 0.04794055, 0.09502834, 0.18996269],
  448. [ 0.02396251, 0.04794055, 0.09631614, 0.19092151, 0.38165151],
  449. [ 0.04750988, 0.09502834, 0.19092151, 0.38341252, 0.7664427 ],
  450. [ 0.09495377, 0.18996269, 0.38165151, 0.7664427, 1.53713523]
  451. ])
  452. Next, consider a minimization problem with several constraints (namely
  453. Example 16.4 from [5]_). The objective function is:
  454. >>> fun = lambda x: (x[0] - 1)**2 + (x[1] - 2.5)**2
  455. There are three constraints defined as:
  456. >>> cons = ({'type': 'ineq', 'fun': lambda x: x[0] - 2 * x[1] + 2},
  457. ... {'type': 'ineq', 'fun': lambda x: -x[0] - 2 * x[1] + 6},
  458. ... {'type': 'ineq', 'fun': lambda x: -x[0] + 2 * x[1] + 2})
  459. And variables must be positive, hence the following bounds:
  460. >>> bnds = ((0, None), (0, None))
  461. The optimization problem is solved using the SLSQP method as:
  462. >>> res = minimize(fun, (2, 0), method='SLSQP', bounds=bnds, constraints=cons)
  463. It should converge to the theoretical solution ``[1.4 ,1.7]``. *SLSQP* also
  464. returns the multipliers that are used in the solution of the problem. These
  465. multipliers, when the problem constraints are linear, can be thought of as the
  466. Karush-Kuhn-Tucker (KKT) multipliers, which are a generalization
  467. of Lagrange multipliers to inequality-constrained optimization problems ([20]_).
  468. Notice that at the solution, the first constraint is active. Let's evaluate the
  469. function at solution:
  470. >>> cons[0]['fun'](res.x)
  471. np.float64(1.4901224698604665e-09)
  472. Also, notice that at optimality there is a non-zero multiplier:
  473. >>> res.multipliers
  474. array([0.8, 0. , 0. ])
  475. This can be understood as the local sensitivity of the optimal value of the
  476. objective function with respect to changes in the first constraint. If we
  477. tighten the constraint by a small amount ``eps``:
  478. >>> eps = 0.01
  479. >>> cons[0]['fun'] = lambda x: x[0] - 2 * x[1] + 2 - eps
  480. we expect the optimal value of the objective function to increase by
  481. approximately ``eps * res.multipliers[0]``:
  482. >>> eps * res.multipliers[0] # Expected change in f0
  483. np.float64(0.008000000027153205)
  484. >>> f0 = res.fun # Keep track of the previous optimal value
  485. >>> res = minimize(fun, (2, 0), method='SLSQP', bounds=bnds, constraints=cons)
  486. >>> f1 = res.fun # New optimal value
  487. >>> f1 - f0
  488. np.float64(0.008019998807885509)
  489. """
  490. x0 = np.atleast_1d(np.asarray(x0))
  491. if x0.ndim != 1:
  492. raise ValueError("'x0' must only have one dimension.")
  493. if x0.dtype.kind in np.typecodes["AllInteger"]:
  494. x0 = np.asarray(x0, dtype=float)
  495. if not isinstance(args, tuple):
  496. args = (args,)
  497. if method is None:
  498. # Select automatically
  499. if constraints:
  500. method = 'SLSQP'
  501. elif bounds is not None:
  502. method = 'L-BFGS-B'
  503. else:
  504. method = 'BFGS'
  505. if callable(method):
  506. meth = "_custom"
  507. else:
  508. meth = method.lower()
  509. if options is None:
  510. options = {}
  511. # check if optional parameters are supported by the selected method
  512. # - jac
  513. if meth in ('nelder-mead', 'powell', 'cobyla', 'cobyqa') and bool(jac):
  514. warn(f'Method {method} does not use gradient information (jac).',
  515. RuntimeWarning, stacklevel=2)
  516. # - hess
  517. if meth not in ('newton-cg', 'dogleg', 'trust-ncg', 'trust-constr',
  518. 'trust-krylov', 'trust-exact', '_custom') and hess is not None:
  519. warn(f'Method {method} does not use Hessian information (hess).',
  520. RuntimeWarning, stacklevel=2)
  521. # - hessp
  522. if meth not in ('newton-cg', 'trust-ncg', 'trust-constr',
  523. 'trust-krylov', '_custom') \
  524. and hessp is not None:
  525. warn(f'Method {method} does not use Hessian-vector product'
  526. ' information (hessp).',
  527. RuntimeWarning, stacklevel=2)
  528. # - constraints or bounds
  529. if (meth not in ('cobyla', 'cobyqa', 'slsqp', 'trust-constr', '_custom') and
  530. np.any(constraints)):
  531. warn(f'Method {method} cannot handle constraints.',
  532. RuntimeWarning, stacklevel=2)
  533. if meth not in (
  534. 'nelder-mead', 'powell', 'l-bfgs-b', 'cobyla', 'cobyqa', 'slsqp',
  535. 'tnc', 'trust-constr', '_custom') and bounds is not None:
  536. warn(f'Method {method} cannot handle bounds.',
  537. RuntimeWarning, stacklevel=2)
  538. # - return_all
  539. if (meth in ('l-bfgs-b', 'tnc', 'cobyla', 'cobyqa', 'slsqp') and
  540. options.get('return_all', False)):
  541. warn(f'Method {method} does not support the return_all option.',
  542. RuntimeWarning, stacklevel=2)
  543. # check gradient vector
  544. if callable(jac):
  545. pass
  546. elif jac is True:
  547. # fun returns func and grad
  548. fun = MemoizeJac(fun)
  549. jac = fun.derivative
  550. elif (jac in FD_METHODS and
  551. meth in ['trust-constr', 'bfgs', 'cg', 'l-bfgs-b', 'tnc', 'slsqp']):
  552. # finite differences with relative step
  553. pass
  554. elif meth in ['trust-constr']:
  555. # default jac calculation for this method
  556. jac = '2-point'
  557. elif jac is None or bool(jac) is False:
  558. # this will cause e.g. LBFGS to use forward difference, absolute step
  559. jac = None
  560. else:
  561. # default if jac option is not understood
  562. jac = None
  563. # set default tolerances
  564. if tol is not None:
  565. options = dict(options)
  566. if meth == 'nelder-mead':
  567. options.setdefault('xatol', tol)
  568. options.setdefault('fatol', tol)
  569. if meth in ('newton-cg', 'powell', 'tnc'):
  570. options.setdefault('xtol', tol)
  571. if meth in ('powell', 'l-bfgs-b', 'tnc', 'slsqp'):
  572. options.setdefault('ftol', tol)
  573. if meth in ('bfgs', 'cg', 'l-bfgs-b', 'tnc', 'dogleg',
  574. 'trust-ncg', 'trust-exact', 'trust-krylov'):
  575. options.setdefault('gtol', tol)
  576. if meth in ('cobyla', '_custom'):
  577. options.setdefault('tol', tol)
  578. if meth == 'cobyqa':
  579. options.setdefault('final_tr_radius', tol)
  580. if meth == 'trust-constr':
  581. options.setdefault('xtol', tol)
  582. options.setdefault('gtol', tol)
  583. options.setdefault('barrier_tol', tol)
  584. if meth == '_custom':
  585. # custom method called before bounds and constraints are 'standardised'
  586. # custom method should be able to accept whatever bounds/constraints
  587. # are provided to it.
  588. return method(fun, x0, args=args, jac=jac, hess=hess, hessp=hessp,
  589. bounds=bounds, constraints=constraints,
  590. callback=callback, **options)
  591. constraints = standardize_constraints(constraints, x0, meth)
  592. remove_vars = False
  593. if bounds is not None:
  594. # convert to new-style bounds so we only have to consider one case
  595. bounds = standardize_bounds(bounds, x0, 'new')
  596. bounds = _validate_bounds(bounds, x0, meth)
  597. if meth in {"tnc", "slsqp", "l-bfgs-b"}:
  598. # These methods can't take the finite-difference derivatives they
  599. # need when a variable is fixed by the bounds. To avoid this issue,
  600. # remove fixed variables from the problem.
  601. # NOTE: if this list is expanded, then be sure to update the
  602. # accompanying tests and test_optimize.eb_data. Consider also if
  603. # default OptimizeResult will need updating.
  604. # determine whether any variables are fixed
  605. i_fixed = (bounds.lb == bounds.ub)
  606. if np.all(i_fixed):
  607. # all the parameters are fixed, a minimizer is not able to do
  608. # anything
  609. return _optimize_result_for_equal_bounds(
  610. fun, bounds, meth, args=args, constraints=constraints
  611. )
  612. # determine whether finite differences are needed for any grad/jac
  613. fd_needed = (not callable(jac))
  614. for con in constraints:
  615. if not callable(con.get('jac', None)):
  616. fd_needed = True
  617. # If finite differences are ever used, remove all fixed variables
  618. # Always remove fixed variables for TNC; see gh-14565
  619. remove_vars = i_fixed.any() and (fd_needed or meth == "tnc")
  620. if remove_vars:
  621. x_fixed = (bounds.lb)[i_fixed]
  622. x0 = x0[~i_fixed]
  623. bounds = _remove_from_bounds(bounds, i_fixed)
  624. fun = _Remove_From_Func(fun, i_fixed, x_fixed)
  625. if callable(callback):
  626. sig = wrapped_inspect_signature(callback)
  627. if set(sig.parameters) == {'intermediate_result'}:
  628. # callback(intermediate_result)
  629. print(callback)
  630. callback = _Patch_Callback_Equal_Variables(
  631. callback, i_fixed, x_fixed
  632. )
  633. else:
  634. # callback(x)
  635. callback = _Remove_From_Func(callback, i_fixed, x_fixed)
  636. if callable(jac):
  637. jac = _Remove_From_Func(jac, i_fixed, x_fixed, remove=1)
  638. # make a copy of the constraints so the user's version doesn't
  639. # get changed. (Shallow copy is ok)
  640. constraints = [con.copy() for con in constraints]
  641. for con in constraints: # yes, guaranteed to be a list
  642. con['fun'] = _Remove_From_Func(con['fun'], i_fixed,
  643. x_fixed, min_dim=1,
  644. remove=0)
  645. if callable(con.get('jac', None)):
  646. con['jac'] = _Remove_From_Func(con['jac'], i_fixed,
  647. x_fixed, min_dim=2,
  648. remove=1)
  649. bounds = standardize_bounds(bounds, x0, meth)
  650. # selects whether to use callback(x) or callback(intermediate_result)
  651. callback = _wrap_callback(callback, meth)
  652. if meth == 'nelder-mead':
  653. res = _minimize_neldermead(fun, x0, args, callback, bounds=bounds,
  654. **options)
  655. elif meth == 'powell':
  656. res = _minimize_powell(fun, x0, args, callback, bounds, **options)
  657. elif meth == 'cg':
  658. res = _minimize_cg(fun, x0, args, jac, callback, **options)
  659. elif meth == 'bfgs':
  660. res = _minimize_bfgs(fun, x0, args, jac, callback, **options)
  661. elif meth == 'newton-cg':
  662. res = _minimize_newtoncg(fun, x0, args, jac, hess, hessp, callback,
  663. **options)
  664. elif meth == 'l-bfgs-b':
  665. res = _minimize_lbfgsb(fun, x0, args, jac, bounds,
  666. callback=callback, **options)
  667. elif meth == 'tnc':
  668. res = _minimize_tnc(fun, x0, args, jac, bounds, callback=callback,
  669. **options)
  670. elif meth == 'cobyla':
  671. res = _minimize_cobyla(fun, x0, args, constraints, callback=callback,
  672. bounds=bounds, **options)
  673. elif meth == 'cobyqa':
  674. res = _minimize_cobyqa(fun, x0, args, bounds, constraints, callback,
  675. **options)
  676. elif meth == 'slsqp':
  677. res = _minimize_slsqp(fun, x0, args, jac, bounds,
  678. constraints, callback=callback, **options)
  679. elif meth == 'trust-constr':
  680. res = _minimize_trustregion_constr(fun, x0, args, jac, hess, hessp,
  681. bounds, constraints,
  682. callback=callback, **options)
  683. elif meth == 'dogleg':
  684. res = _minimize_dogleg(fun, x0, args, jac, hess,
  685. callback=callback, **options)
  686. elif meth == 'trust-ncg':
  687. res = _minimize_trust_ncg(fun, x0, args, jac, hess, hessp,
  688. callback=callback, **options)
  689. elif meth == 'trust-krylov':
  690. res = _minimize_trust_krylov(fun, x0, args, jac, hess, hessp,
  691. callback=callback, **options)
  692. elif meth == 'trust-exact':
  693. res = _minimize_trustregion_exact(fun, x0, args, jac, hess,
  694. callback=callback, **options)
  695. else:
  696. raise ValueError(f'Unknown solver {method}')
  697. if remove_vars:
  698. res.x = _add_to_array(res.x, i_fixed, x_fixed)
  699. res.jac = _add_to_array(res.jac, i_fixed, np.nan)
  700. if "hess_inv" in res:
  701. res.hess_inv = None # unknown
  702. if getattr(callback, 'stop_iteration', False):
  703. res.success = False
  704. res.status = 99
  705. res.message = "`callback` raised `StopIteration`."
  706. return res
  707. def minimize_scalar(fun, bracket=None, bounds=None, args=(),
  708. method=None, tol=None, options=None):
  709. """Local minimization of scalar function of one variable.
  710. Parameters
  711. ----------
  712. fun : callable
  713. Objective function.
  714. Scalar function, must return a scalar.
  715. Suppose the callable has signature ``f0(x, *my_args, **my_kwargs)``, where
  716. ``my_args`` and ``my_kwargs`` are required positional and keyword arguments.
  717. Rather than passing ``f0`` as the callable, wrap it to accept
  718. only ``x``; e.g., pass ``fun=lambda x: f0(x, *my_args, **my_kwargs)`` as the
  719. callable, where ``my_args`` (tuple) and ``my_kwargs`` (dict) have been
  720. gathered before invoking this function.
  721. bracket : sequence, optional
  722. For methods 'brent' and 'golden', `bracket` defines the bracketing
  723. interval and is required.
  724. Either a triple ``(xa, xb, xc)`` satisfying ``xa < xb < xc`` and
  725. ``func(xb) < func(xa) and func(xb) < func(xc)``, or a pair
  726. ``(xa, xb)`` to be used as initial points for a downhill bracket search
  727. (see `scipy.optimize.bracket`).
  728. The minimizer ``res.x`` will not necessarily satisfy
  729. ``xa <= res.x <= xb``.
  730. bounds : sequence, optional
  731. For method 'bounded', `bounds` is mandatory and must have two finite
  732. items corresponding to the optimization bounds.
  733. args : tuple, optional
  734. Extra arguments passed to the objective function.
  735. method : str or callable, optional
  736. Type of solver. Should be one of:
  737. - :ref:`Brent <optimize.minimize_scalar-brent>`
  738. - :ref:`Bounded <optimize.minimize_scalar-bounded>`
  739. - :ref:`Golden <optimize.minimize_scalar-golden>`
  740. - custom - a callable object (added in version 0.14.0), see below
  741. Default is "Bounded" if bounds are provided and "Brent" otherwise.
  742. See the 'Notes' section for details of each solver.
  743. tol : float, optional
  744. Tolerance for termination. For detailed control, use solver-specific
  745. options.
  746. options : dict, optional
  747. A dictionary of solver options.
  748. maxiter : int
  749. Maximum number of iterations to perform.
  750. disp : bool
  751. Set to True to print convergence messages.
  752. See :func:`show_options()` for solver-specific options.
  753. Returns
  754. -------
  755. res : OptimizeResult
  756. The optimization result represented as a ``OptimizeResult`` object.
  757. Important attributes are: ``x`` the solution array, ``success`` a
  758. Boolean flag indicating if the optimizer exited successfully and
  759. ``message`` which describes the cause of the termination. See
  760. `OptimizeResult` for a description of other attributes.
  761. See also
  762. --------
  763. minimize : Interface to minimization algorithms for scalar multivariate
  764. functions
  765. show_options : Additional options accepted by the solvers
  766. Notes
  767. -----
  768. This section describes the available solvers that can be selected by the
  769. 'method' parameter. The default method is the ``"Bounded"`` Brent method if
  770. `bounds` are passed and unbounded ``"Brent"`` otherwise.
  771. Method :ref:`Brent <optimize.minimize_scalar-brent>` uses Brent's
  772. algorithm [1]_ to find a local minimum. The algorithm uses inverse
  773. parabolic interpolation when possible to speed up convergence of
  774. the golden section method.
  775. Method :ref:`Golden <optimize.minimize_scalar-golden>` uses the
  776. golden section search technique [1]_. It uses analog of the bisection
  777. method to decrease the bracketed interval. It is usually
  778. preferable to use the *Brent* method.
  779. Method :ref:`Bounded <optimize.minimize_scalar-bounded>` can
  780. perform bounded minimization [2]_ [3]_. It uses the Brent method to find a
  781. local minimum in the interval x1 < xopt < x2.
  782. Note that the Brent and Golden methods do not guarantee success unless a
  783. valid ``bracket`` triple is provided. If a three-point bracket cannot be
  784. found, consider `scipy.optimize.minimize`. Also, all methods are intended
  785. only for local minimization. When the function of interest has more than
  786. one local minimum, consider :ref:`global_optimization`.
  787. **Custom minimizers**
  788. It may be useful to pass a custom minimization method, for example
  789. when using some library frontend to minimize_scalar. You can simply
  790. pass a callable as the ``method`` parameter.
  791. The callable is called as ``method(fun, args, **kwargs, **options)``
  792. where ``kwargs`` corresponds to any other parameters passed to `minimize`
  793. (such as `bracket`, `tol`, etc.), except the `options` dict, which has
  794. its contents also passed as `method` parameters pair by pair. The method
  795. shall return an `OptimizeResult` object.
  796. The provided `method` callable must be able to accept (and possibly ignore)
  797. arbitrary parameters; the set of parameters accepted by `minimize` may
  798. expand in future versions and then these parameters will be passed to
  799. the method. You can find an example in the scipy.optimize tutorial.
  800. .. versionadded:: 0.11.0
  801. References
  802. ----------
  803. .. [1] Press, W., S.A. Teukolsky, W.T. Vetterling, and B.P. Flannery.
  804. Numerical Recipes in C. Cambridge University Press.
  805. .. [2] Forsythe, G.E., M. A. Malcolm, and C. B. Moler. "Computer Methods
  806. for Mathematical Computations." Prentice-Hall Series in Automatic
  807. Computation 259 (1977).
  808. .. [3] Brent, Richard P. Algorithms for Minimization Without Derivatives.
  809. Courier Corporation, 2013.
  810. Examples
  811. --------
  812. Consider the problem of minimizing the following function.
  813. >>> def f(x):
  814. ... return (x - 2) * x * (x + 2)**2
  815. Using the *Brent* method, we find the local minimum as:
  816. >>> from scipy.optimize import minimize_scalar
  817. >>> res = minimize_scalar(f)
  818. >>> res.fun
  819. -9.9149495908
  820. The minimizer is:
  821. >>> res.x
  822. 1.28077640403
  823. Using the *Bounded* method, we find a local minimum with specified
  824. bounds as:
  825. >>> res = minimize_scalar(f, bounds=(-3, -1), method='bounded')
  826. >>> res.fun # minimum
  827. 3.28365179850e-13
  828. >>> res.x # minimizer
  829. -2.0000002026
  830. """
  831. if not isinstance(args, tuple):
  832. args = (args,)
  833. if callable(method):
  834. meth = "_custom"
  835. elif method is None:
  836. meth = 'brent' if bounds is None else 'bounded'
  837. else:
  838. meth = method.lower()
  839. if options is None:
  840. options = {}
  841. if bounds is not None and meth in {'brent', 'golden'}:
  842. message = f"Use of `bounds` is incompatible with 'method={method}'."
  843. raise ValueError(message)
  844. if tol is not None:
  845. options = dict(options)
  846. if meth == 'bounded' and 'xatol' not in options:
  847. warn("Method 'bounded' does not support relative tolerance in x; "
  848. "defaulting to absolute tolerance.",
  849. RuntimeWarning, stacklevel=2)
  850. options['xatol'] = tol
  851. elif meth == '_custom':
  852. options.setdefault('tol', tol)
  853. else:
  854. options.setdefault('xtol', tol)
  855. # replace boolean "disp" option, if specified, by an integer value.
  856. disp = options.get('disp')
  857. if isinstance(disp, bool):
  858. options['disp'] = 2 * int(disp)
  859. if meth == '_custom':
  860. res = method(fun, args=args, bracket=bracket, bounds=bounds, **options)
  861. elif meth == 'brent':
  862. res = _recover_from_bracket_error(_minimize_scalar_brent,
  863. fun, bracket, args, **options)
  864. elif meth == 'bounded':
  865. if bounds is None:
  866. raise ValueError('The `bounds` parameter is mandatory for '
  867. 'method `bounded`.')
  868. res = _minimize_scalar_bounded(fun, bounds, args, **options)
  869. elif meth == 'golden':
  870. res = _recover_from_bracket_error(_minimize_scalar_golden,
  871. fun, bracket, args, **options)
  872. else:
  873. raise ValueError(f'Unknown solver {method}')
  874. # gh-16196 reported inconsistencies in the output shape of `res.x`. While
  875. # fixing this, future-proof it for when the function is vectorized:
  876. # the shape of `res.x` should match that of `res.fun`.
  877. res.fun = np.asarray(res.fun)[()]
  878. res.x = np.reshape(res.x, res.fun.shape)[()]
  879. return res
  880. def _remove_from_bounds(bounds, i_fixed):
  881. """Removes fixed variables from a `Bounds` instance"""
  882. lb = bounds.lb[~i_fixed]
  883. ub = bounds.ub[~i_fixed]
  884. return Bounds(lb, ub) # don't mutate original Bounds object
  885. class _Patch_Callback_Equal_Variables:
  886. # Patches a callback that accepts an intermediate_result
  887. def __init__(self, callback, i_fixed, x_fixed):
  888. self.callback = callback
  889. self.i_fixed = i_fixed
  890. self.x_fixed = x_fixed
  891. def __call__(self, intermediate_result):
  892. x_in = intermediate_result.x
  893. x_out = np.zeros_like(self.i_fixed, dtype=x_in.dtype)
  894. x_out[self.i_fixed] = self.x_fixed
  895. x_out[~self.i_fixed] = x_in
  896. intermediate_result.x = x_out
  897. return self.callback(intermediate_result)
  898. class _Remove_From_Func:
  899. """Wraps a function such that fixed variables need not be passed in"""
  900. def __init__(self, fun_in, i_fixed, x_fixed, min_dim=None, remove=0):
  901. self.fun_in = fun_in
  902. self.i_fixed = i_fixed
  903. self.x_fixed = x_fixed
  904. self.min_dim = min_dim
  905. self.remove = remove
  906. def __call__(self, x_in, *args, **kwargs):
  907. x_out = np.zeros_like(self.i_fixed, dtype=x_in.dtype)
  908. x_out[self.i_fixed] = self.x_fixed
  909. x_out[~self.i_fixed] = x_in
  910. y_out = self.fun_in(x_out, *args, **kwargs)
  911. y_out = np.array(y_out)
  912. if self.min_dim == 1:
  913. y_out = np.atleast_1d(y_out)
  914. elif self.min_dim == 2:
  915. y_out = np.atleast_2d(y_out)
  916. if self.remove == 1:
  917. y_out = y_out[..., ~self.i_fixed]
  918. elif self.remove == 2:
  919. y_out = y_out[~self.i_fixed, ~self.i_fixed]
  920. return y_out
  921. def _add_to_array(x_in, i_fixed, x_fixed):
  922. """Adds fixed variables back to an array"""
  923. i_free = ~i_fixed
  924. if x_in.ndim == 2:
  925. i_free = i_free[:, None] @ i_free[None, :]
  926. x_out = np.zeros_like(i_free, dtype=x_in.dtype)
  927. x_out[~i_free] = x_fixed
  928. x_out[i_free] = x_in.ravel()
  929. return x_out
  930. def _validate_bounds(bounds, x0, meth):
  931. """Check that bounds are valid."""
  932. msg = "An upper bound is less than the corresponding lower bound."
  933. if np.any(bounds.ub < bounds.lb):
  934. raise ValueError(msg)
  935. msg = "The number of bounds is not compatible with the length of `x0`."
  936. try:
  937. bounds.lb = np.broadcast_to(bounds.lb, x0.shape)
  938. bounds.ub = np.broadcast_to(bounds.ub, x0.shape)
  939. except Exception as e:
  940. raise ValueError(msg) from e
  941. return bounds
  942. def standardize_bounds(bounds, x0, meth):
  943. """Converts bounds to the form required by the solver."""
  944. if meth in {'trust-constr', 'powell', 'nelder-mead', 'cobyla', 'cobyqa',
  945. 'new'}:
  946. if not isinstance(bounds, Bounds):
  947. lb, ub = old_bound_to_new(bounds)
  948. bounds = Bounds(lb, ub)
  949. elif meth in ('l-bfgs-b', 'tnc', 'slsqp', 'old'):
  950. if isinstance(bounds, Bounds):
  951. bounds = new_bounds_to_old(bounds.lb, bounds.ub, x0.shape[0])
  952. return bounds
  953. def standardize_constraints(constraints, x0, meth):
  954. """Converts constraints to the form required by the solver."""
  955. all_constraint_types = (NonlinearConstraint, LinearConstraint, dict)
  956. new_constraint_types = all_constraint_types[:-1]
  957. if constraints is None:
  958. constraints = []
  959. elif isinstance(constraints, all_constraint_types):
  960. constraints = [constraints]
  961. else:
  962. constraints = list(constraints) # ensure it's a mutable sequence
  963. if meth in ['trust-constr', 'cobyqa', 'new', 'cobyla']:
  964. for i, con in enumerate(constraints):
  965. if not isinstance(con, new_constraint_types):
  966. constraints[i] = old_constraint_to_new(i, con)
  967. else:
  968. # iterate over copy, changing original
  969. for i, con in enumerate(list(constraints)):
  970. if isinstance(con, new_constraint_types):
  971. old_constraints = new_constraint_to_old(con, x0)
  972. constraints[i] = old_constraints[0]
  973. constraints.extend(old_constraints[1:]) # appends 1 if present
  974. return constraints
  975. def _optimize_result_for_equal_bounds(
  976. fun, bounds, method, args=(), constraints=()
  977. ):
  978. """
  979. Provides a default OptimizeResult for when a bounded minimization method
  980. has (lb == ub).all().
  981. Parameters
  982. ----------
  983. fun: callable
  984. bounds: Bounds
  985. method: str
  986. constraints: Constraint
  987. """
  988. success = True
  989. message = 'All independent variables were fixed by bounds.'
  990. # bounds is new-style
  991. x0 = bounds.lb
  992. if constraints:
  993. message = ("All independent variables were fixed by bounds at values"
  994. " that satisfy the constraints.")
  995. constraints = standardize_constraints(constraints, x0, 'new')
  996. maxcv = 0
  997. for c in constraints:
  998. pc = PreparedConstraint(c, x0)
  999. violation = pc.violation(x0)
  1000. if np.sum(violation):
  1001. maxcv = max(maxcv, np.max(violation))
  1002. success = False
  1003. message = (f"All independent variables were fixed by bounds, but "
  1004. f"the independent variables do not satisfy the "
  1005. f"constraints exactly. (Maximum violation: {maxcv}).")
  1006. return OptimizeResult(
  1007. x=x0, fun=fun(x0, *args), success=success, message=message, nfev=1,
  1008. njev=0, nhev=0,
  1009. )