optimize.py 103 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710711712713714715716717718719720721722723724725726727728729730731732733734735736737738739740741742743744745746747748749750751752753754755756757758759760761762763764765766767768769770771772773774775776777778779780781782783784785786787788789790791792793794795796797798799800801802803804805806807808809810811812813814815816817818819820821822823824825826827828829830831832833834835836837838839840841842843844845846847848849850851852853854855856857858859860861862863864865866867868869870871872873874875876877878879880881882883884885886887888889890891892893894895896897898899900901902903904905906907908909910911912913914915916917918919920921922923924925926927928929930931932933934935936937938939940941942943944945946947948949950951952953954955956957958959960961962963964965966967968969970971972973974975976977978979980981982983984985986987988989990991992993994995996997998999100010011002100310041005100610071008100910101011101210131014101510161017101810191020102110221023102410251026102710281029103010311032103310341035103610371038103910401041104210431044104510461047104810491050105110521053105410551056105710581059106010611062106310641065106610671068106910701071107210731074107510761077107810791080108110821083108410851086108710881089109010911092109310941095109610971098109911001101110211031104110511061107110811091110111111121113111411151116111711181119112011211122112311241125112611271128112911301131113211331134113511361137113811391140114111421143114411451146114711481149115011511152115311541155115611571158115911601161116211631164116511661167116811691170117111721173117411751176117711781179118011811182118311841185118611871188118911901191119211931194119511961197119811991200120112021203120412051206120712081209121012111212121312141215121612171218121912201221122212231224122512261227122812291230123112321233123412351236123712381239124012411242124312441245124612471248124912501251125212531254125512561257125812591260126112621263126412651266126712681269127012711272127312741275127612771278127912801281128212831284128512861287128812891290129112921293129412951296129712981299130013011302130313041305130613071308130913101311131213131314131513161317131813191320132113221323132413251326132713281329133013311332133313341335133613371338133913401341134213431344134513461347134813491350135113521353135413551356135713581359136013611362136313641365136613671368136913701371137213731374137513761377137813791380138113821383138413851386138713881389139013911392139313941395139613971398139914001401140214031404140514061407140814091410141114121413141414151416141714181419142014211422142314241425142614271428142914301431143214331434143514361437143814391440144114421443144414451446144714481449145014511452145314541455145614571458145914601461146214631464146514661467146814691470147114721473147414751476147714781479148014811482148314841485148614871488148914901491149214931494149514961497149814991500150115021503150415051506150715081509151015111512151315141515151615171518151915201521152215231524152515261527152815291530153115321533153415351536153715381539154015411542154315441545154615471548154915501551155215531554155515561557155815591560156115621563156415651566156715681569157015711572157315741575157615771578157915801581158215831584158515861587158815891590159115921593159415951596159715981599160016011602160316041605160616071608160916101611161216131614161516161617161816191620162116221623162416251626162716281629163016311632163316341635163616371638163916401641164216431644164516461647164816491650165116521653165416551656165716581659166016611662166316641665166616671668166916701671167216731674167516761677167816791680168116821683168416851686168716881689169016911692169316941695169616971698169917001701170217031704170517061707170817091710171117121713171417151716171717181719172017211722172317241725172617271728172917301731173217331734173517361737173817391740174117421743174417451746174717481749175017511752175317541755175617571758175917601761176217631764176517661767176817691770177117721773177417751776177717781779178017811782178317841785178617871788178917901791179217931794179517961797179817991800180118021803180418051806180718081809181018111812181318141815181618171818181918201821182218231824182518261827182818291830183118321833183418351836183718381839184018411842184318441845184618471848184918501851185218531854185518561857185818591860186118621863186418651866186718681869187018711872187318741875187618771878187918801881188218831884188518861887188818891890189118921893189418951896189718981899190019011902190319041905190619071908190919101911191219131914191519161917191819191920192119221923192419251926192719281929193019311932193319341935193619371938193919401941194219431944194519461947194819491950195119521953195419551956195719581959196019611962196319641965196619671968196919701971197219731974197519761977197819791980198119821983198419851986198719881989199019911992199319941995199619971998199920002001200220032004200520062007200820092010201120122013201420152016201720182019202020212022202320242025202620272028202920302031203220332034203520362037203820392040204120422043204420452046204720482049205020512052205320542055205620572058205920602061206220632064206520662067206820692070207120722073207420752076207720782079208020812082208320842085208620872088208920902091209220932094209520962097209820992100210121022103210421052106210721082109211021112112211321142115211621172118211921202121212221232124212521262127212821292130213121322133213421352136213721382139214021412142214321442145214621472148214921502151215221532154215521562157215821592160216121622163216421652166216721682169217021712172217321742175217621772178217921802181218221832184218521862187218821892190219121922193219421952196219721982199220022012202220322042205220622072208220922102211221222132214221522162217221822192220222122222223222422252226222722282229223022312232223322342235223622372238223922402241224222432244224522462247224822492250225122522253225422552256225722582259226022612262226322642265226622672268226922702271227222732274227522762277227822792280228122822283228422852286228722882289229022912292229322942295229622972298229923002301230223032304230523062307230823092310231123122313231423152316231723182319232023212322232323242325232623272328232923302331233223332334233523362337233823392340234123422343234423452346234723482349235023512352235323542355235623572358235923602361236223632364236523662367236823692370237123722373237423752376237723782379238023812382238323842385238623872388238923902391239223932394239523962397239823992400240124022403240424052406240724082409241024112412241324142415241624172418241924202421242224232424242524262427242824292430243124322433243424352436243724382439244024412442244324442445244624472448244924502451245224532454245524562457245824592460246124622463246424652466246724682469247024712472247324742475247624772478247924802481248224832484248524862487248824892490249124922493249424952496249724982499250025012502250325042505250625072508250925102511251225132514251525162517251825192520252125222523252425252526252725282529253025312532253325342535253625372538253925402541254225432544254525462547254825492550255125522553255425552556255725582559256025612562256325642565256625672568256925702571257225732574257525762577257825792580258125822583258425852586258725882589259025912592259325942595259625972598259926002601260226032604260526062607260826092610261126122613261426152616261726182619262026212622262326242625262626272628262926302631263226332634263526362637263826392640264126422643264426452646264726482649265026512652265326542655265626572658265926602661266226632664266526662667266826692670267126722673267426752676267726782679268026812682268326842685268626872688268926902691269226932694269526962697269826992700270127022703270427052706270727082709271027112712271327142715271627172718271927202721272227232724272527262727272827292730273127322733273427352736273727382739274027412742274327442745274627472748274927502751275227532754275527562757275827592760276127622763276427652766276727682769277027712772277327742775277627772778277927802781278227832784278527862787278827892790279127922793279427952796279727982799280028012802280328042805280628072808280928102811281228132814281528162817281828192820282128222823282428252826282728282829283028312832283328342835283628372838283928402841284228432844284528462847284828492850285128522853285428552856285728582859286028612862286328642865286628672868286928702871287228732874287528762877287828792880288128822883288428852886288728882889289028912892289328942895289628972898289929002901290229032904290529062907290829092910291129122913291429152916291729182919292029212922292329242925292629272928292929302931293229332934293529362937293829392940294129422943294429452946294729482949295029512952295329542955295629572958295929602961296229632964296529662967296829692970297129722973297429752976297729782979298029812982298329842985298629872988298929902991299229932994299529962997299829993000300130023003300430053006300730083009301030113012301330143015301630173018301930203021302230233024302530263027302830293030303130323033303430353036303730383039304030413042304330443045304630473048304930503051305230533054305530563057305830593060306130623063306430653066306730683069307030713072307330743075307630773078307930803081308230833084308530863087308830893090309130923093309430953096309730983099310031013102310331043105310631073108310931103111311231133114311531163117311831193120312131223123312431253126312731283129313031313132313331343135
  1. #__docformat__ = "restructuredtext en"
  2. # ******NOTICE***************
  3. # optimize.py module by Travis E. Oliphant
  4. #
  5. # You may copy and use this module as you see fit with no
  6. # guarantee implied provided you keep this notice in all copies.
  7. # *****END NOTICE************
  8. # A collection of optimization algorithms. Version 0.5
  9. # CHANGES
  10. # Added fminbound (July 2001)
  11. # Added brute (Aug. 2002)
  12. # Finished line search satisfying strong Wolfe conditions (Mar. 2004)
  13. # Updated strong Wolfe conditions line search to use
  14. # cubic-interpolation (Mar. 2004)
  15. from __future__ import division, print_function, absolute_import
  16. # Minimization routines
  17. __all__ = ['fmin', 'fmin_powell', 'fmin_bfgs', 'fmin_ncg', 'fmin_cg',
  18. 'fminbound', 'brent', 'golden', 'bracket', 'rosen', 'rosen_der',
  19. 'rosen_hess', 'rosen_hess_prod', 'brute', 'approx_fprime',
  20. 'line_search', 'check_grad', 'OptimizeResult', 'show_options',
  21. 'OptimizeWarning']
  22. __docformat__ = "restructuredtext en"
  23. import warnings
  24. import sys
  25. import numpy
  26. from scipy._lib.six import callable, xrange
  27. from numpy import (atleast_1d, eye, mgrid, argmin, zeros, shape, squeeze,
  28. vectorize, asarray, sqrt, Inf, asfarray, isinf)
  29. import numpy as np
  30. from .linesearch import (line_search_wolfe1, line_search_wolfe2,
  31. line_search_wolfe2 as line_search,
  32. LineSearchWarning)
  33. from scipy._lib._util import getargspec_no_self as _getargspec
  34. # standard status messages of optimizers
  35. _status_message = {'success': 'Optimization terminated successfully.',
  36. 'maxfev': 'Maximum number of function evaluations has '
  37. 'been exceeded.',
  38. 'maxiter': 'Maximum number of iterations has been '
  39. 'exceeded.',
  40. 'pr_loss': 'Desired error not necessarily achieved due '
  41. 'to precision loss.'}
  42. class MemoizeJac(object):
  43. """ Decorator that caches the value gradient of function each time it
  44. is called. """
  45. def __init__(self, fun):
  46. self.fun = fun
  47. self.jac = None
  48. self.x = None
  49. def __call__(self, x, *args):
  50. self.x = numpy.asarray(x).copy()
  51. fg = self.fun(x, *args)
  52. self.jac = fg[1]
  53. return fg[0]
  54. def derivative(self, x, *args):
  55. if self.jac is not None and numpy.alltrue(x == self.x):
  56. return self.jac
  57. else:
  58. self(x, *args)
  59. return self.jac
  60. class OptimizeResult(dict):
  61. """ Represents the optimization result.
  62. Attributes
  63. ----------
  64. x : ndarray
  65. The solution of the optimization.
  66. success : bool
  67. Whether or not the optimizer exited successfully.
  68. status : int
  69. Termination status of the optimizer. Its value depends on the
  70. underlying solver. Refer to `message` for details.
  71. message : str
  72. Description of the cause of the termination.
  73. fun, jac, hess: ndarray
  74. Values of objective function, its Jacobian and its Hessian (if
  75. available). The Hessians may be approximations, see the documentation
  76. of the function in question.
  77. hess_inv : object
  78. Inverse of the objective function's Hessian; may be an approximation.
  79. Not available for all solvers. The type of this attribute may be
  80. either np.ndarray or scipy.sparse.linalg.LinearOperator.
  81. nfev, njev, nhev : int
  82. Number of evaluations of the objective functions and of its
  83. Jacobian and Hessian.
  84. nit : int
  85. Number of iterations performed by the optimizer.
  86. maxcv : float
  87. The maximum constraint violation.
  88. Notes
  89. -----
  90. There may be additional attributes not listed above depending of the
  91. specific solver. Since this class is essentially a subclass of dict
  92. with attribute accessors, one can see which attributes are available
  93. using the `keys()` method.
  94. """
  95. def __getattr__(self, name):
  96. try:
  97. return self[name]
  98. except KeyError:
  99. raise AttributeError(name)
  100. __setattr__ = dict.__setitem__
  101. __delattr__ = dict.__delitem__
  102. def __repr__(self):
  103. if self.keys():
  104. m = max(map(len, list(self.keys()))) + 1
  105. return '\n'.join([k.rjust(m) + ': ' + repr(v)
  106. for k, v in sorted(self.items())])
  107. else:
  108. return self.__class__.__name__ + "()"
  109. def __dir__(self):
  110. return list(self.keys())
  111. class OptimizeWarning(UserWarning):
  112. pass
  113. def _check_unknown_options(unknown_options):
  114. if unknown_options:
  115. msg = ", ".join(map(str, unknown_options.keys()))
  116. # Stack level 4: this is called from _minimize_*, which is
  117. # called from another function in Scipy. Level 4 is the first
  118. # level in user code.
  119. warnings.warn("Unknown solver options: %s" % msg, OptimizeWarning, 4)
  120. def is_array_scalar(x):
  121. """Test whether `x` is either a scalar or an array scalar.
  122. """
  123. return np.size(x) == 1
  124. _epsilon = sqrt(numpy.finfo(float).eps)
  125. def vecnorm(x, ord=2):
  126. if ord == Inf:
  127. return numpy.amax(numpy.abs(x))
  128. elif ord == -Inf:
  129. return numpy.amin(numpy.abs(x))
  130. else:
  131. return numpy.sum(numpy.abs(x)**ord, axis=0)**(1.0 / ord)
  132. def rosen(x):
  133. """
  134. The Rosenbrock function.
  135. The function computed is::
  136. sum(100.0*(x[1:] - x[:-1]**2.0)**2.0 + (1 - x[:-1])**2.0)
  137. Parameters
  138. ----------
  139. x : array_like
  140. 1-D array of points at which the Rosenbrock function is to be computed.
  141. Returns
  142. -------
  143. f : float
  144. The value of the Rosenbrock function.
  145. See Also
  146. --------
  147. rosen_der, rosen_hess, rosen_hess_prod
  148. Examples
  149. --------
  150. >>> from scipy.optimize import rosen
  151. >>> X = 0.1 * np.arange(10)
  152. >>> rosen(X)
  153. 76.56
  154. """
  155. x = asarray(x)
  156. r = numpy.sum(100.0 * (x[1:] - x[:-1]**2.0)**2.0 + (1 - x[:-1])**2.0,
  157. axis=0)
  158. return r
  159. def rosen_der(x):
  160. """
  161. The derivative (i.e. gradient) of the Rosenbrock function.
  162. Parameters
  163. ----------
  164. x : array_like
  165. 1-D array of points at which the derivative is to be computed.
  166. Returns
  167. -------
  168. rosen_der : (N,) ndarray
  169. The gradient of the Rosenbrock function at `x`.
  170. See Also
  171. --------
  172. rosen, rosen_hess, rosen_hess_prod
  173. """
  174. x = asarray(x)
  175. xm = x[1:-1]
  176. xm_m1 = x[:-2]
  177. xm_p1 = x[2:]
  178. der = numpy.zeros_like(x)
  179. der[1:-1] = (200 * (xm - xm_m1**2) -
  180. 400 * (xm_p1 - xm**2) * xm - 2 * (1 - xm))
  181. der[0] = -400 * x[0] * (x[1] - x[0]**2) - 2 * (1 - x[0])
  182. der[-1] = 200 * (x[-1] - x[-2]**2)
  183. return der
  184. def rosen_hess(x):
  185. """
  186. The Hessian matrix of the Rosenbrock function.
  187. Parameters
  188. ----------
  189. x : array_like
  190. 1-D array of points at which the Hessian matrix is to be computed.
  191. Returns
  192. -------
  193. rosen_hess : ndarray
  194. The Hessian matrix of the Rosenbrock function at `x`.
  195. See Also
  196. --------
  197. rosen, rosen_der, rosen_hess_prod
  198. """
  199. x = atleast_1d(x)
  200. H = numpy.diag(-400 * x[:-1], 1) - numpy.diag(400 * x[:-1], -1)
  201. diagonal = numpy.zeros(len(x), dtype=x.dtype)
  202. diagonal[0] = 1200 * x[0]**2 - 400 * x[1] + 2
  203. diagonal[-1] = 200
  204. diagonal[1:-1] = 202 + 1200 * x[1:-1]**2 - 400 * x[2:]
  205. H = H + numpy.diag(diagonal)
  206. return H
  207. def rosen_hess_prod(x, p):
  208. """
  209. Product of the Hessian matrix of the Rosenbrock function with a vector.
  210. Parameters
  211. ----------
  212. x : array_like
  213. 1-D array of points at which the Hessian matrix is to be computed.
  214. p : array_like
  215. 1-D array, the vector to be multiplied by the Hessian matrix.
  216. Returns
  217. -------
  218. rosen_hess_prod : ndarray
  219. The Hessian matrix of the Rosenbrock function at `x` multiplied
  220. by the vector `p`.
  221. See Also
  222. --------
  223. rosen, rosen_der, rosen_hess
  224. """
  225. x = atleast_1d(x)
  226. Hp = numpy.zeros(len(x), dtype=x.dtype)
  227. Hp[0] = (1200 * x[0]**2 - 400 * x[1] + 2) * p[0] - 400 * x[0] * p[1]
  228. Hp[1:-1] = (-400 * x[:-2] * p[:-2] +
  229. (202 + 1200 * x[1:-1]**2 - 400 * x[2:]) * p[1:-1] -
  230. 400 * x[1:-1] * p[2:])
  231. Hp[-1] = -400 * x[-2] * p[-2] + 200*p[-1]
  232. return Hp
  233. def wrap_function(function, args):
  234. ncalls = [0]
  235. if function is None:
  236. return ncalls, None
  237. def function_wrapper(*wrapper_args):
  238. ncalls[0] += 1
  239. return function(*(wrapper_args + args))
  240. return ncalls, function_wrapper
  241. def fmin(func, x0, args=(), xtol=1e-4, ftol=1e-4, maxiter=None, maxfun=None,
  242. full_output=0, disp=1, retall=0, callback=None, initial_simplex=None):
  243. """
  244. Minimize a function using the downhill simplex algorithm.
  245. This algorithm only uses function values, not derivatives or second
  246. derivatives.
  247. Parameters
  248. ----------
  249. func : callable func(x,*args)
  250. The objective function to be minimized.
  251. x0 : ndarray
  252. Initial guess.
  253. args : tuple, optional
  254. Extra arguments passed to func, i.e. ``f(x,*args)``.
  255. xtol : float, optional
  256. Absolute error in xopt between iterations that is acceptable for
  257. convergence.
  258. ftol : number, optional
  259. Absolute error in func(xopt) between iterations that is acceptable for
  260. convergence.
  261. maxiter : int, optional
  262. Maximum number of iterations to perform.
  263. maxfun : number, optional
  264. Maximum number of function evaluations to make.
  265. full_output : bool, optional
  266. Set to True if fopt and warnflag outputs are desired.
  267. disp : bool, optional
  268. Set to True to print convergence messages.
  269. retall : bool, optional
  270. Set to True to return list of solutions at each iteration.
  271. callback : callable, optional
  272. Called after each iteration, as callback(xk), where xk is the
  273. current parameter vector.
  274. initial_simplex : array_like of shape (N + 1, N), optional
  275. Initial simplex. If given, overrides `x0`.
  276. ``initial_simplex[j,:]`` should contain the coordinates of
  277. the j-th vertex of the ``N+1`` vertices in the simplex, where
  278. ``N`` is the dimension.
  279. Returns
  280. -------
  281. xopt : ndarray
  282. Parameter that minimizes function.
  283. fopt : float
  284. Value of function at minimum: ``fopt = func(xopt)``.
  285. iter : int
  286. Number of iterations performed.
  287. funcalls : int
  288. Number of function calls made.
  289. warnflag : int
  290. 1 : Maximum number of function evaluations made.
  291. 2 : Maximum number of iterations reached.
  292. allvecs : list
  293. Solution at each iteration.
  294. See also
  295. --------
  296. minimize: Interface to minimization algorithms for multivariate
  297. functions. See the 'Nelder-Mead' `method` in particular.
  298. Notes
  299. -----
  300. Uses a Nelder-Mead simplex algorithm to find the minimum of function of
  301. one or more variables.
  302. This algorithm has a long history of successful use in applications.
  303. But it will usually be slower than an algorithm that uses first or
  304. second derivative information. In practice it can have poor
  305. performance in high-dimensional problems and is not robust to
  306. minimizing complicated functions. Additionally, there currently is no
  307. complete theory describing when the algorithm will successfully
  308. converge to the minimum, or how fast it will if it does. Both the ftol and
  309. xtol criteria must be met for convergence.
  310. Examples
  311. --------
  312. >>> def f(x):
  313. ... return x**2
  314. >>> from scipy import optimize
  315. >>> minimum = optimize.fmin(f, 1)
  316. Optimization terminated successfully.
  317. Current function value: 0.000000
  318. Iterations: 17
  319. Function evaluations: 34
  320. >>> minimum[0]
  321. -8.8817841970012523e-16
  322. References
  323. ----------
  324. .. [1] Nelder, J.A. and Mead, R. (1965), "A simplex method for function
  325. minimization", The Computer Journal, 7, pp. 308-313
  326. .. [2] Wright, M.H. (1996), "Direct Search Methods: Once Scorned, Now
  327. Respectable", in Numerical Analysis 1995, Proceedings of the
  328. 1995 Dundee Biennial Conference in Numerical Analysis, D.F.
  329. Griffiths and G.A. Watson (Eds.), Addison Wesley Longman,
  330. Harlow, UK, pp. 191-208.
  331. """
  332. opts = {'xatol': xtol,
  333. 'fatol': ftol,
  334. 'maxiter': maxiter,
  335. 'maxfev': maxfun,
  336. 'disp': disp,
  337. 'return_all': retall,
  338. 'initial_simplex': initial_simplex}
  339. res = _minimize_neldermead(func, x0, args, callback=callback, **opts)
  340. if full_output:
  341. retlist = res['x'], res['fun'], res['nit'], res['nfev'], res['status']
  342. if retall:
  343. retlist += (res['allvecs'], )
  344. return retlist
  345. else:
  346. if retall:
  347. return res['x'], res['allvecs']
  348. else:
  349. return res['x']
  350. def _minimize_neldermead(func, x0, args=(), callback=None,
  351. maxiter=None, maxfev=None, disp=False,
  352. return_all=False, initial_simplex=None,
  353. xatol=1e-4, fatol=1e-4, adaptive=False,
  354. **unknown_options):
  355. """
  356. Minimization of scalar function of one or more variables using the
  357. Nelder-Mead algorithm.
  358. Options
  359. -------
  360. disp : bool
  361. Set to True to print convergence messages.
  362. maxiter, maxfev : int
  363. Maximum allowed number of iterations and function evaluations.
  364. Will default to ``N*200``, where ``N`` is the number of
  365. variables, if neither `maxiter` or `maxfev` is set. If both
  366. `maxiter` and `maxfev` are set, minimization will stop at the
  367. first reached.
  368. initial_simplex : array_like of shape (N + 1, N)
  369. Initial simplex. If given, overrides `x0`.
  370. ``initial_simplex[j,:]`` should contain the coordinates of
  371. the j-th vertex of the ``N+1`` vertices in the simplex, where
  372. ``N`` is the dimension.
  373. xatol : float, optional
  374. Absolute error in xopt between iterations that is acceptable for
  375. convergence.
  376. fatol : number, optional
  377. Absolute error in func(xopt) between iterations that is acceptable for
  378. convergence.
  379. adaptive : bool, optional
  380. Adapt algorithm parameters to dimensionality of problem. Useful for
  381. high-dimensional minimization [1]_.
  382. References
  383. ----------
  384. .. [1] Gao, F. and Han, L.
  385. Implementing the Nelder-Mead simplex algorithm with adaptive
  386. parameters. 2012. Computational Optimization and Applications.
  387. 51:1, pp. 259-277
  388. """
  389. if 'ftol' in unknown_options:
  390. warnings.warn("ftol is deprecated for Nelder-Mead,"
  391. " use fatol instead. If you specified both, only"
  392. " fatol is used.",
  393. DeprecationWarning)
  394. if (np.isclose(fatol, 1e-4) and
  395. not np.isclose(unknown_options['ftol'], 1e-4)):
  396. # only ftol was probably specified, use it.
  397. fatol = unknown_options['ftol']
  398. unknown_options.pop('ftol')
  399. if 'xtol' in unknown_options:
  400. warnings.warn("xtol is deprecated for Nelder-Mead,"
  401. " use xatol instead. If you specified both, only"
  402. " xatol is used.",
  403. DeprecationWarning)
  404. if (np.isclose(xatol, 1e-4) and
  405. not np.isclose(unknown_options['xtol'], 1e-4)):
  406. # only xtol was probably specified, use it.
  407. xatol = unknown_options['xtol']
  408. unknown_options.pop('xtol')
  409. _check_unknown_options(unknown_options)
  410. maxfun = maxfev
  411. retall = return_all
  412. fcalls, func = wrap_function(func, args)
  413. if adaptive:
  414. dim = float(len(x0))
  415. rho = 1
  416. chi = 1 + 2/dim
  417. psi = 0.75 - 1/(2*dim)
  418. sigma = 1 - 1/dim
  419. else:
  420. rho = 1
  421. chi = 2
  422. psi = 0.5
  423. sigma = 0.5
  424. nonzdelt = 0.05
  425. zdelt = 0.00025
  426. x0 = asfarray(x0).flatten()
  427. if initial_simplex is None:
  428. N = len(x0)
  429. sim = numpy.zeros((N + 1, N), dtype=x0.dtype)
  430. sim[0] = x0
  431. for k in range(N):
  432. y = numpy.array(x0, copy=True)
  433. if y[k] != 0:
  434. y[k] = (1 + nonzdelt)*y[k]
  435. else:
  436. y[k] = zdelt
  437. sim[k + 1] = y
  438. else:
  439. sim = np.asfarray(initial_simplex).copy()
  440. if sim.ndim != 2 or sim.shape[0] != sim.shape[1] + 1:
  441. raise ValueError("`initial_simplex` should be an array of shape (N+1,N)")
  442. if len(x0) != sim.shape[1]:
  443. raise ValueError("Size of `initial_simplex` is not consistent with `x0`")
  444. N = sim.shape[1]
  445. if retall:
  446. allvecs = [sim[0]]
  447. # If neither are set, then set both to default
  448. if maxiter is None and maxfun is None:
  449. maxiter = N * 200
  450. maxfun = N * 200
  451. elif maxiter is None:
  452. # Convert remaining Nones, to np.inf, unless the other is np.inf, in
  453. # which case use the default to avoid unbounded iteration
  454. if maxfun == np.inf:
  455. maxiter = N * 200
  456. else:
  457. maxiter = np.inf
  458. elif maxfun is None:
  459. if maxiter == np.inf:
  460. maxfun = N * 200
  461. else:
  462. maxfun = np.inf
  463. one2np1 = list(range(1, N + 1))
  464. fsim = numpy.zeros((N + 1,), float)
  465. for k in range(N + 1):
  466. fsim[k] = func(sim[k])
  467. ind = numpy.argsort(fsim)
  468. fsim = numpy.take(fsim, ind, 0)
  469. # sort so sim[0,:] has the lowest function value
  470. sim = numpy.take(sim, ind, 0)
  471. iterations = 1
  472. while (fcalls[0] < maxfun and iterations < maxiter):
  473. if (numpy.max(numpy.ravel(numpy.abs(sim[1:] - sim[0]))) <= xatol and
  474. numpy.max(numpy.abs(fsim[0] - fsim[1:])) <= fatol):
  475. break
  476. xbar = numpy.add.reduce(sim[:-1], 0) / N
  477. xr = (1 + rho) * xbar - rho * sim[-1]
  478. fxr = func(xr)
  479. doshrink = 0
  480. if fxr < fsim[0]:
  481. xe = (1 + rho * chi) * xbar - rho * chi * sim[-1]
  482. fxe = func(xe)
  483. if fxe < fxr:
  484. sim[-1] = xe
  485. fsim[-1] = fxe
  486. else:
  487. sim[-1] = xr
  488. fsim[-1] = fxr
  489. else: # fsim[0] <= fxr
  490. if fxr < fsim[-2]:
  491. sim[-1] = xr
  492. fsim[-1] = fxr
  493. else: # fxr >= fsim[-2]
  494. # Perform contraction
  495. if fxr < fsim[-1]:
  496. xc = (1 + psi * rho) * xbar - psi * rho * sim[-1]
  497. fxc = func(xc)
  498. if fxc <= fxr:
  499. sim[-1] = xc
  500. fsim[-1] = fxc
  501. else:
  502. doshrink = 1
  503. else:
  504. # Perform an inside contraction
  505. xcc = (1 - psi) * xbar + psi * sim[-1]
  506. fxcc = func(xcc)
  507. if fxcc < fsim[-1]:
  508. sim[-1] = xcc
  509. fsim[-1] = fxcc
  510. else:
  511. doshrink = 1
  512. if doshrink:
  513. for j in one2np1:
  514. sim[j] = sim[0] + sigma * (sim[j] - sim[0])
  515. fsim[j] = func(sim[j])
  516. ind = numpy.argsort(fsim)
  517. sim = numpy.take(sim, ind, 0)
  518. fsim = numpy.take(fsim, ind, 0)
  519. if callback is not None:
  520. callback(sim[0])
  521. iterations += 1
  522. if retall:
  523. allvecs.append(sim[0])
  524. x = sim[0]
  525. fval = numpy.min(fsim)
  526. warnflag = 0
  527. if fcalls[0] >= maxfun:
  528. warnflag = 1
  529. msg = _status_message['maxfev']
  530. if disp:
  531. print('Warning: ' + msg)
  532. elif iterations >= maxiter:
  533. warnflag = 2
  534. msg = _status_message['maxiter']
  535. if disp:
  536. print('Warning: ' + msg)
  537. else:
  538. msg = _status_message['success']
  539. if disp:
  540. print(msg)
  541. print(" Current function value: %f" % fval)
  542. print(" Iterations: %d" % iterations)
  543. print(" Function evaluations: %d" % fcalls[0])
  544. result = OptimizeResult(fun=fval, nit=iterations, nfev=fcalls[0],
  545. status=warnflag, success=(warnflag == 0),
  546. message=msg, x=x, final_simplex=(sim, fsim))
  547. if retall:
  548. result['allvecs'] = allvecs
  549. return result
  550. def _approx_fprime_helper(xk, f, epsilon, args=(), f0=None):
  551. """
  552. See ``approx_fprime``. An optional initial function value arg is added.
  553. """
  554. if f0 is None:
  555. f0 = f(*((xk,) + args))
  556. grad = numpy.zeros((len(xk),), float)
  557. ei = numpy.zeros((len(xk),), float)
  558. for k in range(len(xk)):
  559. ei[k] = 1.0
  560. d = epsilon * ei
  561. grad[k] = (f(*((xk + d,) + args)) - f0) / d[k]
  562. ei[k] = 0.0
  563. return grad
  564. def approx_fprime(xk, f, epsilon, *args):
  565. """Finite-difference approximation of the gradient of a scalar function.
  566. Parameters
  567. ----------
  568. xk : array_like
  569. The coordinate vector at which to determine the gradient of `f`.
  570. f : callable
  571. The function of which to determine the gradient (partial derivatives).
  572. Should take `xk` as first argument, other arguments to `f` can be
  573. supplied in ``*args``. Should return a scalar, the value of the
  574. function at `xk`.
  575. epsilon : array_like
  576. Increment to `xk` to use for determining the function gradient.
  577. If a scalar, uses the same finite difference delta for all partial
  578. derivatives. If an array, should contain one value per element of
  579. `xk`.
  580. \\*args : args, optional
  581. Any other arguments that are to be passed to `f`.
  582. Returns
  583. -------
  584. grad : ndarray
  585. The partial derivatives of `f` to `xk`.
  586. See Also
  587. --------
  588. check_grad : Check correctness of gradient function against approx_fprime.
  589. Notes
  590. -----
  591. The function gradient is determined by the forward finite difference
  592. formula::
  593. f(xk[i] + epsilon[i]) - f(xk[i])
  594. f'[i] = ---------------------------------
  595. epsilon[i]
  596. The main use of `approx_fprime` is in scalar function optimizers like
  597. `fmin_bfgs`, to determine numerically the Jacobian of a function.
  598. Examples
  599. --------
  600. >>> from scipy import optimize
  601. >>> def func(x, c0, c1):
  602. ... "Coordinate vector `x` should be an array of size two."
  603. ... return c0 * x[0]**2 + c1*x[1]**2
  604. >>> x = np.ones(2)
  605. >>> c0, c1 = (1, 200)
  606. >>> eps = np.sqrt(np.finfo(float).eps)
  607. >>> optimize.approx_fprime(x, func, [eps, np.sqrt(200) * eps], c0, c1)
  608. array([ 2. , 400.00004198])
  609. """
  610. return _approx_fprime_helper(xk, f, epsilon, args=args)
  611. def check_grad(func, grad, x0, *args, **kwargs):
  612. """Check the correctness of a gradient function by comparing it against a
  613. (forward) finite-difference approximation of the gradient.
  614. Parameters
  615. ----------
  616. func : callable ``func(x0, *args)``
  617. Function whose derivative is to be checked.
  618. grad : callable ``grad(x0, *args)``
  619. Gradient of `func`.
  620. x0 : ndarray
  621. Points to check `grad` against forward difference approximation of grad
  622. using `func`.
  623. args : \\*args, optional
  624. Extra arguments passed to `func` and `grad`.
  625. epsilon : float, optional
  626. Step size used for the finite difference approximation. It defaults to
  627. ``sqrt(numpy.finfo(float).eps)``, which is approximately 1.49e-08.
  628. Returns
  629. -------
  630. err : float
  631. The square root of the sum of squares (i.e. the 2-norm) of the
  632. difference between ``grad(x0, *args)`` and the finite difference
  633. approximation of `grad` using func at the points `x0`.
  634. See Also
  635. --------
  636. approx_fprime
  637. Examples
  638. --------
  639. >>> def func(x):
  640. ... return x[0]**2 - 0.5 * x[1]**3
  641. >>> def grad(x):
  642. ... return [2 * x[0], -1.5 * x[1]**2]
  643. >>> from scipy.optimize import check_grad
  644. >>> check_grad(func, grad, [1.5, -1.5])
  645. 2.9802322387695312e-08
  646. """
  647. step = kwargs.pop('epsilon', _epsilon)
  648. if kwargs:
  649. raise ValueError("Unknown keyword arguments: %r" %
  650. (list(kwargs.keys()),))
  651. return sqrt(sum((grad(x0, *args) -
  652. approx_fprime(x0, func, step, *args))**2))
  653. def approx_fhess_p(x0, p, fprime, epsilon, *args):
  654. f2 = fprime(*((x0 + epsilon*p,) + args))
  655. f1 = fprime(*((x0,) + args))
  656. return (f2 - f1) / epsilon
  657. class _LineSearchError(RuntimeError):
  658. pass
  659. def _line_search_wolfe12(f, fprime, xk, pk, gfk, old_fval, old_old_fval,
  660. **kwargs):
  661. """
  662. Same as line_search_wolfe1, but fall back to line_search_wolfe2 if
  663. suitable step length is not found, and raise an exception if a
  664. suitable step length is not found.
  665. Raises
  666. ------
  667. _LineSearchError
  668. If no suitable step size is found
  669. """
  670. extra_condition = kwargs.pop('extra_condition', None)
  671. ret = line_search_wolfe1(f, fprime, xk, pk, gfk,
  672. old_fval, old_old_fval,
  673. **kwargs)
  674. if ret[0] is not None and extra_condition is not None:
  675. xp1 = xk + ret[0] * pk
  676. if not extra_condition(ret[0], xp1, ret[3], ret[5]):
  677. # Reject step if extra_condition fails
  678. ret = (None,)
  679. if ret[0] is None:
  680. # line search failed: try different one.
  681. with warnings.catch_warnings():
  682. warnings.simplefilter('ignore', LineSearchWarning)
  683. kwargs2 = {}
  684. for key in ('c1', 'c2', 'amax'):
  685. if key in kwargs:
  686. kwargs2[key] = kwargs[key]
  687. ret = line_search_wolfe2(f, fprime, xk, pk, gfk,
  688. old_fval, old_old_fval,
  689. extra_condition=extra_condition,
  690. **kwargs2)
  691. if ret[0] is None:
  692. raise _LineSearchError()
  693. return ret
  694. def fmin_bfgs(f, x0, fprime=None, args=(), gtol=1e-5, norm=Inf,
  695. epsilon=_epsilon, maxiter=None, full_output=0, disp=1,
  696. retall=0, callback=None):
  697. """
  698. Minimize a function using the BFGS algorithm.
  699. Parameters
  700. ----------
  701. f : callable f(x,*args)
  702. Objective function to be minimized.
  703. x0 : ndarray
  704. Initial guess.
  705. fprime : callable f'(x,*args), optional
  706. Gradient of f.
  707. args : tuple, optional
  708. Extra arguments passed to f and fprime.
  709. gtol : float, optional
  710. Gradient norm must be less than gtol before successful termination.
  711. norm : float, optional
  712. Order of norm (Inf is max, -Inf is min)
  713. epsilon : int or ndarray, optional
  714. If fprime is approximated, use this value for the step size.
  715. callback : callable, optional
  716. An optional user-supplied function to call after each
  717. iteration. Called as callback(xk), where xk is the
  718. current parameter vector.
  719. maxiter : int, optional
  720. Maximum number of iterations to perform.
  721. full_output : bool, optional
  722. If True,return fopt, func_calls, grad_calls, and warnflag
  723. in addition to xopt.
  724. disp : bool, optional
  725. Print convergence message if True.
  726. retall : bool, optional
  727. Return a list of results at each iteration if True.
  728. Returns
  729. -------
  730. xopt : ndarray
  731. Parameters which minimize f, i.e. f(xopt) == fopt.
  732. fopt : float
  733. Minimum value.
  734. gopt : ndarray
  735. Value of gradient at minimum, f'(xopt), which should be near 0.
  736. Bopt : ndarray
  737. Value of 1/f''(xopt), i.e. the inverse hessian matrix.
  738. func_calls : int
  739. Number of function_calls made.
  740. grad_calls : int
  741. Number of gradient calls made.
  742. warnflag : integer
  743. 1 : Maximum number of iterations exceeded.
  744. 2 : Gradient and/or function calls not changing.
  745. allvecs : list
  746. The value of xopt at each iteration. Only returned if retall is True.
  747. See also
  748. --------
  749. minimize: Interface to minimization algorithms for multivariate
  750. functions. See the 'BFGS' `method` in particular.
  751. Notes
  752. -----
  753. Optimize the function, f, whose gradient is given by fprime
  754. using the quasi-Newton method of Broyden, Fletcher, Goldfarb,
  755. and Shanno (BFGS)
  756. References
  757. ----------
  758. Wright, and Nocedal 'Numerical Optimization', 1999, pg. 198.
  759. """
  760. opts = {'gtol': gtol,
  761. 'norm': norm,
  762. 'eps': epsilon,
  763. 'disp': disp,
  764. 'maxiter': maxiter,
  765. 'return_all': retall}
  766. res = _minimize_bfgs(f, x0, args, fprime, callback=callback, **opts)
  767. if full_output:
  768. retlist = (res['x'], res['fun'], res['jac'], res['hess_inv'],
  769. res['nfev'], res['njev'], res['status'])
  770. if retall:
  771. retlist += (res['allvecs'], )
  772. return retlist
  773. else:
  774. if retall:
  775. return res['x'], res['allvecs']
  776. else:
  777. return res['x']
  778. def _minimize_bfgs(fun, x0, args=(), jac=None, callback=None,
  779. gtol=1e-5, norm=Inf, eps=_epsilon, maxiter=None,
  780. disp=False, return_all=False,
  781. **unknown_options):
  782. """
  783. Minimization of scalar function of one or more variables using the
  784. BFGS algorithm.
  785. Options
  786. -------
  787. disp : bool
  788. Set to True to print convergence messages.
  789. maxiter : int
  790. Maximum number of iterations to perform.
  791. gtol : float
  792. Gradient norm must be less than `gtol` before successful
  793. termination.
  794. norm : float
  795. Order of norm (Inf is max, -Inf is min).
  796. eps : float or ndarray
  797. If `jac` is approximated, use this value for the step size.
  798. """
  799. _check_unknown_options(unknown_options)
  800. f = fun
  801. fprime = jac
  802. epsilon = eps
  803. retall = return_all
  804. x0 = asarray(x0).flatten()
  805. if x0.ndim == 0:
  806. x0.shape = (1,)
  807. if maxiter is None:
  808. maxiter = len(x0) * 200
  809. func_calls, f = wrap_function(f, args)
  810. if fprime is None:
  811. grad_calls, myfprime = wrap_function(approx_fprime, (f, epsilon))
  812. else:
  813. grad_calls, myfprime = wrap_function(fprime, args)
  814. gfk = myfprime(x0)
  815. k = 0
  816. N = len(x0)
  817. I = numpy.eye(N, dtype=int)
  818. Hk = I
  819. # Sets the initial step guess to dx ~ 1
  820. old_fval = f(x0)
  821. old_old_fval = old_fval + np.linalg.norm(gfk) / 2
  822. xk = x0
  823. if retall:
  824. allvecs = [x0]
  825. warnflag = 0
  826. gnorm = vecnorm(gfk, ord=norm)
  827. while (gnorm > gtol) and (k < maxiter):
  828. pk = -numpy.dot(Hk, gfk)
  829. try:
  830. alpha_k, fc, gc, old_fval, old_old_fval, gfkp1 = \
  831. _line_search_wolfe12(f, myfprime, xk, pk, gfk,
  832. old_fval, old_old_fval, amin=1e-100, amax=1e100)
  833. except _LineSearchError:
  834. # Line search failed to find a better solution.
  835. warnflag = 2
  836. break
  837. xkp1 = xk + alpha_k * pk
  838. if retall:
  839. allvecs.append(xkp1)
  840. sk = xkp1 - xk
  841. xk = xkp1
  842. if gfkp1 is None:
  843. gfkp1 = myfprime(xkp1)
  844. yk = gfkp1 - gfk
  845. gfk = gfkp1
  846. if callback is not None:
  847. callback(xk)
  848. k += 1
  849. gnorm = vecnorm(gfk, ord=norm)
  850. if (gnorm <= gtol):
  851. break
  852. if not numpy.isfinite(old_fval):
  853. # We correctly found +-Inf as optimal value, or something went
  854. # wrong.
  855. warnflag = 2
  856. break
  857. try: # this was handled in numeric, let it remaines for more safety
  858. rhok = 1.0 / (numpy.dot(yk, sk))
  859. except ZeroDivisionError:
  860. rhok = 1000.0
  861. if disp:
  862. print("Divide-by-zero encountered: rhok assumed large")
  863. if isinf(rhok): # this is patch for numpy
  864. rhok = 1000.0
  865. if disp:
  866. print("Divide-by-zero encountered: rhok assumed large")
  867. A1 = I - sk[:, numpy.newaxis] * yk[numpy.newaxis, :] * rhok
  868. A2 = I - yk[:, numpy.newaxis] * sk[numpy.newaxis, :] * rhok
  869. Hk = numpy.dot(A1, numpy.dot(Hk, A2)) + (rhok * sk[:, numpy.newaxis] *
  870. sk[numpy.newaxis, :])
  871. fval = old_fval
  872. if np.isnan(fval):
  873. # This can happen if the first call to f returned NaN;
  874. # the loop is then never entered.
  875. warnflag = 2
  876. if warnflag == 2:
  877. msg = _status_message['pr_loss']
  878. elif k >= maxiter:
  879. warnflag = 1
  880. msg = _status_message['maxiter']
  881. else:
  882. msg = _status_message['success']
  883. if disp:
  884. print("%s%s" % ("Warning: " if warnflag != 0 else "", msg))
  885. print(" Current function value: %f" % fval)
  886. print(" Iterations: %d" % k)
  887. print(" Function evaluations: %d" % func_calls[0])
  888. print(" Gradient evaluations: %d" % grad_calls[0])
  889. result = OptimizeResult(fun=fval, jac=gfk, hess_inv=Hk, nfev=func_calls[0],
  890. njev=grad_calls[0], status=warnflag,
  891. success=(warnflag == 0), message=msg, x=xk,
  892. nit=k)
  893. if retall:
  894. result['allvecs'] = allvecs
  895. return result
  896. def fmin_cg(f, x0, fprime=None, args=(), gtol=1e-5, norm=Inf, epsilon=_epsilon,
  897. maxiter=None, full_output=0, disp=1, retall=0, callback=None):
  898. """
  899. Minimize a function using a nonlinear conjugate gradient algorithm.
  900. Parameters
  901. ----------
  902. f : callable, ``f(x, *args)``
  903. Objective function to be minimized. Here `x` must be a 1-D array of
  904. the variables that are to be changed in the search for a minimum, and
  905. `args` are the other (fixed) parameters of `f`.
  906. x0 : ndarray
  907. A user-supplied initial estimate of `xopt`, the optimal value of `x`.
  908. It must be a 1-D array of values.
  909. fprime : callable, ``fprime(x, *args)``, optional
  910. A function that returns the gradient of `f` at `x`. Here `x` and `args`
  911. are as described above for `f`. The returned value must be a 1-D array.
  912. Defaults to None, in which case the gradient is approximated
  913. numerically (see `epsilon`, below).
  914. args : tuple, optional
  915. Parameter values passed to `f` and `fprime`. Must be supplied whenever
  916. additional fixed parameters are needed to completely specify the
  917. functions `f` and `fprime`.
  918. gtol : float, optional
  919. Stop when the norm of the gradient is less than `gtol`.
  920. norm : float, optional
  921. Order to use for the norm of the gradient
  922. (``-np.Inf`` is min, ``np.Inf`` is max).
  923. epsilon : float or ndarray, optional
  924. Step size(s) to use when `fprime` is approximated numerically. Can be a
  925. scalar or a 1-D array. Defaults to ``sqrt(eps)``, with eps the
  926. floating point machine precision. Usually ``sqrt(eps)`` is about
  927. 1.5e-8.
  928. maxiter : int, optional
  929. Maximum number of iterations to perform. Default is ``200 * len(x0)``.
  930. full_output : bool, optional
  931. If True, return `fopt`, `func_calls`, `grad_calls`, and `warnflag` in
  932. addition to `xopt`. See the Returns section below for additional
  933. information on optional return values.
  934. disp : bool, optional
  935. If True, return a convergence message, followed by `xopt`.
  936. retall : bool, optional
  937. If True, add to the returned values the results of each iteration.
  938. callback : callable, optional
  939. An optional user-supplied function, called after each iteration.
  940. Called as ``callback(xk)``, where ``xk`` is the current value of `x0`.
  941. Returns
  942. -------
  943. xopt : ndarray
  944. Parameters which minimize f, i.e. ``f(xopt) == fopt``.
  945. fopt : float, optional
  946. Minimum value found, f(xopt). Only returned if `full_output` is True.
  947. func_calls : int, optional
  948. The number of function_calls made. Only returned if `full_output`
  949. is True.
  950. grad_calls : int, optional
  951. The number of gradient calls made. Only returned if `full_output` is
  952. True.
  953. warnflag : int, optional
  954. Integer value with warning status, only returned if `full_output` is
  955. True.
  956. 0 : Success.
  957. 1 : The maximum number of iterations was exceeded.
  958. 2 : Gradient and/or function calls were not changing. May indicate
  959. that precision was lost, i.e., the routine did not converge.
  960. allvecs : list of ndarray, optional
  961. List of arrays, containing the results at each iteration.
  962. Only returned if `retall` is True.
  963. See Also
  964. --------
  965. minimize : common interface to all `scipy.optimize` algorithms for
  966. unconstrained and constrained minimization of multivariate
  967. functions. It provides an alternative way to call
  968. ``fmin_cg``, by specifying ``method='CG'``.
  969. Notes
  970. -----
  971. This conjugate gradient algorithm is based on that of Polak and Ribiere
  972. [1]_.
  973. Conjugate gradient methods tend to work better when:
  974. 1. `f` has a unique global minimizing point, and no local minima or
  975. other stationary points,
  976. 2. `f` is, at least locally, reasonably well approximated by a
  977. quadratic function of the variables,
  978. 3. `f` is continuous and has a continuous gradient,
  979. 4. `fprime` is not too large, e.g., has a norm less than 1000,
  980. 5. The initial guess, `x0`, is reasonably close to `f` 's global
  981. minimizing point, `xopt`.
  982. References
  983. ----------
  984. .. [1] Wright & Nocedal, "Numerical Optimization", 1999, pp. 120-122.
  985. Examples
  986. --------
  987. Example 1: seek the minimum value of the expression
  988. ``a*u**2 + b*u*v + c*v**2 + d*u + e*v + f`` for given values
  989. of the parameters and an initial guess ``(u, v) = (0, 0)``.
  990. >>> args = (2, 3, 7, 8, 9, 10) # parameter values
  991. >>> def f(x, *args):
  992. ... u, v = x
  993. ... a, b, c, d, e, f = args
  994. ... return a*u**2 + b*u*v + c*v**2 + d*u + e*v + f
  995. >>> def gradf(x, *args):
  996. ... u, v = x
  997. ... a, b, c, d, e, f = args
  998. ... gu = 2*a*u + b*v + d # u-component of the gradient
  999. ... gv = b*u + 2*c*v + e # v-component of the gradient
  1000. ... return np.asarray((gu, gv))
  1001. >>> x0 = np.asarray((0, 0)) # Initial guess.
  1002. >>> from scipy import optimize
  1003. >>> res1 = optimize.fmin_cg(f, x0, fprime=gradf, args=args)
  1004. Optimization terminated successfully.
  1005. Current function value: 1.617021
  1006. Iterations: 4
  1007. Function evaluations: 8
  1008. Gradient evaluations: 8
  1009. >>> res1
  1010. array([-1.80851064, -0.25531915])
  1011. Example 2: solve the same problem using the `minimize` function.
  1012. (This `myopts` dictionary shows all of the available options,
  1013. although in practice only non-default values would be needed.
  1014. The returned value will be a dictionary.)
  1015. >>> opts = {'maxiter' : None, # default value.
  1016. ... 'disp' : True, # non-default value.
  1017. ... 'gtol' : 1e-5, # default value.
  1018. ... 'norm' : np.inf, # default value.
  1019. ... 'eps' : 1.4901161193847656e-08} # default value.
  1020. >>> res2 = optimize.minimize(f, x0, jac=gradf, args=args,
  1021. ... method='CG', options=opts)
  1022. Optimization terminated successfully.
  1023. Current function value: 1.617021
  1024. Iterations: 4
  1025. Function evaluations: 8
  1026. Gradient evaluations: 8
  1027. >>> res2.x # minimum found
  1028. array([-1.80851064, -0.25531915])
  1029. """
  1030. opts = {'gtol': gtol,
  1031. 'norm': norm,
  1032. 'eps': epsilon,
  1033. 'disp': disp,
  1034. 'maxiter': maxiter,
  1035. 'return_all': retall}
  1036. res = _minimize_cg(f, x0, args, fprime, callback=callback, **opts)
  1037. if full_output:
  1038. retlist = res['x'], res['fun'], res['nfev'], res['njev'], res['status']
  1039. if retall:
  1040. retlist += (res['allvecs'], )
  1041. return retlist
  1042. else:
  1043. if retall:
  1044. return res['x'], res['allvecs']
  1045. else:
  1046. return res['x']
  1047. def _minimize_cg(fun, x0, args=(), jac=None, callback=None,
  1048. gtol=1e-5, norm=Inf, eps=_epsilon, maxiter=None,
  1049. disp=False, return_all=False,
  1050. **unknown_options):
  1051. """
  1052. Minimization of scalar function of one or more variables using the
  1053. conjugate gradient algorithm.
  1054. Options
  1055. -------
  1056. disp : bool
  1057. Set to True to print convergence messages.
  1058. maxiter : int
  1059. Maximum number of iterations to perform.
  1060. gtol : float
  1061. Gradient norm must be less than `gtol` before successful
  1062. termination.
  1063. norm : float
  1064. Order of norm (Inf is max, -Inf is min).
  1065. eps : float or ndarray
  1066. If `jac` is approximated, use this value for the step size.
  1067. """
  1068. _check_unknown_options(unknown_options)
  1069. f = fun
  1070. fprime = jac
  1071. epsilon = eps
  1072. retall = return_all
  1073. x0 = asarray(x0).flatten()
  1074. if maxiter is None:
  1075. maxiter = len(x0) * 200
  1076. func_calls, f = wrap_function(f, args)
  1077. if fprime is None:
  1078. grad_calls, myfprime = wrap_function(approx_fprime, (f, epsilon))
  1079. else:
  1080. grad_calls, myfprime = wrap_function(fprime, args)
  1081. gfk = myfprime(x0)
  1082. k = 0
  1083. xk = x0
  1084. # Sets the initial step guess to dx ~ 1
  1085. old_fval = f(xk)
  1086. old_old_fval = old_fval + np.linalg.norm(gfk) / 2
  1087. if retall:
  1088. allvecs = [xk]
  1089. warnflag = 0
  1090. pk = -gfk
  1091. gnorm = vecnorm(gfk, ord=norm)
  1092. sigma_3 = 0.01
  1093. while (gnorm > gtol) and (k < maxiter):
  1094. deltak = numpy.dot(gfk, gfk)
  1095. cached_step = [None]
  1096. def polak_ribiere_powell_step(alpha, gfkp1=None):
  1097. xkp1 = xk + alpha * pk
  1098. if gfkp1 is None:
  1099. gfkp1 = myfprime(xkp1)
  1100. yk = gfkp1 - gfk
  1101. beta_k = max(0, numpy.dot(yk, gfkp1) / deltak)
  1102. pkp1 = -gfkp1 + beta_k * pk
  1103. gnorm = vecnorm(gfkp1, ord=norm)
  1104. return (alpha, xkp1, pkp1, gfkp1, gnorm)
  1105. def descent_condition(alpha, xkp1, fp1, gfkp1):
  1106. # Polak-Ribiere+ needs an explicit check of a sufficient
  1107. # descent condition, which is not guaranteed by strong Wolfe.
  1108. #
  1109. # See Gilbert & Nocedal, "Global convergence properties of
  1110. # conjugate gradient methods for optimization",
  1111. # SIAM J. Optimization 2, 21 (1992).
  1112. cached_step[:] = polak_ribiere_powell_step(alpha, gfkp1)
  1113. alpha, xk, pk, gfk, gnorm = cached_step
  1114. # Accept step if it leads to convergence.
  1115. if gnorm <= gtol:
  1116. return True
  1117. # Accept step if sufficient descent condition applies.
  1118. return numpy.dot(pk, gfk) <= -sigma_3 * numpy.dot(gfk, gfk)
  1119. try:
  1120. alpha_k, fc, gc, old_fval, old_old_fval, gfkp1 = \
  1121. _line_search_wolfe12(f, myfprime, xk, pk, gfk, old_fval,
  1122. old_old_fval, c2=0.4, amin=1e-100, amax=1e100,
  1123. extra_condition=descent_condition)
  1124. except _LineSearchError:
  1125. # Line search failed to find a better solution.
  1126. warnflag = 2
  1127. break
  1128. # Reuse already computed results if possible
  1129. if alpha_k == cached_step[0]:
  1130. alpha_k, xk, pk, gfk, gnorm = cached_step
  1131. else:
  1132. alpha_k, xk, pk, gfk, gnorm = polak_ribiere_powell_step(alpha_k, gfkp1)
  1133. if retall:
  1134. allvecs.append(xk)
  1135. if callback is not None:
  1136. callback(xk)
  1137. k += 1
  1138. fval = old_fval
  1139. if warnflag == 2:
  1140. msg = _status_message['pr_loss']
  1141. elif k >= maxiter:
  1142. warnflag = 1
  1143. msg = _status_message['maxiter']
  1144. else:
  1145. msg = _status_message['success']
  1146. if disp:
  1147. print("%s%s" % ("Warning: " if warnflag != 0 else "", msg))
  1148. print(" Current function value: %f" % fval)
  1149. print(" Iterations: %d" % k)
  1150. print(" Function evaluations: %d" % func_calls[0])
  1151. print(" Gradient evaluations: %d" % grad_calls[0])
  1152. result = OptimizeResult(fun=fval, jac=gfk, nfev=func_calls[0],
  1153. njev=grad_calls[0], status=warnflag,
  1154. success=(warnflag == 0), message=msg, x=xk,
  1155. nit=k)
  1156. if retall:
  1157. result['allvecs'] = allvecs
  1158. return result
  1159. def fmin_ncg(f, x0, fprime, fhess_p=None, fhess=None, args=(), avextol=1e-5,
  1160. epsilon=_epsilon, maxiter=None, full_output=0, disp=1, retall=0,
  1161. callback=None):
  1162. """
  1163. Unconstrained minimization of a function using the Newton-CG method.
  1164. Parameters
  1165. ----------
  1166. f : callable ``f(x, *args)``
  1167. Objective function to be minimized.
  1168. x0 : ndarray
  1169. Initial guess.
  1170. fprime : callable ``f'(x, *args)``
  1171. Gradient of f.
  1172. fhess_p : callable ``fhess_p(x, p, *args)``, optional
  1173. Function which computes the Hessian of f times an
  1174. arbitrary vector, p.
  1175. fhess : callable ``fhess(x, *args)``, optional
  1176. Function to compute the Hessian matrix of f.
  1177. args : tuple, optional
  1178. Extra arguments passed to f, fprime, fhess_p, and fhess
  1179. (the same set of extra arguments is supplied to all of
  1180. these functions).
  1181. epsilon : float or ndarray, optional
  1182. If fhess is approximated, use this value for the step size.
  1183. callback : callable, optional
  1184. An optional user-supplied function which is called after
  1185. each iteration. Called as callback(xk), where xk is the
  1186. current parameter vector.
  1187. avextol : float, optional
  1188. Convergence is assumed when the average relative error in
  1189. the minimizer falls below this amount.
  1190. maxiter : int, optional
  1191. Maximum number of iterations to perform.
  1192. full_output : bool, optional
  1193. If True, return the optional outputs.
  1194. disp : bool, optional
  1195. If True, print convergence message.
  1196. retall : bool, optional
  1197. If True, return a list of results at each iteration.
  1198. Returns
  1199. -------
  1200. xopt : ndarray
  1201. Parameters which minimize f, i.e. ``f(xopt) == fopt``.
  1202. fopt : float
  1203. Value of the function at xopt, i.e. ``fopt = f(xopt)``.
  1204. fcalls : int
  1205. Number of function calls made.
  1206. gcalls : int
  1207. Number of gradient calls made.
  1208. hcalls : int
  1209. Number of hessian calls made.
  1210. warnflag : int
  1211. Warnings generated by the algorithm.
  1212. 1 : Maximum number of iterations exceeded.
  1213. allvecs : list
  1214. The result at each iteration, if retall is True (see below).
  1215. See also
  1216. --------
  1217. minimize: Interface to minimization algorithms for multivariate
  1218. functions. See the 'Newton-CG' `method` in particular.
  1219. Notes
  1220. -----
  1221. Only one of `fhess_p` or `fhess` need to be given. If `fhess`
  1222. is provided, then `fhess_p` will be ignored. If neither `fhess`
  1223. nor `fhess_p` is provided, then the hessian product will be
  1224. approximated using finite differences on `fprime`. `fhess_p`
  1225. must compute the hessian times an arbitrary vector. If it is not
  1226. given, finite-differences on `fprime` are used to compute
  1227. it.
  1228. Newton-CG methods are also called truncated Newton methods. This
  1229. function differs from scipy.optimize.fmin_tnc because
  1230. 1. scipy.optimize.fmin_ncg is written purely in python using numpy
  1231. and scipy while scipy.optimize.fmin_tnc calls a C function.
  1232. 2. scipy.optimize.fmin_ncg is only for unconstrained minimization
  1233. while scipy.optimize.fmin_tnc is for unconstrained minimization
  1234. or box constrained minimization. (Box constraints give
  1235. lower and upper bounds for each variable separately.)
  1236. References
  1237. ----------
  1238. Wright & Nocedal, 'Numerical Optimization', 1999, pg. 140.
  1239. """
  1240. opts = {'xtol': avextol,
  1241. 'eps': epsilon,
  1242. 'maxiter': maxiter,
  1243. 'disp': disp,
  1244. 'return_all': retall}
  1245. res = _minimize_newtoncg(f, x0, args, fprime, fhess, fhess_p,
  1246. callback=callback, **opts)
  1247. if full_output:
  1248. retlist = (res['x'], res['fun'], res['nfev'], res['njev'],
  1249. res['nhev'], res['status'])
  1250. if retall:
  1251. retlist += (res['allvecs'], )
  1252. return retlist
  1253. else:
  1254. if retall:
  1255. return res['x'], res['allvecs']
  1256. else:
  1257. return res['x']
  1258. def _minimize_newtoncg(fun, x0, args=(), jac=None, hess=None, hessp=None,
  1259. callback=None, xtol=1e-5, eps=_epsilon, maxiter=None,
  1260. disp=False, return_all=False,
  1261. **unknown_options):
  1262. """
  1263. Minimization of scalar function of one or more variables using the
  1264. Newton-CG algorithm.
  1265. Note that the `jac` parameter (Jacobian) is required.
  1266. Options
  1267. -------
  1268. disp : bool
  1269. Set to True to print convergence messages.
  1270. xtol : float
  1271. Average relative error in solution `xopt` acceptable for
  1272. convergence.
  1273. maxiter : int
  1274. Maximum number of iterations to perform.
  1275. eps : float or ndarray
  1276. If `jac` is approximated, use this value for the step size.
  1277. """
  1278. _check_unknown_options(unknown_options)
  1279. if jac is None:
  1280. raise ValueError('Jacobian is required for Newton-CG method')
  1281. f = fun
  1282. fprime = jac
  1283. fhess_p = hessp
  1284. fhess = hess
  1285. avextol = xtol
  1286. epsilon = eps
  1287. retall = return_all
  1288. def terminate(warnflag, msg):
  1289. if disp:
  1290. print(msg)
  1291. print(" Current function value: %f" % old_fval)
  1292. print(" Iterations: %d" % k)
  1293. print(" Function evaluations: %d" % fcalls[0])
  1294. print(" Gradient evaluations: %d" % gcalls[0])
  1295. print(" Hessian evaluations: %d" % hcalls)
  1296. fval = old_fval
  1297. result = OptimizeResult(fun=fval, jac=gfk, nfev=fcalls[0],
  1298. njev=gcalls[0], nhev=hcalls, status=warnflag,
  1299. success=(warnflag == 0), message=msg, x=xk,
  1300. nit=k)
  1301. if retall:
  1302. result['allvecs'] = allvecs
  1303. return result
  1304. x0 = asarray(x0).flatten()
  1305. fcalls, f = wrap_function(f, args)
  1306. gcalls, fprime = wrap_function(fprime, args)
  1307. hcalls = 0
  1308. if maxiter is None:
  1309. maxiter = len(x0)*200
  1310. cg_maxiter = 20*len(x0)
  1311. xtol = len(x0) * avextol
  1312. update = [2 * xtol]
  1313. xk = x0
  1314. if retall:
  1315. allvecs = [xk]
  1316. k = 0
  1317. gfk = None
  1318. old_fval = f(x0)
  1319. old_old_fval = None
  1320. float64eps = numpy.finfo(numpy.float64).eps
  1321. while numpy.add.reduce(numpy.abs(update)) > xtol:
  1322. if k >= maxiter:
  1323. msg = "Warning: " + _status_message['maxiter']
  1324. return terminate(1, msg)
  1325. # Compute a search direction pk by applying the CG method to
  1326. # del2 f(xk) p = - grad f(xk) starting from 0.
  1327. b = -fprime(xk)
  1328. maggrad = numpy.add.reduce(numpy.abs(b))
  1329. eta = numpy.min([0.5, numpy.sqrt(maggrad)])
  1330. termcond = eta * maggrad
  1331. xsupi = zeros(len(x0), dtype=x0.dtype)
  1332. ri = -b
  1333. psupi = -ri
  1334. i = 0
  1335. dri0 = numpy.dot(ri, ri)
  1336. if fhess is not None: # you want to compute hessian once.
  1337. A = fhess(*(xk,) + args)
  1338. hcalls = hcalls + 1
  1339. for k2 in xrange(cg_maxiter):
  1340. if numpy.add.reduce(numpy.abs(ri)) <= termcond:
  1341. break
  1342. if fhess is None:
  1343. if fhess_p is None:
  1344. Ap = approx_fhess_p(xk, psupi, fprime, epsilon)
  1345. else:
  1346. Ap = fhess_p(xk, psupi, *args)
  1347. hcalls = hcalls + 1
  1348. else:
  1349. Ap = numpy.dot(A, psupi)
  1350. # check curvature
  1351. Ap = asarray(Ap).squeeze() # get rid of matrices...
  1352. curv = numpy.dot(psupi, Ap)
  1353. if 0 <= curv <= 3 * float64eps:
  1354. break
  1355. elif curv < 0:
  1356. if (i > 0):
  1357. break
  1358. else:
  1359. # fall back to steepest descent direction
  1360. xsupi = dri0 / (-curv) * b
  1361. break
  1362. alphai = dri0 / curv
  1363. xsupi = xsupi + alphai * psupi
  1364. ri = ri + alphai * Ap
  1365. dri1 = numpy.dot(ri, ri)
  1366. betai = dri1 / dri0
  1367. psupi = -ri + betai * psupi
  1368. i = i + 1
  1369. dri0 = dri1 # update numpy.dot(ri,ri) for next time.
  1370. else:
  1371. # curvature keeps increasing, bail out
  1372. msg = ("Warning: CG iterations didn't converge. The Hessian is not "
  1373. "positive definite.")
  1374. return terminate(3, msg)
  1375. pk = xsupi # search direction is solution to system.
  1376. gfk = -b # gradient at xk
  1377. try:
  1378. alphak, fc, gc, old_fval, old_old_fval, gfkp1 = \
  1379. _line_search_wolfe12(f, fprime, xk, pk, gfk,
  1380. old_fval, old_old_fval)
  1381. except _LineSearchError:
  1382. # Line search failed to find a better solution.
  1383. msg = "Warning: " + _status_message['pr_loss']
  1384. return terminate(2, msg)
  1385. update = alphak * pk
  1386. xk = xk + update # upcast if necessary
  1387. if callback is not None:
  1388. callback(xk)
  1389. if retall:
  1390. allvecs.append(xk)
  1391. k += 1
  1392. else:
  1393. msg = _status_message['success']
  1394. return terminate(0, msg)
  1395. def fminbound(func, x1, x2, args=(), xtol=1e-5, maxfun=500,
  1396. full_output=0, disp=1):
  1397. """Bounded minimization for scalar functions.
  1398. Parameters
  1399. ----------
  1400. func : callable f(x,*args)
  1401. Objective function to be minimized (must accept and return scalars).
  1402. x1, x2 : float or array scalar
  1403. The optimization bounds.
  1404. args : tuple, optional
  1405. Extra arguments passed to function.
  1406. xtol : float, optional
  1407. The convergence tolerance.
  1408. maxfun : int, optional
  1409. Maximum number of function evaluations allowed.
  1410. full_output : bool, optional
  1411. If True, return optional outputs.
  1412. disp : int, optional
  1413. If non-zero, print messages.
  1414. 0 : no message printing.
  1415. 1 : non-convergence notification messages only.
  1416. 2 : print a message on convergence too.
  1417. 3 : print iteration results.
  1418. Returns
  1419. -------
  1420. xopt : ndarray
  1421. Parameters (over given interval) which minimize the
  1422. objective function.
  1423. fval : number
  1424. The function value at the minimum point.
  1425. ierr : int
  1426. An error flag (0 if converged, 1 if maximum number of
  1427. function calls reached).
  1428. numfunc : int
  1429. The number of function calls made.
  1430. See also
  1431. --------
  1432. minimize_scalar: Interface to minimization algorithms for scalar
  1433. univariate functions. See the 'Bounded' `method` in particular.
  1434. Notes
  1435. -----
  1436. Finds a local minimizer of the scalar function `func` in the
  1437. interval x1 < xopt < x2 using Brent's method. (See `brent`
  1438. for auto-bracketing).
  1439. Examples
  1440. --------
  1441. `fminbound` finds the minimum of the function in the given range.
  1442. The following examples illustrate the same
  1443. >>> def f(x):
  1444. ... return x**2
  1445. >>> from scipy import optimize
  1446. >>> minimum = optimize.fminbound(f, -1, 2)
  1447. >>> minimum
  1448. 0.0
  1449. >>> minimum = optimize.fminbound(f, 1, 2)
  1450. >>> minimum
  1451. 1.0000059608609866
  1452. """
  1453. options = {'xatol': xtol,
  1454. 'maxiter': maxfun,
  1455. 'disp': disp}
  1456. res = _minimize_scalar_bounded(func, (x1, x2), args, **options)
  1457. if full_output:
  1458. return res['x'], res['fun'], res['status'], res['nfev']
  1459. else:
  1460. return res['x']
  1461. def _minimize_scalar_bounded(func, bounds, args=(),
  1462. xatol=1e-5, maxiter=500, disp=0,
  1463. **unknown_options):
  1464. """
  1465. Options
  1466. -------
  1467. maxiter : int
  1468. Maximum number of iterations to perform.
  1469. disp: int, optional
  1470. If non-zero, print messages.
  1471. 0 : no message printing.
  1472. 1 : non-convergence notification messages only.
  1473. 2 : print a message on convergence too.
  1474. 3 : print iteration results.
  1475. xatol : float
  1476. Absolute error in solution `xopt` acceptable for convergence.
  1477. """
  1478. _check_unknown_options(unknown_options)
  1479. maxfun = maxiter
  1480. # Test bounds are of correct form
  1481. if len(bounds) != 2:
  1482. raise ValueError('bounds must have two elements.')
  1483. x1, x2 = bounds
  1484. if not (is_array_scalar(x1) and is_array_scalar(x2)):
  1485. raise ValueError("Optimisation bounds must be scalars"
  1486. " or array scalars.")
  1487. if x1 > x2:
  1488. raise ValueError("The lower bound exceeds the upper bound.")
  1489. flag = 0
  1490. header = ' Func-count x f(x) Procedure'
  1491. step = ' initial'
  1492. sqrt_eps = sqrt(2.2e-16)
  1493. golden_mean = 0.5 * (3.0 - sqrt(5.0))
  1494. a, b = x1, x2
  1495. fulc = a + golden_mean * (b - a)
  1496. nfc, xf = fulc, fulc
  1497. rat = e = 0.0
  1498. x = xf
  1499. fx = func(x, *args)
  1500. num = 1
  1501. fmin_data = (1, xf, fx)
  1502. ffulc = fnfc = fx
  1503. xm = 0.5 * (a + b)
  1504. tol1 = sqrt_eps * numpy.abs(xf) + xatol / 3.0
  1505. tol2 = 2.0 * tol1
  1506. if disp > 2:
  1507. print(" ")
  1508. print(header)
  1509. print("%5.0f %12.6g %12.6g %s" % (fmin_data + (step,)))
  1510. while (numpy.abs(xf - xm) > (tol2 - 0.5 * (b - a))):
  1511. golden = 1
  1512. # Check for parabolic fit
  1513. if numpy.abs(e) > tol1:
  1514. golden = 0
  1515. r = (xf - nfc) * (fx - ffulc)
  1516. q = (xf - fulc) * (fx - fnfc)
  1517. p = (xf - fulc) * q - (xf - nfc) * r
  1518. q = 2.0 * (q - r)
  1519. if q > 0.0:
  1520. p = -p
  1521. q = numpy.abs(q)
  1522. r = e
  1523. e = rat
  1524. # Check for acceptability of parabola
  1525. if ((numpy.abs(p) < numpy.abs(0.5*q*r)) and (p > q*(a - xf)) and
  1526. (p < q * (b - xf))):
  1527. rat = (p + 0.0) / q
  1528. x = xf + rat
  1529. step = ' parabolic'
  1530. if ((x - a) < tol2) or ((b - x) < tol2):
  1531. si = numpy.sign(xm - xf) + ((xm - xf) == 0)
  1532. rat = tol1 * si
  1533. else: # do a golden section step
  1534. golden = 1
  1535. if golden: # Do a golden-section step
  1536. if xf >= xm:
  1537. e = a - xf
  1538. else:
  1539. e = b - xf
  1540. rat = golden_mean*e
  1541. step = ' golden'
  1542. si = numpy.sign(rat) + (rat == 0)
  1543. x = xf + si * numpy.max([numpy.abs(rat), tol1])
  1544. fu = func(x, *args)
  1545. num += 1
  1546. fmin_data = (num, x, fu)
  1547. if disp > 2:
  1548. print("%5.0f %12.6g %12.6g %s" % (fmin_data + (step,)))
  1549. if fu <= fx:
  1550. if x >= xf:
  1551. a = xf
  1552. else:
  1553. b = xf
  1554. fulc, ffulc = nfc, fnfc
  1555. nfc, fnfc = xf, fx
  1556. xf, fx = x, fu
  1557. else:
  1558. if x < xf:
  1559. a = x
  1560. else:
  1561. b = x
  1562. if (fu <= fnfc) or (nfc == xf):
  1563. fulc, ffulc = nfc, fnfc
  1564. nfc, fnfc = x, fu
  1565. elif (fu <= ffulc) or (fulc == xf) or (fulc == nfc):
  1566. fulc, ffulc = x, fu
  1567. xm = 0.5 * (a + b)
  1568. tol1 = sqrt_eps * numpy.abs(xf) + xatol / 3.0
  1569. tol2 = 2.0 * tol1
  1570. if num >= maxfun:
  1571. flag = 1
  1572. break
  1573. fval = fx
  1574. if disp > 0:
  1575. _endprint(x, flag, fval, maxfun, xatol, disp)
  1576. result = OptimizeResult(fun=fval, status=flag, success=(flag == 0),
  1577. message={0: 'Solution found.',
  1578. 1: 'Maximum number of function calls '
  1579. 'reached.'}.get(flag, ''),
  1580. x=xf, nfev=num)
  1581. return result
  1582. class Brent:
  1583. #need to rethink design of __init__
  1584. def __init__(self, func, args=(), tol=1.48e-8, maxiter=500,
  1585. full_output=0):
  1586. self.func = func
  1587. self.args = args
  1588. self.tol = tol
  1589. self.maxiter = maxiter
  1590. self._mintol = 1.0e-11
  1591. self._cg = 0.3819660
  1592. self.xmin = None
  1593. self.fval = None
  1594. self.iter = 0
  1595. self.funcalls = 0
  1596. # need to rethink design of set_bracket (new options, etc)
  1597. def set_bracket(self, brack=None):
  1598. self.brack = brack
  1599. def get_bracket_info(self):
  1600. #set up
  1601. func = self.func
  1602. args = self.args
  1603. brack = self.brack
  1604. ### BEGIN core bracket_info code ###
  1605. ### carefully DOCUMENT any CHANGES in core ##
  1606. if brack is None:
  1607. xa, xb, xc, fa, fb, fc, funcalls = bracket(func, args=args)
  1608. elif len(brack) == 2:
  1609. xa, xb, xc, fa, fb, fc, funcalls = bracket(func, xa=brack[0],
  1610. xb=brack[1], args=args)
  1611. elif len(brack) == 3:
  1612. xa, xb, xc = brack
  1613. if (xa > xc): # swap so xa < xc can be assumed
  1614. xc, xa = xa, xc
  1615. if not ((xa < xb) and (xb < xc)):
  1616. raise ValueError("Not a bracketing interval.")
  1617. fa = func(*((xa,) + args))
  1618. fb = func(*((xb,) + args))
  1619. fc = func(*((xc,) + args))
  1620. if not ((fb < fa) and (fb < fc)):
  1621. raise ValueError("Not a bracketing interval.")
  1622. funcalls = 3
  1623. else:
  1624. raise ValueError("Bracketing interval must be "
  1625. "length 2 or 3 sequence.")
  1626. ### END core bracket_info code ###
  1627. return xa, xb, xc, fa, fb, fc, funcalls
  1628. def optimize(self):
  1629. # set up for optimization
  1630. func = self.func
  1631. xa, xb, xc, fa, fb, fc, funcalls = self.get_bracket_info()
  1632. _mintol = self._mintol
  1633. _cg = self._cg
  1634. #################################
  1635. #BEGIN CORE ALGORITHM
  1636. #################################
  1637. x = w = v = xb
  1638. fw = fv = fx = func(*((x,) + self.args))
  1639. if (xa < xc):
  1640. a = xa
  1641. b = xc
  1642. else:
  1643. a = xc
  1644. b = xa
  1645. deltax = 0.0
  1646. funcalls += 1
  1647. iter = 0
  1648. while (iter < self.maxiter):
  1649. tol1 = self.tol * numpy.abs(x) + _mintol
  1650. tol2 = 2.0 * tol1
  1651. xmid = 0.5 * (a + b)
  1652. # check for convergence
  1653. if numpy.abs(x - xmid) < (tol2 - 0.5 * (b - a)):
  1654. break
  1655. # XXX In the first iteration, rat is only bound in the true case
  1656. # of this conditional. This used to cause an UnboundLocalError
  1657. # (gh-4140). It should be set before the if (but to what?).
  1658. if (numpy.abs(deltax) <= tol1):
  1659. if (x >= xmid):
  1660. deltax = a - x # do a golden section step
  1661. else:
  1662. deltax = b - x
  1663. rat = _cg * deltax
  1664. else: # do a parabolic step
  1665. tmp1 = (x - w) * (fx - fv)
  1666. tmp2 = (x - v) * (fx - fw)
  1667. p = (x - v) * tmp2 - (x - w) * tmp1
  1668. tmp2 = 2.0 * (tmp2 - tmp1)
  1669. if (tmp2 > 0.0):
  1670. p = -p
  1671. tmp2 = numpy.abs(tmp2)
  1672. dx_temp = deltax
  1673. deltax = rat
  1674. # check parabolic fit
  1675. if ((p > tmp2 * (a - x)) and (p < tmp2 * (b - x)) and
  1676. (numpy.abs(p) < numpy.abs(0.5 * tmp2 * dx_temp))):
  1677. rat = p * 1.0 / tmp2 # if parabolic step is useful.
  1678. u = x + rat
  1679. if ((u - a) < tol2 or (b - u) < tol2):
  1680. if xmid - x >= 0:
  1681. rat = tol1
  1682. else:
  1683. rat = -tol1
  1684. else:
  1685. if (x >= xmid):
  1686. deltax = a - x # if it's not do a golden section step
  1687. else:
  1688. deltax = b - x
  1689. rat = _cg * deltax
  1690. if (numpy.abs(rat) < tol1): # update by at least tol1
  1691. if rat >= 0:
  1692. u = x + tol1
  1693. else:
  1694. u = x - tol1
  1695. else:
  1696. u = x + rat
  1697. fu = func(*((u,) + self.args)) # calculate new output value
  1698. funcalls += 1
  1699. if (fu > fx): # if it's bigger than current
  1700. if (u < x):
  1701. a = u
  1702. else:
  1703. b = u
  1704. if (fu <= fw) or (w == x):
  1705. v = w
  1706. w = u
  1707. fv = fw
  1708. fw = fu
  1709. elif (fu <= fv) or (v == x) or (v == w):
  1710. v = u
  1711. fv = fu
  1712. else:
  1713. if (u >= x):
  1714. a = x
  1715. else:
  1716. b = x
  1717. v = w
  1718. w = x
  1719. x = u
  1720. fv = fw
  1721. fw = fx
  1722. fx = fu
  1723. iter += 1
  1724. #################################
  1725. #END CORE ALGORITHM
  1726. #################################
  1727. self.xmin = x
  1728. self.fval = fx
  1729. self.iter = iter
  1730. self.funcalls = funcalls
  1731. def get_result(self, full_output=False):
  1732. if full_output:
  1733. return self.xmin, self.fval, self.iter, self.funcalls
  1734. else:
  1735. return self.xmin
  1736. def brent(func, args=(), brack=None, tol=1.48e-8, full_output=0, maxiter=500):
  1737. """
  1738. Given a function of one-variable and a possible bracket, return
  1739. the local minimum of the function isolated to a fractional precision
  1740. of tol.
  1741. Parameters
  1742. ----------
  1743. func : callable f(x,*args)
  1744. Objective function.
  1745. args : tuple, optional
  1746. Additional arguments (if present).
  1747. brack : tuple, optional
  1748. Either a triple (xa,xb,xc) where xa<xb<xc and func(xb) <
  1749. func(xa), func(xc) or a pair (xa,xb) which are used as a
  1750. starting interval for a downhill bracket search (see
  1751. `bracket`). Providing the pair (xa,xb) does not always mean
  1752. the obtained solution will satisfy xa<=x<=xb.
  1753. tol : float, optional
  1754. Stop if between iteration change is less than `tol`.
  1755. full_output : bool, optional
  1756. If True, return all output args (xmin, fval, iter,
  1757. funcalls).
  1758. maxiter : int, optional
  1759. Maximum number of iterations in solution.
  1760. Returns
  1761. -------
  1762. xmin : ndarray
  1763. Optimum point.
  1764. fval : float
  1765. Optimum value.
  1766. iter : int
  1767. Number of iterations.
  1768. funcalls : int
  1769. Number of objective function evaluations made.
  1770. See also
  1771. --------
  1772. minimize_scalar: Interface to minimization algorithms for scalar
  1773. univariate functions. See the 'Brent' `method` in particular.
  1774. Notes
  1775. -----
  1776. Uses inverse parabolic interpolation when possible to speed up
  1777. convergence of golden section method.
  1778. Does not ensure that the minimum lies in the range specified by
  1779. `brack`. See `fminbound`.
  1780. Examples
  1781. --------
  1782. We illustrate the behaviour of the function when `brack` is of
  1783. size 2 and 3 respectively. In the case where `brack` is of the
  1784. form (xa,xb), we can see for the given values, the output need
  1785. not necessarily lie in the range (xa,xb).
  1786. >>> def f(x):
  1787. ... return x**2
  1788. >>> from scipy import optimize
  1789. >>> minimum = optimize.brent(f,brack=(1,2))
  1790. >>> minimum
  1791. 0.0
  1792. >>> minimum = optimize.brent(f,brack=(-1,0.5,2))
  1793. >>> minimum
  1794. -2.7755575615628914e-17
  1795. """
  1796. options = {'xtol': tol,
  1797. 'maxiter': maxiter}
  1798. res = _minimize_scalar_brent(func, brack, args, **options)
  1799. if full_output:
  1800. return res['x'], res['fun'], res['nit'], res['nfev']
  1801. else:
  1802. return res['x']
  1803. def _minimize_scalar_brent(func, brack=None, args=(),
  1804. xtol=1.48e-8, maxiter=500,
  1805. **unknown_options):
  1806. """
  1807. Options
  1808. -------
  1809. maxiter : int
  1810. Maximum number of iterations to perform.
  1811. xtol : float
  1812. Relative error in solution `xopt` acceptable for convergence.
  1813. Notes
  1814. -----
  1815. Uses inverse parabolic interpolation when possible to speed up
  1816. convergence of golden section method.
  1817. """
  1818. _check_unknown_options(unknown_options)
  1819. tol = xtol
  1820. if tol < 0:
  1821. raise ValueError('tolerance should be >= 0, got %r' % tol)
  1822. brent = Brent(func=func, args=args, tol=tol,
  1823. full_output=True, maxiter=maxiter)
  1824. brent.set_bracket(brack)
  1825. brent.optimize()
  1826. x, fval, nit, nfev = brent.get_result(full_output=True)
  1827. return OptimizeResult(fun=fval, x=x, nit=nit, nfev=nfev,
  1828. success=nit < maxiter)
  1829. def golden(func, args=(), brack=None, tol=_epsilon,
  1830. full_output=0, maxiter=5000):
  1831. """
  1832. Return the minimum of a function of one variable using golden section
  1833. method.
  1834. Given a function of one variable and a possible bracketing interval,
  1835. return the minimum of the function isolated to a fractional precision of
  1836. tol.
  1837. Parameters
  1838. ----------
  1839. func : callable func(x,*args)
  1840. Objective function to minimize.
  1841. args : tuple, optional
  1842. Additional arguments (if present), passed to func.
  1843. brack : tuple, optional
  1844. Triple (a,b,c), where (a<b<c) and func(b) <
  1845. func(a),func(c). If bracket consists of two numbers (a,
  1846. c), then they are assumed to be a starting interval for a
  1847. downhill bracket search (see `bracket`); it doesn't always
  1848. mean that obtained solution will satisfy a<=x<=c.
  1849. tol : float, optional
  1850. x tolerance stop criterion
  1851. full_output : bool, optional
  1852. If True, return optional outputs.
  1853. maxiter : int
  1854. Maximum number of iterations to perform.
  1855. See also
  1856. --------
  1857. minimize_scalar: Interface to minimization algorithms for scalar
  1858. univariate functions. See the 'Golden' `method` in particular.
  1859. Notes
  1860. -----
  1861. Uses analog of bisection method to decrease the bracketed
  1862. interval.
  1863. Examples
  1864. --------
  1865. We illustrate the behaviour of the function when `brack` is of
  1866. size 2 and 3 respectively. In the case where `brack` is of the
  1867. form (xa,xb), we can see for the given values, the output need
  1868. not necessarily lie in the range ``(xa, xb)``.
  1869. >>> def f(x):
  1870. ... return x**2
  1871. >>> from scipy import optimize
  1872. >>> minimum = optimize.golden(f, brack=(1, 2))
  1873. >>> minimum
  1874. 1.5717277788484873e-162
  1875. >>> minimum = optimize.golden(f, brack=(-1, 0.5, 2))
  1876. >>> minimum
  1877. -1.5717277788484873e-162
  1878. """
  1879. options = {'xtol': tol, 'maxiter': maxiter}
  1880. res = _minimize_scalar_golden(func, brack, args, **options)
  1881. if full_output:
  1882. return res['x'], res['fun'], res['nfev']
  1883. else:
  1884. return res['x']
  1885. def _minimize_scalar_golden(func, brack=None, args=(),
  1886. xtol=_epsilon, maxiter=5000, **unknown_options):
  1887. """
  1888. Options
  1889. -------
  1890. maxiter : int
  1891. Maximum number of iterations to perform.
  1892. xtol : float
  1893. Relative error in solution `xopt` acceptable for convergence.
  1894. """
  1895. _check_unknown_options(unknown_options)
  1896. tol = xtol
  1897. if brack is None:
  1898. xa, xb, xc, fa, fb, fc, funcalls = bracket(func, args=args)
  1899. elif len(brack) == 2:
  1900. xa, xb, xc, fa, fb, fc, funcalls = bracket(func, xa=brack[0],
  1901. xb=brack[1], args=args)
  1902. elif len(brack) == 3:
  1903. xa, xb, xc = brack
  1904. if (xa > xc): # swap so xa < xc can be assumed
  1905. xc, xa = xa, xc
  1906. if not ((xa < xb) and (xb < xc)):
  1907. raise ValueError("Not a bracketing interval.")
  1908. fa = func(*((xa,) + args))
  1909. fb = func(*((xb,) + args))
  1910. fc = func(*((xc,) + args))
  1911. if not ((fb < fa) and (fb < fc)):
  1912. raise ValueError("Not a bracketing interval.")
  1913. funcalls = 3
  1914. else:
  1915. raise ValueError("Bracketing interval must be length 2 or 3 sequence.")
  1916. _gR = 0.61803399 # golden ratio conjugate: 2.0/(1.0+sqrt(5.0))
  1917. _gC = 1.0 - _gR
  1918. x3 = xc
  1919. x0 = xa
  1920. if (numpy.abs(xc - xb) > numpy.abs(xb - xa)):
  1921. x1 = xb
  1922. x2 = xb + _gC * (xc - xb)
  1923. else:
  1924. x2 = xb
  1925. x1 = xb - _gC * (xb - xa)
  1926. f1 = func(*((x1,) + args))
  1927. f2 = func(*((x2,) + args))
  1928. funcalls += 2
  1929. nit = 0
  1930. for i in xrange(maxiter):
  1931. if numpy.abs(x3 - x0) <= tol * (numpy.abs(x1) + numpy.abs(x2)):
  1932. break
  1933. if (f2 < f1):
  1934. x0 = x1
  1935. x1 = x2
  1936. x2 = _gR * x1 + _gC * x3
  1937. f1 = f2
  1938. f2 = func(*((x2,) + args))
  1939. else:
  1940. x3 = x2
  1941. x2 = x1
  1942. x1 = _gR * x2 + _gC * x0
  1943. f2 = f1
  1944. f1 = func(*((x1,) + args))
  1945. funcalls += 1
  1946. nit += 1
  1947. if (f1 < f2):
  1948. xmin = x1
  1949. fval = f1
  1950. else:
  1951. xmin = x2
  1952. fval = f2
  1953. return OptimizeResult(fun=fval, nfev=funcalls, x=xmin, nit=nit,
  1954. success=nit < maxiter)
  1955. def bracket(func, xa=0.0, xb=1.0, args=(), grow_limit=110.0, maxiter=1000):
  1956. """
  1957. Bracket the minimum of the function.
  1958. Given a function and distinct initial points, search in the
  1959. downhill direction (as defined by the initital points) and return
  1960. new points xa, xb, xc that bracket the minimum of the function
  1961. f(xa) > f(xb) < f(xc). It doesn't always mean that obtained
  1962. solution will satisfy xa<=x<=xb
  1963. Parameters
  1964. ----------
  1965. func : callable f(x,*args)
  1966. Objective function to minimize.
  1967. xa, xb : float, optional
  1968. Bracketing interval. Defaults `xa` to 0.0, and `xb` to 1.0.
  1969. args : tuple, optional
  1970. Additional arguments (if present), passed to `func`.
  1971. grow_limit : float, optional
  1972. Maximum grow limit. Defaults to 110.0
  1973. maxiter : int, optional
  1974. Maximum number of iterations to perform. Defaults to 1000.
  1975. Returns
  1976. -------
  1977. xa, xb, xc : float
  1978. Bracket.
  1979. fa, fb, fc : float
  1980. Objective function values in bracket.
  1981. funcalls : int
  1982. Number of function evaluations made.
  1983. """
  1984. _gold = 1.618034 # golden ratio: (1.0+sqrt(5.0))/2.0
  1985. _verysmall_num = 1e-21
  1986. fa = func(*(xa,) + args)
  1987. fb = func(*(xb,) + args)
  1988. if (fa < fb): # Switch so fa > fb
  1989. xa, xb = xb, xa
  1990. fa, fb = fb, fa
  1991. xc = xb + _gold * (xb - xa)
  1992. fc = func(*((xc,) + args))
  1993. funcalls = 3
  1994. iter = 0
  1995. while (fc < fb):
  1996. tmp1 = (xb - xa) * (fb - fc)
  1997. tmp2 = (xb - xc) * (fb - fa)
  1998. val = tmp2 - tmp1
  1999. if numpy.abs(val) < _verysmall_num:
  2000. denom = 2.0 * _verysmall_num
  2001. else:
  2002. denom = 2.0 * val
  2003. w = xb - ((xb - xc) * tmp2 - (xb - xa) * tmp1) / denom
  2004. wlim = xb + grow_limit * (xc - xb)
  2005. if iter > maxiter:
  2006. raise RuntimeError("Too many iterations.")
  2007. iter += 1
  2008. if (w - xc) * (xb - w) > 0.0:
  2009. fw = func(*((w,) + args))
  2010. funcalls += 1
  2011. if (fw < fc):
  2012. xa = xb
  2013. xb = w
  2014. fa = fb
  2015. fb = fw
  2016. return xa, xb, xc, fa, fb, fc, funcalls
  2017. elif (fw > fb):
  2018. xc = w
  2019. fc = fw
  2020. return xa, xb, xc, fa, fb, fc, funcalls
  2021. w = xc + _gold * (xc - xb)
  2022. fw = func(*((w,) + args))
  2023. funcalls += 1
  2024. elif (w - wlim)*(wlim - xc) >= 0.0:
  2025. w = wlim
  2026. fw = func(*((w,) + args))
  2027. funcalls += 1
  2028. elif (w - wlim)*(xc - w) > 0.0:
  2029. fw = func(*((w,) + args))
  2030. funcalls += 1
  2031. if (fw < fc):
  2032. xb = xc
  2033. xc = w
  2034. w = xc + _gold * (xc - xb)
  2035. fb = fc
  2036. fc = fw
  2037. fw = func(*((w,) + args))
  2038. funcalls += 1
  2039. else:
  2040. w = xc + _gold * (xc - xb)
  2041. fw = func(*((w,) + args))
  2042. funcalls += 1
  2043. xa = xb
  2044. xb = xc
  2045. xc = w
  2046. fa = fb
  2047. fb = fc
  2048. fc = fw
  2049. return xa, xb, xc, fa, fb, fc, funcalls
  2050. def _linesearch_powell(func, p, xi, tol=1e-3):
  2051. """Line-search algorithm using fminbound.
  2052. Find the minimium of the function ``func(x0+ alpha*direc)``.
  2053. """
  2054. def myfunc(alpha):
  2055. return func(p + alpha*xi)
  2056. alpha_min, fret, iter, num = brent(myfunc, full_output=1, tol=tol)
  2057. xi = alpha_min*xi
  2058. return squeeze(fret), p + xi, xi
  2059. def fmin_powell(func, x0, args=(), xtol=1e-4, ftol=1e-4, maxiter=None,
  2060. maxfun=None, full_output=0, disp=1, retall=0, callback=None,
  2061. direc=None):
  2062. """
  2063. Minimize a function using modified Powell's method. This method
  2064. only uses function values, not derivatives.
  2065. Parameters
  2066. ----------
  2067. func : callable f(x,*args)
  2068. Objective function to be minimized.
  2069. x0 : ndarray
  2070. Initial guess.
  2071. args : tuple, optional
  2072. Extra arguments passed to func.
  2073. callback : callable, optional
  2074. An optional user-supplied function, called after each
  2075. iteration. Called as ``callback(xk)``, where ``xk`` is the
  2076. current parameter vector.
  2077. direc : ndarray, optional
  2078. Initial direction set.
  2079. xtol : float, optional
  2080. Line-search error tolerance.
  2081. ftol : float, optional
  2082. Relative error in ``func(xopt)`` acceptable for convergence.
  2083. maxiter : int, optional
  2084. Maximum number of iterations to perform.
  2085. maxfun : int, optional
  2086. Maximum number of function evaluations to make.
  2087. full_output : bool, optional
  2088. If True, fopt, xi, direc, iter, funcalls, and
  2089. warnflag are returned.
  2090. disp : bool, optional
  2091. If True, print convergence messages.
  2092. retall : bool, optional
  2093. If True, return a list of the solution at each iteration.
  2094. Returns
  2095. -------
  2096. xopt : ndarray
  2097. Parameter which minimizes `func`.
  2098. fopt : number
  2099. Value of function at minimum: ``fopt = func(xopt)``.
  2100. direc : ndarray
  2101. Current direction set.
  2102. iter : int
  2103. Number of iterations.
  2104. funcalls : int
  2105. Number of function calls made.
  2106. warnflag : int
  2107. Integer warning flag:
  2108. 1 : Maximum number of function evaluations.
  2109. 2 : Maximum number of iterations.
  2110. allvecs : list
  2111. List of solutions at each iteration.
  2112. See also
  2113. --------
  2114. minimize: Interface to unconstrained minimization algorithms for
  2115. multivariate functions. See the 'Powell' `method` in particular.
  2116. Notes
  2117. -----
  2118. Uses a modification of Powell's method to find the minimum of
  2119. a function of N variables. Powell's method is a conjugate
  2120. direction method.
  2121. The algorithm has two loops. The outer loop
  2122. merely iterates over the inner loop. The inner loop minimizes
  2123. over each current direction in the direction set. At the end
  2124. of the inner loop, if certain conditions are met, the direction
  2125. that gave the largest decrease is dropped and replaced with
  2126. the difference between the current estimated x and the estimated
  2127. x from the beginning of the inner-loop.
  2128. The technical conditions for replacing the direction of greatest
  2129. increase amount to checking that
  2130. 1. No further gain can be made along the direction of greatest increase
  2131. from that iteration.
  2132. 2. The direction of greatest increase accounted for a large sufficient
  2133. fraction of the decrease in the function value from that iteration of
  2134. the inner loop.
  2135. Examples
  2136. --------
  2137. >>> def f(x):
  2138. ... return x**2
  2139. >>> from scipy import optimize
  2140. >>> minimum = optimize.fmin_powell(f, -1)
  2141. Optimization terminated successfully.
  2142. Current function value: 0.000000
  2143. Iterations: 2
  2144. Function evaluations: 18
  2145. >>> minimum
  2146. array(0.0)
  2147. References
  2148. ----------
  2149. Powell M.J.D. (1964) An efficient method for finding the minimum of a
  2150. function of several variables without calculating derivatives,
  2151. Computer Journal, 7 (2):155-162.
  2152. Press W., Teukolsky S.A., Vetterling W.T., and Flannery B.P.:
  2153. Numerical Recipes (any edition), Cambridge University Press
  2154. """
  2155. opts = {'xtol': xtol,
  2156. 'ftol': ftol,
  2157. 'maxiter': maxiter,
  2158. 'maxfev': maxfun,
  2159. 'disp': disp,
  2160. 'direc': direc,
  2161. 'return_all': retall}
  2162. res = _minimize_powell(func, x0, args, callback=callback, **opts)
  2163. if full_output:
  2164. retlist = (res['x'], res['fun'], res['direc'], res['nit'],
  2165. res['nfev'], res['status'])
  2166. if retall:
  2167. retlist += (res['allvecs'], )
  2168. return retlist
  2169. else:
  2170. if retall:
  2171. return res['x'], res['allvecs']
  2172. else:
  2173. return res['x']
  2174. def _minimize_powell(func, x0, args=(), callback=None,
  2175. xtol=1e-4, ftol=1e-4, maxiter=None, maxfev=None,
  2176. disp=False, direc=None, return_all=False,
  2177. **unknown_options):
  2178. """
  2179. Minimization of scalar function of one or more variables using the
  2180. modified Powell algorithm.
  2181. Options
  2182. -------
  2183. disp : bool
  2184. Set to True to print convergence messages.
  2185. xtol : float
  2186. Relative error in solution `xopt` acceptable for convergence.
  2187. ftol : float
  2188. Relative error in ``fun(xopt)`` acceptable for convergence.
  2189. maxiter, maxfev : int
  2190. Maximum allowed number of iterations and function evaluations.
  2191. Will default to ``N*1000``, where ``N`` is the number of
  2192. variables, if neither `maxiter` or `maxfev` is set. If both
  2193. `maxiter` and `maxfev` are set, minimization will stop at the
  2194. first reached.
  2195. direc : ndarray
  2196. Initial set of direction vectors for the Powell method.
  2197. """
  2198. _check_unknown_options(unknown_options)
  2199. maxfun = maxfev
  2200. retall = return_all
  2201. # we need to use a mutable object here that we can update in the
  2202. # wrapper function
  2203. fcalls, func = wrap_function(func, args)
  2204. x = asarray(x0).flatten()
  2205. if retall:
  2206. allvecs = [x]
  2207. N = len(x)
  2208. # If neither are set, then set both to default
  2209. if maxiter is None and maxfun is None:
  2210. maxiter = N * 1000
  2211. maxfun = N * 1000
  2212. elif maxiter is None:
  2213. # Convert remaining Nones, to np.inf, unless the other is np.inf, in
  2214. # which case use the default to avoid unbounded iteration
  2215. if maxfun == np.inf:
  2216. maxiter = N * 1000
  2217. else:
  2218. maxiter = np.inf
  2219. elif maxfun is None:
  2220. if maxiter == np.inf:
  2221. maxfun = N * 1000
  2222. else:
  2223. maxfun = np.inf
  2224. if direc is None:
  2225. direc = eye(N, dtype=float)
  2226. else:
  2227. direc = asarray(direc, dtype=float)
  2228. fval = squeeze(func(x))
  2229. x1 = x.copy()
  2230. iter = 0
  2231. ilist = list(range(N))
  2232. while True:
  2233. fx = fval
  2234. bigind = 0
  2235. delta = 0.0
  2236. for i in ilist:
  2237. direc1 = direc[i]
  2238. fx2 = fval
  2239. fval, x, direc1 = _linesearch_powell(func, x, direc1,
  2240. tol=xtol * 100)
  2241. if (fx2 - fval) > delta:
  2242. delta = fx2 - fval
  2243. bigind = i
  2244. iter += 1
  2245. if callback is not None:
  2246. callback(x)
  2247. if retall:
  2248. allvecs.append(x)
  2249. bnd = ftol * (numpy.abs(fx) + numpy.abs(fval)) + 1e-20
  2250. if 2.0 * (fx - fval) <= bnd:
  2251. break
  2252. if fcalls[0] >= maxfun:
  2253. break
  2254. if iter >= maxiter:
  2255. break
  2256. # Construct the extrapolated point
  2257. direc1 = x - x1
  2258. x2 = 2*x - x1
  2259. x1 = x.copy()
  2260. fx2 = squeeze(func(x2))
  2261. if (fx > fx2):
  2262. t = 2.0*(fx + fx2 - 2.0*fval)
  2263. temp = (fx - fval - delta)
  2264. t *= temp*temp
  2265. temp = fx - fx2
  2266. t -= delta*temp*temp
  2267. if t < 0.0:
  2268. fval, x, direc1 = _linesearch_powell(func, x, direc1,
  2269. tol=xtol*100)
  2270. direc[bigind] = direc[-1]
  2271. direc[-1] = direc1
  2272. warnflag = 0
  2273. if fcalls[0] >= maxfun:
  2274. warnflag = 1
  2275. msg = _status_message['maxfev']
  2276. if disp:
  2277. print("Warning: " + msg)
  2278. elif iter >= maxiter:
  2279. warnflag = 2
  2280. msg = _status_message['maxiter']
  2281. if disp:
  2282. print("Warning: " + msg)
  2283. else:
  2284. msg = _status_message['success']
  2285. if disp:
  2286. print(msg)
  2287. print(" Current function value: %f" % fval)
  2288. print(" Iterations: %d" % iter)
  2289. print(" Function evaluations: %d" % fcalls[0])
  2290. x = squeeze(x)
  2291. result = OptimizeResult(fun=fval, direc=direc, nit=iter, nfev=fcalls[0],
  2292. status=warnflag, success=(warnflag == 0),
  2293. message=msg, x=x)
  2294. if retall:
  2295. result['allvecs'] = allvecs
  2296. return result
  2297. def _endprint(x, flag, fval, maxfun, xtol, disp):
  2298. if flag == 0:
  2299. if disp > 1:
  2300. print("\nOptimization terminated successfully;\n"
  2301. "The returned value satisfies the termination criteria\n"
  2302. "(using xtol = ", xtol, ")")
  2303. if flag == 1:
  2304. if disp:
  2305. print("\nMaximum number of function evaluations exceeded --- "
  2306. "increase maxfun argument.\n")
  2307. return
  2308. def brute(func, ranges, args=(), Ns=20, full_output=0, finish=fmin,
  2309. disp=False):
  2310. """Minimize a function over a given range by brute force.
  2311. Uses the "brute force" method, i.e. computes the function's value
  2312. at each point of a multidimensional grid of points, to find the global
  2313. minimum of the function.
  2314. The function is evaluated everywhere in the range with the datatype of the
  2315. first call to the function, as enforced by the ``vectorize`` NumPy
  2316. function. The value and type of the function evaluation returned when
  2317. ``full_output=True`` are affected in addition by the ``finish`` argument
  2318. (see Notes).
  2319. The brute force approach is inefficient because the number of grid points
  2320. increases exponentially - the number of grid points to evaluate is
  2321. ``Ns ** len(x)``. Consequently, even with coarse grid spacing, even
  2322. moderately sized problems can take a long time to run, and/or run into
  2323. memory limitations.
  2324. Parameters
  2325. ----------
  2326. func : callable
  2327. The objective function to be minimized. Must be in the
  2328. form ``f(x, *args)``, where ``x`` is the argument in
  2329. the form of a 1-D array and ``args`` is a tuple of any
  2330. additional fixed parameters needed to completely specify
  2331. the function.
  2332. ranges : tuple
  2333. Each component of the `ranges` tuple must be either a
  2334. "slice object" or a range tuple of the form ``(low, high)``.
  2335. The program uses these to create the grid of points on which
  2336. the objective function will be computed. See `Note 2` for
  2337. more detail.
  2338. args : tuple, optional
  2339. Any additional fixed parameters needed to completely specify
  2340. the function.
  2341. Ns : int, optional
  2342. Number of grid points along the axes, if not otherwise
  2343. specified. See `Note2`.
  2344. full_output : bool, optional
  2345. If True, return the evaluation grid and the objective function's
  2346. values on it.
  2347. finish : callable, optional
  2348. An optimization function that is called with the result of brute force
  2349. minimization as initial guess. `finish` should take `func` and
  2350. the initial guess as positional arguments, and take `args` as
  2351. keyword arguments. It may additionally take `full_output`
  2352. and/or `disp` as keyword arguments. Use None if no "polishing"
  2353. function is to be used. See Notes for more details.
  2354. disp : bool, optional
  2355. Set to True to print convergence messages.
  2356. Returns
  2357. -------
  2358. x0 : ndarray
  2359. A 1-D array containing the coordinates of a point at which the
  2360. objective function had its minimum value. (See `Note 1` for
  2361. which point is returned.)
  2362. fval : float
  2363. Function value at the point `x0`. (Returned when `full_output` is
  2364. True.)
  2365. grid : tuple
  2366. Representation of the evaluation grid. It has the same
  2367. length as `x0`. (Returned when `full_output` is True.)
  2368. Jout : ndarray
  2369. Function values at each point of the evaluation
  2370. grid, `i.e.`, ``Jout = func(*grid)``. (Returned
  2371. when `full_output` is True.)
  2372. See Also
  2373. --------
  2374. basinhopping, differential_evolution
  2375. Notes
  2376. -----
  2377. *Note 1*: The program finds the gridpoint at which the lowest value
  2378. of the objective function occurs. If `finish` is None, that is the
  2379. point returned. When the global minimum occurs within (or not very far
  2380. outside) the grid's boundaries, and the grid is fine enough, that
  2381. point will be in the neighborhood of the global minimum.
  2382. However, users often employ some other optimization program to
  2383. "polish" the gridpoint values, `i.e.`, to seek a more precise
  2384. (local) minimum near `brute's` best gridpoint.
  2385. The `brute` function's `finish` option provides a convenient way to do
  2386. that. Any polishing program used must take `brute's` output as its
  2387. initial guess as a positional argument, and take `brute's` input values
  2388. for `args` as keyword arguments, otherwise an error will be raised.
  2389. It may additionally take `full_output` and/or `disp` as keyword arguments.
  2390. `brute` assumes that the `finish` function returns either an
  2391. `OptimizeResult` object or a tuple in the form:
  2392. ``(xmin, Jmin, ... , statuscode)``, where ``xmin`` is the minimizing
  2393. value of the argument, ``Jmin`` is the minimum value of the objective
  2394. function, "..." may be some other returned values (which are not used
  2395. by `brute`), and ``statuscode`` is the status code of the `finish` program.
  2396. Note that when `finish` is not None, the values returned are those
  2397. of the `finish` program, *not* the gridpoint ones. Consequently,
  2398. while `brute` confines its search to the input grid points,
  2399. the `finish` program's results usually will not coincide with any
  2400. gridpoint, and may fall outside the grid's boundary. Thus, if a
  2401. minimum only needs to be found over the provided grid points, make
  2402. sure to pass in `finish=None`.
  2403. *Note 2*: The grid of points is a `numpy.mgrid` object.
  2404. For `brute` the `ranges` and `Ns` inputs have the following effect.
  2405. Each component of the `ranges` tuple can be either a slice object or a
  2406. two-tuple giving a range of values, such as (0, 5). If the component is a
  2407. slice object, `brute` uses it directly. If the component is a two-tuple
  2408. range, `brute` internally converts it to a slice object that interpolates
  2409. `Ns` points from its low-value to its high-value, inclusive.
  2410. Examples
  2411. --------
  2412. We illustrate the use of `brute` to seek the global minimum of a function
  2413. of two variables that is given as the sum of a positive-definite
  2414. quadratic and two deep "Gaussian-shaped" craters. Specifically, define
  2415. the objective function `f` as the sum of three other functions,
  2416. ``f = f1 + f2 + f3``. We suppose each of these has a signature
  2417. ``(z, *params)``, where ``z = (x, y)``, and ``params`` and the functions
  2418. are as defined below.
  2419. >>> params = (2, 3, 7, 8, 9, 10, 44, -1, 2, 26, 1, -2, 0.5)
  2420. >>> def f1(z, *params):
  2421. ... x, y = z
  2422. ... a, b, c, d, e, f, g, h, i, j, k, l, scale = params
  2423. ... return (a * x**2 + b * x * y + c * y**2 + d*x + e*y + f)
  2424. >>> def f2(z, *params):
  2425. ... x, y = z
  2426. ... a, b, c, d, e, f, g, h, i, j, k, l, scale = params
  2427. ... return (-g*np.exp(-((x-h)**2 + (y-i)**2) / scale))
  2428. >>> def f3(z, *params):
  2429. ... x, y = z
  2430. ... a, b, c, d, e, f, g, h, i, j, k, l, scale = params
  2431. ... return (-j*np.exp(-((x-k)**2 + (y-l)**2) / scale))
  2432. >>> def f(z, *params):
  2433. ... return f1(z, *params) + f2(z, *params) + f3(z, *params)
  2434. Thus, the objective function may have local minima near the minimum
  2435. of each of the three functions of which it is composed. To
  2436. use `fmin` to polish its gridpoint result, we may then continue as
  2437. follows:
  2438. >>> rranges = (slice(-4, 4, 0.25), slice(-4, 4, 0.25))
  2439. >>> from scipy import optimize
  2440. >>> resbrute = optimize.brute(f, rranges, args=params, full_output=True,
  2441. ... finish=optimize.fmin)
  2442. >>> resbrute[0] # global minimum
  2443. array([-1.05665192, 1.80834843])
  2444. >>> resbrute[1] # function value at global minimum
  2445. -3.4085818767
  2446. Note that if `finish` had been set to None, we would have gotten the
  2447. gridpoint [-1.0 1.75] where the rounded function value is -2.892.
  2448. """
  2449. N = len(ranges)
  2450. if N > 40:
  2451. raise ValueError("Brute Force not possible with more "
  2452. "than 40 variables.")
  2453. lrange = list(ranges)
  2454. for k in range(N):
  2455. if type(lrange[k]) is not type(slice(None)):
  2456. if len(lrange[k]) < 3:
  2457. lrange[k] = tuple(lrange[k]) + (complex(Ns),)
  2458. lrange[k] = slice(*lrange[k])
  2459. if (N == 1):
  2460. lrange = lrange[0]
  2461. def _scalarfunc(*params):
  2462. params = asarray(params).flatten()
  2463. return func(params, *args)
  2464. vecfunc = vectorize(_scalarfunc)
  2465. grid = mgrid[lrange]
  2466. if (N == 1):
  2467. grid = (grid,)
  2468. Jout = vecfunc(*grid)
  2469. Nshape = shape(Jout)
  2470. indx = argmin(Jout.ravel(), axis=-1)
  2471. Nindx = zeros(N, int)
  2472. xmin = zeros(N, float)
  2473. for k in range(N - 1, -1, -1):
  2474. thisN = Nshape[k]
  2475. Nindx[k] = indx % Nshape[k]
  2476. indx = indx // thisN
  2477. for k in range(N):
  2478. xmin[k] = grid[k][tuple(Nindx)]
  2479. Jmin = Jout[tuple(Nindx)]
  2480. if (N == 1):
  2481. grid = grid[0]
  2482. xmin = xmin[0]
  2483. if callable(finish):
  2484. # set up kwargs for `finish` function
  2485. finish_args = _getargspec(finish).args
  2486. finish_kwargs = dict()
  2487. if 'full_output' in finish_args:
  2488. finish_kwargs['full_output'] = 1
  2489. if 'disp' in finish_args:
  2490. finish_kwargs['disp'] = disp
  2491. elif 'options' in finish_args:
  2492. # pass 'disp' as `options`
  2493. # (e.g. if `finish` is `minimize`)
  2494. finish_kwargs['options'] = {'disp': disp}
  2495. # run minimizer
  2496. res = finish(func, xmin, args=args, **finish_kwargs)
  2497. if isinstance(res, OptimizeResult):
  2498. xmin = res.x
  2499. Jmin = res.fun
  2500. success = res.success
  2501. else:
  2502. xmin = res[0]
  2503. Jmin = res[1]
  2504. success = res[-1] == 0
  2505. if not success:
  2506. if disp:
  2507. print("Warning: Either final optimization did not succeed "
  2508. "or `finish` does not return `statuscode` as its last "
  2509. "argument.")
  2510. if full_output:
  2511. return xmin, Jmin, grid, Jout
  2512. else:
  2513. return xmin
  2514. def show_options(solver=None, method=None, disp=True):
  2515. """
  2516. Show documentation for additional options of optimization solvers.
  2517. These are method-specific options that can be supplied through the
  2518. ``options`` dict.
  2519. Parameters
  2520. ----------
  2521. solver : str
  2522. Type of optimization solver. One of 'minimize', 'minimize_scalar',
  2523. 'root', or 'linprog'.
  2524. method : str, optional
  2525. If not given, shows all methods of the specified solver. Otherwise,
  2526. show only the options for the specified method. Valid values
  2527. corresponds to methods' names of respective solver (e.g. 'BFGS' for
  2528. 'minimize').
  2529. disp : bool, optional
  2530. Whether to print the result rather than returning it.
  2531. Returns
  2532. -------
  2533. text
  2534. Either None (for disp=False) or the text string (disp=True)
  2535. Notes
  2536. -----
  2537. The solver-specific methods are:
  2538. `scipy.optimize.minimize`
  2539. - :ref:`Nelder-Mead <optimize.minimize-neldermead>`
  2540. - :ref:`Powell <optimize.minimize-powell>`
  2541. - :ref:`CG <optimize.minimize-cg>`
  2542. - :ref:`BFGS <optimize.minimize-bfgs>`
  2543. - :ref:`Newton-CG <optimize.minimize-newtoncg>`
  2544. - :ref:`L-BFGS-B <optimize.minimize-lbfgsb>`
  2545. - :ref:`TNC <optimize.minimize-tnc>`
  2546. - :ref:`COBYLA <optimize.minimize-cobyla>`
  2547. - :ref:`SLSQP <optimize.minimize-slsqp>`
  2548. - :ref:`dogleg <optimize.minimize-dogleg>`
  2549. - :ref:`trust-ncg <optimize.minimize-trustncg>`
  2550. `scipy.optimize.root`
  2551. - :ref:`hybr <optimize.root-hybr>`
  2552. - :ref:`lm <optimize.root-lm>`
  2553. - :ref:`broyden1 <optimize.root-broyden1>`
  2554. - :ref:`broyden2 <optimize.root-broyden2>`
  2555. - :ref:`anderson <optimize.root-anderson>`
  2556. - :ref:`linearmixing <optimize.root-linearmixing>`
  2557. - :ref:`diagbroyden <optimize.root-diagbroyden>`
  2558. - :ref:`excitingmixing <optimize.root-excitingmixing>`
  2559. - :ref:`krylov <optimize.root-krylov>`
  2560. - :ref:`df-sane <optimize.root-dfsane>`
  2561. `scipy.optimize.minimize_scalar`
  2562. - :ref:`brent <optimize.minimize_scalar-brent>`
  2563. - :ref:`golden <optimize.minimize_scalar-golden>`
  2564. - :ref:`bounded <optimize.minimize_scalar-bounded>`
  2565. `scipy.optimize.linprog`
  2566. - :ref:`simplex <optimize.linprog-simplex>`
  2567. - :ref:`interior-point <optimize.linprog-interior-point>`
  2568. """
  2569. import textwrap
  2570. doc_routines = {
  2571. 'minimize': (
  2572. ('bfgs', 'scipy.optimize.optimize._minimize_bfgs'),
  2573. ('cg', 'scipy.optimize.optimize._minimize_cg'),
  2574. ('cobyla', 'scipy.optimize.cobyla._minimize_cobyla'),
  2575. ('dogleg', 'scipy.optimize._trustregion_dogleg._minimize_dogleg'),
  2576. ('l-bfgs-b', 'scipy.optimize.lbfgsb._minimize_lbfgsb'),
  2577. ('nelder-mead', 'scipy.optimize.optimize._minimize_neldermead'),
  2578. ('newton-cg', 'scipy.optimize.optimize._minimize_newtoncg'),
  2579. ('powell', 'scipy.optimize.optimize._minimize_powell'),
  2580. ('slsqp', 'scipy.optimize.slsqp._minimize_slsqp'),
  2581. ('tnc', 'scipy.optimize.tnc._minimize_tnc'),
  2582. ('trust-ncg', 'scipy.optimize._trustregion_ncg._minimize_trust_ncg'),
  2583. ),
  2584. 'root': (
  2585. ('hybr', 'scipy.optimize.minpack._root_hybr'),
  2586. ('lm', 'scipy.optimize._root._root_leastsq'),
  2587. ('broyden1', 'scipy.optimize._root._root_broyden1_doc'),
  2588. ('broyden2', 'scipy.optimize._root._root_broyden2_doc'),
  2589. ('anderson', 'scipy.optimize._root._root_anderson_doc'),
  2590. ('diagbroyden', 'scipy.optimize._root._root_diagbroyden_doc'),
  2591. ('excitingmixing', 'scipy.optimize._root._root_excitingmixing_doc'),
  2592. ('linearmixing', 'scipy.optimize._root._root_linearmixing_doc'),
  2593. ('krylov', 'scipy.optimize._root._root_krylov_doc'),
  2594. ('df-sane', 'scipy.optimize._spectral._root_df_sane'),
  2595. ),
  2596. 'root_scalar': (
  2597. ('bisect', 'scipy.optimize._root_scalar._root_scalar_bisect_doc'),
  2598. ('brentq', 'scipy.optimize._root_scalar._root_scalar_brentq_doc'),
  2599. ('brenth', 'scipy.optimize._root_scalar._root_scalar_brenth_doc'),
  2600. ('ridder', 'scipy.optimize._root_scalar._root_scalar_ridder_doc'),
  2601. ('toms748', 'scipy.optimize._root_scalar._root_scalar_toms748_doc'),
  2602. ('secant', 'scipy.optimize._root_scalar._root_scalar_secant_doc'),
  2603. ('newton', 'scipy.optimize._root_scalar._root_scalar_newton_doc'),
  2604. ('halley', 'scipy.optimize._root_scalar._root_scalar_halley_doc'),
  2605. ),
  2606. 'linprog': (
  2607. ('simplex', 'scipy.optimize._linprog._linprog_simplex'),
  2608. ('interior-point', 'scipy.optimize._linprog._linprog_ip'),
  2609. ),
  2610. 'minimize_scalar': (
  2611. ('brent', 'scipy.optimize.optimize._minimize_scalar_brent'),
  2612. ('bounded', 'scipy.optimize.optimize._minimize_scalar_bounded'),
  2613. ('golden', 'scipy.optimize.optimize._minimize_scalar_golden'),
  2614. ),
  2615. }
  2616. if solver is None:
  2617. text = ["\n\n\n========\n", "minimize\n", "========\n"]
  2618. text.append(show_options('minimize', disp=False))
  2619. text.extend(["\n\n===============\n", "minimize_scalar\n",
  2620. "===============\n"])
  2621. text.append(show_options('minimize_scalar', disp=False))
  2622. text.extend(["\n\n\n====\n", "root\n",
  2623. "====\n"])
  2624. text.append(show_options('root', disp=False))
  2625. text.extend(['\n\n\n=======\n', 'linprog\n',
  2626. '=======\n'])
  2627. text.append(show_options('linprog', disp=False))
  2628. text = "".join(text)
  2629. else:
  2630. solver = solver.lower()
  2631. if solver not in doc_routines:
  2632. raise ValueError('Unknown solver %r' % (solver,))
  2633. if method is None:
  2634. text = []
  2635. for name, _ in doc_routines[solver]:
  2636. text.extend(["\n\n" + name, "\n" + "="*len(name) + "\n\n"])
  2637. text.append(show_options(solver, name, disp=False))
  2638. text = "".join(text)
  2639. else:
  2640. method = method.lower()
  2641. methods = dict(doc_routines[solver])
  2642. if method not in methods:
  2643. raise ValueError("Unknown method %r" % (method,))
  2644. name = methods[method]
  2645. # Import function object
  2646. parts = name.split('.')
  2647. mod_name = ".".join(parts[:-1])
  2648. __import__(mod_name)
  2649. obj = getattr(sys.modules[mod_name], parts[-1])
  2650. # Get doc
  2651. doc = obj.__doc__
  2652. if doc is not None:
  2653. text = textwrap.dedent(doc).strip()
  2654. else:
  2655. text = ""
  2656. if disp:
  2657. print(text)
  2658. return
  2659. else:
  2660. return text
  2661. def main():
  2662. import time
  2663. times = []
  2664. algor = []
  2665. x0 = [0.8, 1.2, 0.7]
  2666. print("Nelder-Mead Simplex")
  2667. print("===================")
  2668. start = time.time()
  2669. x = fmin(rosen, x0)
  2670. print(x)
  2671. times.append(time.time() - start)
  2672. algor.append('Nelder-Mead Simplex\t')
  2673. print()
  2674. print("Powell Direction Set Method")
  2675. print("===========================")
  2676. start = time.time()
  2677. x = fmin_powell(rosen, x0)
  2678. print(x)
  2679. times.append(time.time() - start)
  2680. algor.append('Powell Direction Set Method.')
  2681. print()
  2682. print("Nonlinear CG")
  2683. print("============")
  2684. start = time.time()
  2685. x = fmin_cg(rosen, x0, fprime=rosen_der, maxiter=200)
  2686. print(x)
  2687. times.append(time.time() - start)
  2688. algor.append('Nonlinear CG \t')
  2689. print()
  2690. print("BFGS Quasi-Newton")
  2691. print("=================")
  2692. start = time.time()
  2693. x = fmin_bfgs(rosen, x0, fprime=rosen_der, maxiter=80)
  2694. print(x)
  2695. times.append(time.time() - start)
  2696. algor.append('BFGS Quasi-Newton\t')
  2697. print()
  2698. print("BFGS approximate gradient")
  2699. print("=========================")
  2700. start = time.time()
  2701. x = fmin_bfgs(rosen, x0, gtol=1e-4, maxiter=100)
  2702. print(x)
  2703. times.append(time.time() - start)
  2704. algor.append('BFGS without gradient\t')
  2705. print()
  2706. print("Newton-CG with Hessian product")
  2707. print("==============================")
  2708. start = time.time()
  2709. x = fmin_ncg(rosen, x0, rosen_der, fhess_p=rosen_hess_prod, maxiter=80)
  2710. print(x)
  2711. times.append(time.time() - start)
  2712. algor.append('Newton-CG with hessian product')
  2713. print()
  2714. print("Newton-CG with full Hessian")
  2715. print("===========================")
  2716. start = time.time()
  2717. x = fmin_ncg(rosen, x0, rosen_der, fhess=rosen_hess, maxiter=80)
  2718. print(x)
  2719. times.append(time.time() - start)
  2720. algor.append('Newton-CG with full hessian')
  2721. print()
  2722. print("\nMinimizing the Rosenbrock function of order 3\n")
  2723. print(" Algorithm \t\t\t Seconds")
  2724. print("===========\t\t\t =========")
  2725. for k in range(len(algor)):
  2726. print(algor[k], "\t -- ", times[k])
  2727. if __name__ == "__main__":
  2728. main()