_zeros_py.py 55 KB

12345678910111213141516171819202122232425262728293031323334353637383940414243444546474849505152535455565758596061626364656667686970717273747576777879808182838485868788899091929394959697989910010110210310410510610710810911011111211311411511611711811912012112212312412512612712812913013113213313413513613713813914014114214314414514614714814915015115215315415515615715815916016116216316416516616716816917017117217317417517617717817918018118218318418518618718818919019119219319419519619719819920020120220320420520620720820921021121221321421521621721821922022122222322422522622722822923023123223323423523623723823924024124224324424524624724824925025125225325425525625725825926026126226326426526626726826927027127227327427527627727827928028128228328428528628728828929029129229329429529629729829930030130230330430530630730830931031131231331431531631731831932032132232332432532632732832933033133233333433533633733833934034134234334434534634734834935035135235335435535635735835936036136236336436536636736836937037137237337437537637737837938038138238338438538638738838939039139239339439539639739839940040140240340440540640740840941041141241341441541641741841942042142242342442542642742842943043143243343443543643743843944044144244344444544644744844945045145245345445545645745845946046146246346446546646746846947047147247347447547647747847948048148248348448548648748848949049149249349449549649749849950050150250350450550650750850951051151251351451551651751851952052152252352452552652752852953053153253353453553653753853954054154254354454554654754854955055155255355455555655755855956056156256356456556656756856957057157257357457557657757857958058158258358458558658758858959059159259359459559659759859960060160260360460560660760860961061161261361461561661761861962062162262362462562662762862963063163263363463563663763863964064164264364464564664764864965065165265365465565665765865966066166266366466566666766866967067167267367467567667767867968068168268368468568668768868969069169269369469569669769869970070170270370470570670770870971071171271371471571671771871972072172272372472572672772872973073173273373473573673773873974074174274374474574674774874975075175275375475575675775875976076176276376476576676776876977077177277377477577677777877978078178278378478578678778878979079179279379479579679779879980080180280380480580680780880981081181281381481581681781881982082182282382482582682782882983083183283383483583683783883984084184284384484584684784884985085185285385485585685785885986086186286386486586686786886987087187287387487587687787887988088188288388488588688788888989089189289389489589689789889990090190290390490590690790890991091191291391491591691791891992092192292392492592692792892993093193293393493593693793893994094194294394494594694794894995095195295395495595695795895996096196296396496596696796896997097197297397497597697797897998098198298398498598698798898999099199299399499599699799899910001001100210031004100510061007100810091010101110121013101410151016101710181019102010211022102310241025102610271028102910301031103210331034103510361037103810391040104110421043104410451046104710481049105010511052105310541055105610571058105910601061106210631064106510661067106810691070107110721073107410751076107710781079108010811082108310841085108610871088108910901091109210931094109510961097109810991100110111021103110411051106110711081109111011111112111311141115111611171118111911201121112211231124112511261127112811291130113111321133113411351136113711381139114011411142114311441145114611471148114911501151115211531154115511561157115811591160116111621163116411651166116711681169117011711172117311741175117611771178117911801181118211831184118511861187118811891190119111921193119411951196119711981199120012011202120312041205120612071208120912101211121212131214121512161217121812191220122112221223122412251226122712281229123012311232123312341235123612371238123912401241124212431244124512461247124812491250125112521253125412551256125712581259126012611262126312641265126612671268126912701271127212731274127512761277127812791280128112821283128412851286128712881289129012911292129312941295129612971298129913001301130213031304130513061307130813091310131113121313131413151316131713181319132013211322132313241325132613271328132913301331133213331334133513361337133813391340134113421343134413451346134713481349135013511352135313541355135613571358135913601361136213631364136513661367136813691370137113721373137413751376137713781379138013811382138313841385138613871388138913901391139213931394139513961397139813991400140114021403140414051406140714081409141014111412141314141415141614171418141914201421142214231424142514261427142814291430143114321433143414351436143714381439144014411442144314441445144614471448144914501451145214531454145514561457145814591460146114621463146414651466146714681469147014711472147314741475
  1. import warnings
  2. from collections import namedtuple
  3. import operator
  4. from . import _zeros
  5. from ._optimize import OptimizeResult
  6. import numpy as np
  7. _iter = 100
  8. _xtol = 2e-12
  9. _rtol = 4 * np.finfo(float).eps
  10. __all__ = ['newton', 'bisect', 'ridder', 'brentq', 'brenth', 'toms748',
  11. 'RootResults']
  12. # Must agree with CONVERGED, SIGNERR, CONVERR, ... in zeros.h
  13. _ECONVERGED = 0
  14. _ESIGNERR = -1 # used in _chandrupatla
  15. _ECONVERR = -2
  16. _EVALUEERR = -3
  17. _ECALLBACK = -4
  18. _EINPROGRESS = 1
  19. CONVERGED = 'converged'
  20. SIGNERR = 'sign error'
  21. CONVERR = 'convergence error'
  22. VALUEERR = 'value error'
  23. INPROGRESS = 'No error'
  24. flag_map = {_ECONVERGED: CONVERGED, _ESIGNERR: SIGNERR, _ECONVERR: CONVERR,
  25. _EVALUEERR: VALUEERR, _EINPROGRESS: INPROGRESS}
  26. class RootResults(OptimizeResult):
  27. """Represents the root finding result.
  28. Attributes
  29. ----------
  30. root : float
  31. Estimated root location.
  32. iterations : int
  33. Number of iterations needed to find the root.
  34. function_calls : int
  35. Number of times the function was called.
  36. converged : bool
  37. True if the routine converged.
  38. flag : str
  39. Description of the cause of termination.
  40. method : str
  41. Root finding method used.
  42. """
  43. def __init__(self, root, iterations, function_calls, flag, method):
  44. self.root = root
  45. self.iterations = iterations
  46. self.function_calls = function_calls
  47. self.converged = flag == _ECONVERGED
  48. if flag in flag_map:
  49. self.flag = flag_map[flag]
  50. else:
  51. self.flag = flag
  52. self.method = method
  53. def results_c(full_output, r, method):
  54. if full_output:
  55. x, funcalls, iterations, flag = r
  56. results = RootResults(root=x,
  57. iterations=iterations,
  58. function_calls=funcalls,
  59. flag=flag, method=method)
  60. return x, results
  61. else:
  62. return r
  63. def _results_select(full_output, r, method):
  64. """Select from a tuple of (root, funccalls, iterations, flag)"""
  65. x, funcalls, iterations, flag = r
  66. if full_output:
  67. results = RootResults(root=x,
  68. iterations=iterations,
  69. function_calls=funcalls,
  70. flag=flag, method=method)
  71. return x, results
  72. return x
  73. def _wrap_nan_raise(f):
  74. def f_raise(x, *args):
  75. fx = f(x, *args)
  76. f_raise._function_calls += 1
  77. if np.isnan(fx):
  78. msg = (f'The function value at x={x} is NaN; '
  79. 'solver cannot continue.')
  80. err = ValueError(msg)
  81. err._x = x
  82. err._function_calls = f_raise._function_calls
  83. raise err
  84. return fx
  85. f_raise._function_calls = 0
  86. return f_raise
  87. def newton(func, x0, fprime=None, args=(), tol=1.48e-8, maxiter=50,
  88. fprime2=None, x1=None, rtol=0.0,
  89. full_output=False, disp=True):
  90. """
  91. Find a root of a real or complex function using the Newton-Raphson
  92. (or secant or Halley's) method.
  93. Find a root of the scalar-valued function `func` given a nearby scalar
  94. starting point `x0`.
  95. The Newton-Raphson method is used if the derivative `fprime` of `func`
  96. is provided, otherwise the secant method is used. If the second order
  97. derivative `fprime2` of `func` is also provided, then Halley's method is
  98. used.
  99. If `x0` is a sequence with more than one item, `newton` returns an array:
  100. the roots of the function from each (scalar) starting point in `x0`.
  101. In this case, `func` must be vectorized to return a sequence or array of
  102. the same shape as its first argument. If `fprime` (`fprime2`) is given,
  103. then its return must also have the same shape: each element is the first
  104. (second) derivative of `func` with respect to its only variable evaluated
  105. at each element of its first argument.
  106. `newton` is for finding roots of a scalar-valued functions of a single
  107. variable. For problems involving several variables, see `root`.
  108. Parameters
  109. ----------
  110. func : callable
  111. The function whose root is wanted. It must be a function of a
  112. single variable of the form ``f(x,a,b,c,...)``, where ``a,b,c,...``
  113. are extra arguments that can be passed in the `args` parameter.
  114. x0 : float, sequence, or ndarray
  115. An initial estimate of the root that should be somewhere near the
  116. actual root. If not scalar, then `func` must be vectorized and return
  117. a sequence or array of the same shape as its first argument.
  118. fprime : callable, optional
  119. The derivative of the function when available and convenient. If it
  120. is None (default), then the secant method is used.
  121. args : tuple, optional
  122. Extra arguments to be used in the function call.
  123. tol : float, optional
  124. The allowable error of the root's value. If `func` is complex-valued,
  125. a larger `tol` is recommended as both the real and imaginary parts
  126. of `x` contribute to ``|x - x0|``.
  127. maxiter : int, optional
  128. Maximum number of iterations.
  129. fprime2 : callable, optional
  130. The second order derivative of the function when available and
  131. convenient. If it is None (default), then the normal Newton-Raphson
  132. or the secant method is used. If it is not None, then Halley's method
  133. is used.
  134. x1 : float, optional
  135. Another estimate of the root that should be somewhere near the
  136. actual root. Used if `fprime` is not provided.
  137. rtol : float, optional
  138. Tolerance (relative) for termination.
  139. full_output : bool, optional
  140. If `full_output` is False (default), the root is returned.
  141. If True and `x0` is scalar, the return value is ``(x, r)``, where ``x``
  142. is the root and ``r`` is a `RootResults` object.
  143. If True and `x0` is non-scalar, the return value is ``(x, converged,
  144. zero_der)`` (see Returns section for details).
  145. disp : bool, optional
  146. If True, raise a RuntimeError if the algorithm didn't converge, with
  147. the error message containing the number of iterations and current
  148. function value. Otherwise, the convergence status is recorded in a
  149. `RootResults` return object.
  150. Ignored if `x0` is not scalar.
  151. *Note: this has little to do with displaying, however,
  152. the `disp` keyword cannot be renamed for backwards compatibility.*
  153. Returns
  154. -------
  155. root : float, sequence, or ndarray
  156. Estimated location where function is zero.
  157. r : `RootResults`, optional
  158. Present if ``full_output=True`` and `x0` is scalar.
  159. Object containing information about the convergence. In particular,
  160. ``r.converged`` is True if the routine converged.
  161. converged : ndarray of bool, optional
  162. Present if ``full_output=True`` and `x0` is non-scalar.
  163. For vector functions, indicates which elements converged successfully.
  164. zero_der : ndarray of bool, optional
  165. Present if ``full_output=True`` and `x0` is non-scalar.
  166. For vector functions, indicates which elements had a zero derivative.
  167. See Also
  168. --------
  169. root_scalar : interface to root solvers for scalar functions
  170. root : interface to root solvers for multi-input, multi-output functions
  171. Notes
  172. -----
  173. The convergence rate of the Newton-Raphson method is quadratic,
  174. the Halley method is cubic, and the secant method is
  175. sub-quadratic. This means that if the function is well-behaved
  176. the actual error in the estimated root after the nth iteration
  177. is approximately the square (cube for Halley) of the error
  178. after the (n-1)th step. However, the stopping criterion used
  179. here is the step size and there is no guarantee that a root
  180. has been found. Consequently, the result should be verified.
  181. Safer algorithms are brentq, brenth, ridder, and bisect,
  182. but they all require that the root first be bracketed in an
  183. interval where the function changes sign. The brentq algorithm
  184. is recommended for general use in one dimensional problems
  185. when such an interval has been found.
  186. When `newton` is used with arrays, it is best suited for the following
  187. types of problems:
  188. * The initial guesses, `x0`, are all relatively the same distance from
  189. the roots.
  190. * Some or all of the extra arguments, `args`, are also arrays so that a
  191. class of similar problems can be solved together.
  192. * The size of the initial guesses, `x0`, is larger than O(100) elements.
  193. Otherwise, a naive loop may perform as well or better than a vector.
  194. Examples
  195. --------
  196. >>> import numpy as np
  197. >>> import matplotlib.pyplot as plt
  198. >>> from scipy import optimize
  199. >>> def f(x):
  200. ... return (x**3 - 1) # only one real root at x = 1
  201. ``fprime`` is not provided, use the secant method:
  202. >>> root = optimize.newton(f, 1.5)
  203. >>> root
  204. 1.0000000000000016
  205. >>> root = optimize.newton(f, 1.5, fprime2=lambda x: 6 * x)
  206. >>> root
  207. 1.0000000000000016
  208. Only ``fprime`` is provided, use the Newton-Raphson method:
  209. >>> root = optimize.newton(f, 1.5, fprime=lambda x: 3 * x**2)
  210. >>> root
  211. 1.0
  212. Both ``fprime2`` and ``fprime`` are provided, use Halley's method:
  213. >>> root = optimize.newton(f, 1.5, fprime=lambda x: 3 * x**2,
  214. ... fprime2=lambda x: 6 * x)
  215. >>> root
  216. 1.0
  217. When we want to find roots for a set of related starting values and/or
  218. function parameters, we can provide both of those as an array of inputs:
  219. >>> f = lambda x, a: x**3 - a
  220. >>> fder = lambda x, a: 3 * x**2
  221. >>> rng = np.random.default_rng()
  222. >>> x = rng.standard_normal(100)
  223. >>> a = np.arange(-50, 50)
  224. >>> vec_res = optimize.newton(f, x, fprime=fder, args=(a, ), maxiter=200)
  225. The above is the equivalent of solving for each value in ``(x, a)``
  226. separately in a for-loop, just faster:
  227. >>> loop_res = [optimize.newton(f, x0, fprime=fder, args=(a0,),
  228. ... maxiter=200)
  229. ... for x0, a0 in zip(x, a)]
  230. >>> np.allclose(vec_res, loop_res)
  231. True
  232. Plot the results found for all values of ``a``:
  233. >>> analytical_result = np.sign(a) * np.abs(a)**(1/3)
  234. >>> fig, ax = plt.subplots()
  235. >>> ax.plot(a, analytical_result, 'o')
  236. >>> ax.plot(a, vec_res, '.')
  237. >>> ax.set_xlabel('$a$')
  238. >>> ax.set_ylabel('$x$ where $f(x, a)=0$')
  239. >>> plt.show()
  240. """
  241. if tol <= 0:
  242. raise ValueError(f"tol too small ({tol:g} <= 0)")
  243. maxiter = operator.index(maxiter)
  244. if maxiter < 1:
  245. raise ValueError("maxiter must be greater than 0")
  246. if np.size(x0) > 1:
  247. return _array_newton(func, x0, fprime, args, tol, maxiter, fprime2,
  248. full_output)
  249. # Convert to float (don't use float(x0); this works also for complex x0)
  250. # Use np.asarray because we want x0 to be a numpy object, not a Python
  251. # object. e.g. np.complex(1+1j) > 0 is possible, but (1 + 1j) > 0 raises
  252. # a TypeError
  253. x0 = np.asarray(x0)[()] * 1.0
  254. p0 = x0
  255. funcalls = 0
  256. if fprime is not None:
  257. # Newton-Raphson method
  258. method = "newton"
  259. for itr in range(maxiter):
  260. # first evaluate fval
  261. fval = func(p0, *args)
  262. funcalls += 1
  263. # If fval is 0, a root has been found, then terminate
  264. if fval == 0:
  265. return _results_select(
  266. full_output, (p0, funcalls, itr, _ECONVERGED), method)
  267. fder = fprime(p0, *args)
  268. funcalls += 1
  269. if fder == 0:
  270. msg = "Derivative was zero."
  271. if disp:
  272. msg += (
  273. f" Failed to converge after {itr + 1} iterations,"
  274. f" value is {p0}."
  275. )
  276. raise RuntimeError(msg)
  277. warnings.warn(msg, RuntimeWarning, stacklevel=2)
  278. return _results_select(
  279. full_output, (p0, funcalls, itr + 1, _ECONVERR), method)
  280. newton_step = fval / fder
  281. if fprime2:
  282. fder2 = fprime2(p0, *args)
  283. funcalls += 1
  284. method = "halley"
  285. # Halley's method:
  286. # newton_step /= (1.0 - 0.5 * newton_step * fder2 / fder)
  287. # Only do it if denominator stays close enough to 1
  288. # Rationale: If 1-adj < 0, then Halley sends x in the
  289. # opposite direction to Newton. Doesn't happen if x is close
  290. # enough to root.
  291. adj = newton_step * fder2 / fder / 2
  292. if np.abs(adj) < 1:
  293. newton_step /= 1.0 - adj
  294. p = p0 - newton_step
  295. if np.isclose(p, p0, rtol=rtol, atol=tol):
  296. return _results_select(
  297. full_output, (p, funcalls, itr + 1, _ECONVERGED), method)
  298. p0 = p
  299. else:
  300. # Secant method
  301. method = "secant"
  302. if x1 is not None:
  303. if x1 == x0:
  304. raise ValueError("x1 and x0 must be different")
  305. p1 = x1
  306. else:
  307. eps = 1e-4
  308. p1 = x0 * (1 + eps)
  309. p1 += (eps if p1 >= 0 else -eps)
  310. q0 = func(p0, *args)
  311. funcalls += 1
  312. q1 = func(p1, *args)
  313. funcalls += 1
  314. if abs(q1) < abs(q0):
  315. p0, p1, q0, q1 = p1, p0, q1, q0
  316. for itr in range(maxiter):
  317. if q1 == q0:
  318. if p1 != p0:
  319. msg = f"Tolerance of {p1 - p0} reached."
  320. if disp:
  321. msg += (
  322. f" Failed to converge after {itr + 1} iterations,"
  323. f" value is {p1}."
  324. )
  325. raise RuntimeError(msg)
  326. warnings.warn(msg, RuntimeWarning, stacklevel=2)
  327. p = (p1 + p0) / 2.0
  328. return _results_select(
  329. full_output, (p, funcalls, itr + 1, _ECONVERR), method)
  330. else:
  331. if abs(q1) > abs(q0):
  332. p = (-q0 / q1 * p1 + p0) / (1 - q0 / q1)
  333. else:
  334. p = (-q1 / q0 * p0 + p1) / (1 - q1 / q0)
  335. if np.isclose(p, p1, rtol=rtol, atol=tol):
  336. return _results_select(
  337. full_output, (p, funcalls, itr + 1, _ECONVERGED), method)
  338. p0, q0 = p1, q1
  339. p1 = p
  340. q1 = func(p1, *args)
  341. funcalls += 1
  342. if disp:
  343. msg = f"Failed to converge after {itr + 1} iterations, value is {p}."
  344. raise RuntimeError(msg)
  345. return _results_select(full_output, (p, funcalls, itr + 1, _ECONVERR), method)
  346. def _array_newton(func, x0, fprime, args, tol, maxiter, fprime2, full_output):
  347. """
  348. A vectorized version of Newton, Halley, and secant methods for arrays.
  349. Do not use this method directly. This method is called from `newton`
  350. when ``np.size(x0) > 1`` is ``True``. For docstring, see `newton`.
  351. """
  352. # Explicitly copy `x0` as `p` will be modified inplace, but the
  353. # user's array should not be altered.
  354. p = np.array(x0, copy=True)
  355. failures = np.ones_like(p, dtype=bool)
  356. nz_der = np.ones_like(failures)
  357. if fprime is not None:
  358. # Newton-Raphson method
  359. for iteration in range(maxiter):
  360. # first evaluate fval
  361. fval = np.asarray(func(p, *args))
  362. # If all fval are 0, all roots have been found, then terminate
  363. if not fval.any():
  364. failures = fval.astype(bool)
  365. break
  366. fder = np.asarray(fprime(p, *args))
  367. nz_der = (fder != 0)
  368. # stop iterating if all derivatives are zero
  369. if not nz_der.any():
  370. break
  371. # Newton step
  372. dp = fval[nz_der] / fder[nz_der]
  373. if fprime2 is not None:
  374. fder2 = np.asarray(fprime2(p, *args))
  375. dp = dp / (1.0 - 0.5 * dp * fder2[nz_der] / fder[nz_der])
  376. # only update nonzero derivatives
  377. p = np.asarray(p, dtype=np.result_type(p, dp, np.float64))
  378. p[nz_der] -= dp
  379. failures[nz_der] = np.abs(dp) >= tol # items not yet converged
  380. # stop iterating if there aren't any failures, not incl zero der
  381. if not failures[nz_der].any():
  382. break
  383. else:
  384. # Secant method
  385. dx = np.finfo(float).eps**0.33
  386. p1 = p * (1 + dx) + np.where(p >= 0, dx, -dx)
  387. q0 = np.asarray(func(p, *args))
  388. q1 = np.asarray(func(p1, *args))
  389. active = np.ones_like(p, dtype=bool)
  390. for iteration in range(maxiter):
  391. nz_der = (q1 != q0)
  392. # stop iterating if all derivatives are zero
  393. if not nz_der.any():
  394. p = (p1 + p) / 2.0
  395. break
  396. # Secant Step
  397. dp = (q1 * (p1 - p))[nz_der] / (q1 - q0)[nz_der]
  398. # only update nonzero derivatives
  399. p = np.asarray(p, dtype=np.result_type(p, p1, dp, np.float64))
  400. p[nz_der] = p1[nz_der] - dp
  401. active_zero_der = ~nz_der & active
  402. p[active_zero_der] = (p1 + p)[active_zero_der] / 2.0
  403. active &= nz_der # don't assign zero derivatives again
  404. failures[nz_der] = np.abs(dp) >= tol # not yet converged
  405. # stop iterating if there aren't any failures, not incl zero der
  406. if not failures[nz_der].any():
  407. break
  408. p1, p = p, p1
  409. q0 = q1
  410. q1 = np.asarray(func(p1, *args))
  411. zero_der = ~nz_der & failures # don't include converged with zero-ders
  412. if zero_der.any():
  413. # Secant warnings
  414. if fprime is None:
  415. nonzero_dp = (p1 != p)
  416. # non-zero dp, but infinite newton step
  417. zero_der_nz_dp = (zero_der & nonzero_dp)
  418. if zero_der_nz_dp.any():
  419. rms = np.sqrt(
  420. sum((p1[zero_der_nz_dp] - p[zero_der_nz_dp]) ** 2)
  421. )
  422. warnings.warn(f'RMS of {rms:g} reached', RuntimeWarning, stacklevel=3)
  423. # Newton or Halley warnings
  424. else:
  425. all_or_some = 'all' if zero_der.all() else 'some'
  426. msg = f'{all_or_some:s} derivatives were zero'
  427. warnings.warn(msg, RuntimeWarning, stacklevel=3)
  428. elif failures.any():
  429. all_or_some = 'all' if failures.all() else 'some'
  430. msg = f'{all_or_some:s} failed to converge after {maxiter:d} iterations'
  431. if failures.all():
  432. raise RuntimeError(msg)
  433. warnings.warn(msg, RuntimeWarning, stacklevel=3)
  434. if full_output:
  435. result = namedtuple('result', ('root', 'converged', 'zero_der'))
  436. p = result(p, ~failures, zero_der)
  437. return p
  438. def bisect(f, a, b, args=(),
  439. xtol=_xtol, rtol=_rtol, maxiter=_iter,
  440. full_output=False, disp=True):
  441. """
  442. Find root of a function within an interval using bisection.
  443. Basic bisection routine to find a root of the function `f` between the
  444. arguments `a` and `b`. `f(a)` and `f(b)` cannot have the same signs.
  445. Slow but sure.
  446. Parameters
  447. ----------
  448. f : function
  449. Python function returning a number. `f` must be continuous, and
  450. f(a) and f(b) must have opposite signs.
  451. a : scalar
  452. One end of the bracketing interval [a,b].
  453. b : scalar
  454. The other end of the bracketing interval [a,b].
  455. xtol : number, optional
  456. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  457. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  458. parameter must be positive.
  459. rtol : number, optional
  460. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  461. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  462. parameter cannot be smaller than its default value of
  463. ``4*np.finfo(float).eps``.
  464. maxiter : int, optional
  465. If convergence is not achieved in `maxiter` iterations, an error is
  466. raised. Must be >= 0.
  467. args : tuple, optional
  468. Containing extra arguments for the function `f`.
  469. `f` is called by ``apply(f, (x)+args)``.
  470. full_output : bool, optional
  471. If `full_output` is False, the root is returned. If `full_output` is
  472. True, the return value is ``(x, r)``, where x is the root, and r is
  473. a `RootResults` object.
  474. disp : bool, optional
  475. If True, raise RuntimeError if the algorithm didn't converge.
  476. Otherwise, the convergence status is recorded in a `RootResults`
  477. return object.
  478. Returns
  479. -------
  480. root : float
  481. Root of `f` between `a` and `b`.
  482. r : `RootResults` (present if ``full_output = True``)
  483. Object containing information about the convergence. In particular,
  484. ``r.converged`` is True if the routine converged.
  485. Notes
  486. -----
  487. As mentioned in the parameter documentation, the computed root ``x0`` will
  488. satisfy ``np.isclose(x, x0, atol=xtol, rtol=rtol)``, where ``x`` is the
  489. exact root. In equation form, this terminating condition is ``abs(x - x0)
  490. <= xtol + rtol * abs(x0)``.
  491. The default value ``xtol=2e-12`` may lead to surprising behavior if one
  492. expects `bisect` to always compute roots with relative error near machine
  493. precision. Care should be taken to select `xtol` for the use case at hand.
  494. Setting ``xtol=5e-324``, the smallest subnormal number, will ensure the
  495. highest level of accuracy. Larger values of `xtol` may be useful for saving
  496. function evaluations when a root is at or near zero in applications where
  497. the tiny absolute differences available between floating point numbers near
  498. zero are not meaningful.
  499. Examples
  500. --------
  501. >>> def f(x):
  502. ... return (x**2 - 1)
  503. >>> from scipy import optimize
  504. >>> root = optimize.bisect(f, 0, 2)
  505. >>> root
  506. 1.0
  507. >>> root = optimize.bisect(f, -2, 0)
  508. >>> root
  509. -1.0
  510. See Also
  511. --------
  512. brentq, brenth, bisect, newton
  513. fixed_point : scalar fixed-point finder
  514. fsolve : n-dimensional root-finding
  515. elementwise.find_root : efficient elementwise 1-D root-finder
  516. """
  517. if not isinstance(args, tuple):
  518. args = (args,)
  519. maxiter = operator.index(maxiter)
  520. if xtol <= 0:
  521. raise ValueError(f"xtol too small ({xtol:g} <= 0)")
  522. if rtol < _rtol:
  523. raise ValueError(f"rtol too small ({rtol:g} < {_rtol:g})")
  524. f = _wrap_nan_raise(f)
  525. r = _zeros._bisect(f, a, b, xtol, rtol, maxiter, args, full_output, disp)
  526. return results_c(full_output, r, "bisect")
  527. def ridder(f, a, b, args=(),
  528. xtol=_xtol, rtol=_rtol, maxiter=_iter,
  529. full_output=False, disp=True):
  530. """
  531. Find a root of a function in an interval using Ridder's method.
  532. Parameters
  533. ----------
  534. f : function
  535. Python function returning a number. f must be continuous, and f(a) and
  536. f(b) must have opposite signs.
  537. a : scalar
  538. One end of the bracketing interval [a,b].
  539. b : scalar
  540. The other end of the bracketing interval [a,b].
  541. xtol : number, optional
  542. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  543. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  544. parameter must be positive.
  545. rtol : number, optional
  546. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  547. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  548. parameter cannot be smaller than its default value of
  549. ``4*np.finfo(float).eps``.
  550. maxiter : int, optional
  551. If convergence is not achieved in `maxiter` iterations, an error is
  552. raised. Must be >= 0.
  553. args : tuple, optional
  554. Containing extra arguments for the function `f`.
  555. `f` is called by ``apply(f, (x)+args)``.
  556. full_output : bool, optional
  557. If `full_output` is False, the root is returned. If `full_output` is
  558. True, the return value is ``(x, r)``, where `x` is the root, and `r` is
  559. a `RootResults` object.
  560. disp : bool, optional
  561. If True, raise RuntimeError if the algorithm didn't converge.
  562. Otherwise, the convergence status is recorded in any `RootResults`
  563. return object.
  564. Returns
  565. -------
  566. root : float
  567. Root of `f` between `a` and `b`.
  568. r : `RootResults` (present if ``full_output = True``)
  569. Object containing information about the convergence.
  570. In particular, ``r.converged`` is True if the routine converged.
  571. See Also
  572. --------
  573. brentq, brenth, bisect, newton : 1-D root-finding
  574. fixed_point : scalar fixed-point finder
  575. elementwise.find_root : efficient elementwise 1-D root-finder
  576. Notes
  577. -----
  578. Uses [Ridders1979]_ method to find a root of the function `f` between the
  579. arguments `a` and `b`. Ridders' method is faster than bisection, but not
  580. generally as fast as the Brent routines. [Ridders1979]_ provides the
  581. classic description and source of the algorithm. A description can also be
  582. found in any recent edition of Numerical Recipes.
  583. The routine used here diverges slightly from standard presentations in
  584. order to be a bit more careful of tolerance.
  585. As mentioned in the parameter documentation, the computed root ``x0`` will
  586. satisfy ``np.isclose(x, x0, atol=xtol, rtol=rtol)``, where ``x`` is the
  587. exact root. In equation form, this terminating condition is ``abs(x - x0)
  588. <= xtol + rtol * abs(x0)``.
  589. The default value ``xtol=2e-12`` may lead to surprising behavior if one
  590. expects `ridder` to always compute roots with relative error near machine
  591. precision. Care should be taken to select `xtol` for the use case at hand.
  592. Setting ``xtol=5e-324``, the smallest subnormal number, will ensure the
  593. highest level of accuracy. Larger values of `xtol` may be useful for saving
  594. function evaluations when a root is at or near zero in applications where
  595. the tiny absolute differences available between floating point numbers near
  596. zero are not meaningful.
  597. References
  598. ----------
  599. .. [Ridders1979]
  600. Ridders, C. F. J. "A New Algorithm for Computing a
  601. Single Root of a Real Continuous Function."
  602. IEEE Trans. Circuits Systems 26, 979-980, 1979.
  603. Examples
  604. --------
  605. >>> def f(x):
  606. ... return (x**2 - 1)
  607. >>> from scipy import optimize
  608. >>> root = optimize.ridder(f, 0, 2)
  609. >>> root
  610. 1.0
  611. >>> root = optimize.ridder(f, -2, 0)
  612. >>> root
  613. -1.0
  614. """
  615. if not isinstance(args, tuple):
  616. args = (args,)
  617. maxiter = operator.index(maxiter)
  618. if xtol <= 0:
  619. raise ValueError(f"xtol too small ({xtol:g} <= 0)")
  620. if rtol < _rtol:
  621. raise ValueError(f"rtol too small ({rtol:g} < {_rtol:g})")
  622. f = _wrap_nan_raise(f)
  623. r = _zeros._ridder(f, a, b, xtol, rtol, maxiter, args, full_output, disp)
  624. return results_c(full_output, r, "ridder")
  625. def brentq(f, a, b, args=(),
  626. xtol=_xtol, rtol=_rtol, maxiter=_iter,
  627. full_output=False, disp=True):
  628. """
  629. Find a root of a function in a bracketing interval using Brent's method.
  630. Uses the classic Brent's method to find a root of the function `f` on
  631. the sign changing interval [a , b]. Generally considered the best of the
  632. rootfinding routines here. It is a safe version of the secant method that
  633. uses inverse quadratic extrapolation. Brent's method combines root
  634. bracketing, interval bisection, and inverse quadratic interpolation. It is
  635. sometimes known as the van Wijngaarden-Dekker-Brent method. Brent (1973)
  636. claims convergence is guaranteed for functions computable within [a,b].
  637. [Brent1973]_ provides the classic description of the algorithm. Another
  638. description can be found in a recent edition of Numerical Recipes, including
  639. [PressEtal1992]_. A third description is at
  640. http://mathworld.wolfram.com/BrentsMethod.html. It should be easy to
  641. understand the algorithm just by reading our code. Our code diverges a bit
  642. from standard presentations: we choose a different formula for the
  643. extrapolation step.
  644. Parameters
  645. ----------
  646. f : function
  647. Python function returning a number. The function :math:`f`
  648. must be continuous, and :math:`f(a)` and :math:`f(b)` must
  649. have opposite signs.
  650. a : scalar
  651. One end of the bracketing interval :math:`[a, b]`.
  652. b : scalar
  653. The other end of the bracketing interval :math:`[a, b]`.
  654. xtol : number, optional
  655. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  656. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  657. parameter must be positive. For nice functions, Brent's
  658. method will often satisfy the above condition with ``xtol/2``
  659. and ``rtol/2``. [Brent1973]_
  660. rtol : number, optional
  661. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  662. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  663. parameter cannot be smaller than its default value of
  664. ``4*np.finfo(float).eps``. For nice functions, Brent's
  665. method will often satisfy the above condition with ``xtol/2``
  666. and ``rtol/2``. [Brent1973]_
  667. maxiter : int, optional
  668. If convergence is not achieved in `maxiter` iterations, an error is
  669. raised. Must be >= 0.
  670. args : tuple, optional
  671. Containing extra arguments for the function `f`.
  672. `f` is called by ``apply(f, (x)+args)``.
  673. full_output : bool, optional
  674. If `full_output` is False, the root is returned. If `full_output` is
  675. True, the return value is ``(x, r)``, where `x` is the root, and `r` is
  676. a `RootResults` object.
  677. disp : bool, optional
  678. If True, raise RuntimeError if the algorithm didn't converge.
  679. Otherwise, the convergence status is recorded in any `RootResults`
  680. return object.
  681. Returns
  682. -------
  683. root : float
  684. Root of `f` between `a` and `b`.
  685. r : `RootResults` (present if ``full_output = True``)
  686. Object containing information about the convergence. In particular,
  687. ``r.converged`` is True if the routine converged.
  688. See Also
  689. --------
  690. fmin, fmin_powell, fmin_cg, fmin_bfgs, fmin_ncg : multivariate local optimizers
  691. leastsq : nonlinear least squares minimizer
  692. fmin_l_bfgs_b, fmin_tnc, fmin_cobyla : constrained multivariate optimizers
  693. basinhopping, differential_evolution, brute : global optimizers
  694. fminbound, brent, golden, bracket : local scalar minimizers
  695. fsolve : N-D root-finding
  696. brenth, ridder, bisect, newton : 1-D root-finding
  697. fixed_point : scalar fixed-point finder
  698. elementwise.find_root : efficient elementwise 1-D root-finder
  699. Notes
  700. -----
  701. `f` must be continuous. f(a) and f(b) must have opposite signs.
  702. As mentioned in the parameter documentation, the computed root ``x0`` will
  703. satisfy ``np.isclose(x, x0, atol=xtol, rtol=rtol)``, where ``x`` is the
  704. exact root. In equation form, this terminating condition is ``abs(x - x0)
  705. <= xtol + rtol * abs(x0)``.
  706. The default value ``xtol=2e-12`` may lead to surprising behavior if one
  707. expects `brentq` to always compute roots with relative error near machine
  708. precision. Care should be taken to select `xtol` for the use case at hand.
  709. Setting ``xtol=5e-324``, the smallest subnormal number, will ensure the
  710. highest level of accuracy. Larger values of `xtol` may be useful for saving
  711. function evaluations when a root is at or near zero in applications where
  712. the tiny absolute differences available between floating point numbers near
  713. zero are not meaningful.
  714. References
  715. ----------
  716. .. [Brent1973]
  717. Brent, R. P.,
  718. *Algorithms for Minimization Without Derivatives*.
  719. Englewood Cliffs, NJ: Prentice-Hall, 1973. Ch. 3-4.
  720. .. [PressEtal1992]
  721. Press, W. H.; Flannery, B. P.; Teukolsky, S. A.; and Vetterling, W. T.
  722. *Numerical Recipes in FORTRAN: The Art of Scientific Computing*, 2nd ed.
  723. Cambridge, England: Cambridge University Press, pp. 352-355, 1992.
  724. Section 9.3: "Van Wijngaarden-Dekker-Brent Method."
  725. Examples
  726. --------
  727. >>> def f(x):
  728. ... return (x**2 - 1)
  729. >>> from scipy import optimize
  730. >>> root = optimize.brentq(f, -2, 0)
  731. >>> root
  732. -1.0
  733. >>> root = optimize.brentq(f, 0, 2)
  734. >>> root
  735. 1.0
  736. """
  737. if not isinstance(args, tuple):
  738. args = (args,)
  739. maxiter = operator.index(maxiter)
  740. if xtol <= 0:
  741. raise ValueError(f"xtol too small ({xtol:g} <= 0)")
  742. if rtol < _rtol:
  743. raise ValueError(f"rtol too small ({rtol:g} < {_rtol:g})")
  744. f = _wrap_nan_raise(f)
  745. r = _zeros._brentq(f, a, b, xtol, rtol, maxiter, args, full_output, disp)
  746. return results_c(full_output, r, "brentq")
  747. def brenth(f, a, b, args=(),
  748. xtol=_xtol, rtol=_rtol, maxiter=_iter,
  749. full_output=False, disp=True):
  750. """Find a root of a function in a bracketing interval using Brent's
  751. method with hyperbolic extrapolation.
  752. A variation on the classic Brent routine to find a root of the function f
  753. between the arguments a and b that uses hyperbolic extrapolation instead of
  754. inverse quadratic extrapolation. Bus & Dekker (1975) guarantee convergence
  755. for this method, claiming that the upper bound of function evaluations here
  756. is 4 or 5 times that of bisection.
  757. f(a) and f(b) cannot have the same signs. Generally, on a par with the
  758. brent routine, but not as heavily tested. It is a safe version of the
  759. secant method that uses hyperbolic extrapolation.
  760. The version here is by Chuck Harris, and implements Algorithm M of
  761. [BusAndDekker1975]_, where further details (convergence properties,
  762. additional remarks and such) can be found
  763. Parameters
  764. ----------
  765. f : function
  766. Python function returning a number. f must be continuous, and f(a) and
  767. f(b) must have opposite signs.
  768. a : scalar
  769. One end of the bracketing interval [a,b].
  770. b : scalar
  771. The other end of the bracketing interval [a,b].
  772. xtol : number, optional
  773. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  774. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  775. parameter must be positive. As with `brentq`, for nice
  776. functions the method will often satisfy the above condition
  777. with ``xtol/2`` and ``rtol/2``.
  778. rtol : number, optional
  779. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  780. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  781. parameter cannot be smaller than its default value of
  782. ``4*np.finfo(float).eps``. As with `brentq`, for nice functions
  783. the method will often satisfy the above condition with
  784. ``xtol/2`` and ``rtol/2``.
  785. maxiter : int, optional
  786. If convergence is not achieved in `maxiter` iterations, an error is
  787. raised. Must be >= 0.
  788. args : tuple, optional
  789. Containing extra arguments for the function `f`.
  790. `f` is called by ``apply(f, (x)+args)``.
  791. full_output : bool, optional
  792. If `full_output` is False, the root is returned. If `full_output` is
  793. True, the return value is ``(x, r)``, where `x` is the root, and `r` is
  794. a `RootResults` object.
  795. disp : bool, optional
  796. If True, raise RuntimeError if the algorithm didn't converge.
  797. Otherwise, the convergence status is recorded in any `RootResults`
  798. return object.
  799. Returns
  800. -------
  801. root : float
  802. Root of `f` between `a` and `b`.
  803. r : `RootResults` (present if ``full_output = True``)
  804. Object containing information about the convergence. In particular,
  805. ``r.converged`` is True if the routine converged.
  806. See Also
  807. --------
  808. fmin, fmin_powell, fmin_cg, fmin_bfgs, fmin_ncg : multivariate local optimizers
  809. leastsq : nonlinear least squares minimizer
  810. fmin_l_bfgs_b, fmin_tnc, fmin_cobyla : constrained multivariate optimizers
  811. basinhopping, differential_evolution, brute : global optimizers
  812. fminbound, brent, golden, bracket : local scalar minimizers
  813. fsolve : N-D root-finding
  814. brentq, ridder, bisect, newton : 1-D root-finding
  815. fixed_point : scalar fixed-point finder
  816. elementwise.find_root : efficient elementwise 1-D root-finder
  817. Notes
  818. -----
  819. As mentioned in the parameter documentation, the computed root ``x0`` will
  820. satisfy ``np.isclose(x, x0, atol=xtol, rtol=rtol)``, where ``x`` is the
  821. exact root. In equation form, this terminating condition is ``abs(x - x0)
  822. <= xtol + rtol * abs(x0)``.
  823. The default value ``xtol=2e-12`` may lead to surprising behavior if one
  824. expects `brenth` to always compute roots with relative error near machine
  825. precision. Care should be taken to select `xtol` for the use case at hand.
  826. Setting ``xtol=5e-324``, the smallest subnormal number, will ensure the
  827. highest level of accuracy. Larger values of `xtol` may be useful for saving
  828. function evaluations when a root is at or near zero in applications where
  829. the tiny absolute differences available between floating point numbers near
  830. zero are not meaningful.
  831. References
  832. ----------
  833. .. [BusAndDekker1975]
  834. Bus, J. C. P., Dekker, T. J.,
  835. "Two Efficient Algorithms with Guaranteed Convergence for Finding a Zero
  836. of a Function", ACM Transactions on Mathematical Software, Vol. 1, Issue
  837. 4, Dec. 1975, pp. 330-345. Section 3: "Algorithm M".
  838. :doi:`10.1145/355656.355659`
  839. Examples
  840. --------
  841. >>> def f(x):
  842. ... return (x**2 - 1)
  843. >>> from scipy import optimize
  844. >>> root = optimize.brenth(f, -2, 0)
  845. >>> root
  846. -1.0
  847. >>> root = optimize.brenth(f, 0, 2)
  848. >>> root
  849. 1.0
  850. """
  851. if not isinstance(args, tuple):
  852. args = (args,)
  853. maxiter = operator.index(maxiter)
  854. if xtol <= 0:
  855. raise ValueError(f"xtol too small ({xtol:g} <= 0)")
  856. if rtol < _rtol:
  857. raise ValueError(f"rtol too small ({rtol:g} < {_rtol:g})")
  858. f = _wrap_nan_raise(f)
  859. r = _zeros._brenth(f, a, b, xtol, rtol, maxiter, args, full_output, disp)
  860. return results_c(full_output, r, "brenth")
  861. ################################
  862. # TOMS "Algorithm 748: Enclosing Zeros of Continuous Functions", by
  863. # Alefeld, G. E. and Potra, F. A. and Shi, Yixun,
  864. # See [1]
  865. def _notclose(fs, rtol=_rtol, atol=_xtol):
  866. # Ensure not None, not 0, all finite, and not very close to each other
  867. notclosefvals = (
  868. all(fs) and all(np.isfinite(fs)) and
  869. not any(any(np.isclose(_f, fs[i + 1:], rtol=rtol, atol=atol))
  870. for i, _f in enumerate(fs[:-1])))
  871. return notclosefvals
  872. def _secant(xvals, fvals):
  873. """Perform a secant step, taking a little care"""
  874. # Secant has many "mathematically" equivalent formulations
  875. # x2 = x0 - (x1 - x0)/(f1 - f0) * f0
  876. # = x1 - (x1 - x0)/(f1 - f0) * f1
  877. # = (-x1 * f0 + x0 * f1) / (f1 - f0)
  878. # = (-f0 / f1 * x1 + x0) / (1 - f0 / f1)
  879. # = (-f1 / f0 * x0 + x1) / (1 - f1 / f0)
  880. x0, x1 = xvals[:2]
  881. f0, f1 = fvals[:2]
  882. if f0 == f1:
  883. return np.nan
  884. if np.abs(f1) > np.abs(f0):
  885. x2 = (-f0 / f1 * x1 + x0) / (1 - f0 / f1)
  886. else:
  887. x2 = (-f1 / f0 * x0 + x1) / (1 - f1 / f0)
  888. return x2
  889. def _update_bracket(ab, fab, c, fc):
  890. """Update a bracket given (c, fc), return the discarded endpoints."""
  891. fa, fb = fab
  892. idx = (0 if np.sign(fa) * np.sign(fc) > 0 else 1)
  893. rx, rfx = ab[idx], fab[idx]
  894. fab[idx] = fc
  895. ab[idx] = c
  896. return rx, rfx
  897. def _compute_divided_differences(xvals, fvals, N=None, full=True,
  898. forward=True):
  899. """Return a matrix of divided differences for the xvals, fvals pairs
  900. DD[i, j] = f[x_{i-j}, ..., x_i] for 0 <= j <= i
  901. If full is False, just return the main diagonal(or last row):
  902. f[a], f[a, b] and f[a, b, c].
  903. If forward is False, return f[c], f[b, c], f[a, b, c]."""
  904. if full:
  905. if forward:
  906. xvals = np.asarray(xvals)
  907. else:
  908. xvals = np.array(xvals)[::-1]
  909. M = len(xvals)
  910. N = M if N is None else min(N, M)
  911. DD = np.zeros([M, N])
  912. DD[:, 0] = fvals[:]
  913. for i in range(1, N):
  914. DD[i:, i] = (np.diff(DD[i - 1:, i - 1]) /
  915. (xvals[i:] - xvals[:M - i]))
  916. return DD
  917. xvals = np.asarray(xvals)
  918. dd = np.array(fvals)
  919. row = np.array(fvals)
  920. idx2Use = (0 if forward else -1)
  921. dd[0] = fvals[idx2Use]
  922. for i in range(1, len(xvals)):
  923. denom = xvals[i:i + len(row) - 1] - xvals[:len(row) - 1]
  924. row = np.diff(row)[:] / denom
  925. dd[i] = row[idx2Use]
  926. return dd
  927. def _interpolated_poly(xvals, fvals, x):
  928. """Compute p(x) for the polynomial passing through the specified locations.
  929. Use Neville's algorithm to compute p(x) where p is the minimal degree
  930. polynomial passing through the points xvals, fvals"""
  931. xvals = np.asarray(xvals)
  932. N = len(xvals)
  933. Q = np.zeros([N, N])
  934. D = np.zeros([N, N])
  935. Q[:, 0] = fvals[:]
  936. D[:, 0] = fvals[:]
  937. for k in range(1, N):
  938. alpha = D[k:, k - 1] - Q[k - 1:N - 1, k - 1]
  939. diffik = xvals[0:N - k] - xvals[k:N]
  940. Q[k:, k] = (xvals[k:] - x) / diffik * alpha
  941. D[k:, k] = (xvals[:N - k] - x) / diffik * alpha
  942. # Expect Q[-1, 1:] to be small relative to Q[-1, 0] as x approaches a root
  943. return np.sum(Q[-1, 1:]) + Q[-1, 0]
  944. def _inverse_poly_zero(a, b, c, d, fa, fb, fc, fd):
  945. """Inverse cubic interpolation f-values -> x-values
  946. Given four points (fa, a), (fb, b), (fc, c), (fd, d) with
  947. fa, fb, fc, fd all distinct, find poly IP(y) through the 4 points
  948. and compute x=IP(0).
  949. """
  950. return _interpolated_poly([fa, fb, fc, fd], [a, b, c, d], 0)
  951. def _newton_quadratic(ab, fab, d, fd, k):
  952. """Apply Newton-Raphson like steps, using divided differences to approximate f'
  953. ab is a real interval [a, b] containing a root,
  954. fab holds the real values of f(a), f(b)
  955. d is a real number outside [ab, b]
  956. k is the number of steps to apply
  957. """
  958. a, b = ab
  959. fa, fb = fab
  960. _, B, A = _compute_divided_differences([a, b, d], [fa, fb, fd],
  961. forward=True, full=False)
  962. # _P is the quadratic polynomial through the 3 points
  963. def _P(x):
  964. # Horner evaluation of fa + B * (x - a) + A * (x - a) * (x - b)
  965. return (A * (x - b) + B) * (x - a) + fa
  966. if A == 0:
  967. r = a - fa / B
  968. else:
  969. r = (a if np.sign(A) * np.sign(fa) > 0 else b)
  970. # Apply k Newton-Raphson steps to _P(x), starting from x=r
  971. for i in range(k):
  972. r1 = r - _P(r) / (B + A * (2 * r - a - b))
  973. if not (ab[0] < r1 < ab[1]):
  974. if (ab[0] < r < ab[1]):
  975. return r
  976. r = sum(ab) / 2.0
  977. break
  978. r = r1
  979. return r
  980. class TOMS748Solver:
  981. """Solve f(x, *args) == 0 using Algorithm748 of Alefeld, Potro & Shi.
  982. """
  983. _MU = 0.5
  984. _K_MIN = 1
  985. _K_MAX = 100 # A very high value for real usage. Expect 1, 2, maybe 3.
  986. def __init__(self):
  987. self.f = None
  988. self.args = None
  989. self.function_calls = 0
  990. self.iterations = 0
  991. self.k = 2
  992. # ab=[a,b] is a global interval containing a root
  993. self.ab = [np.nan, np.nan]
  994. # fab is function values at a, b
  995. self.fab = [np.nan, np.nan]
  996. self.d = None
  997. self.fd = None
  998. self.e = None
  999. self.fe = None
  1000. self.disp = False
  1001. self.xtol = _xtol
  1002. self.rtol = _rtol
  1003. self.maxiter = _iter
  1004. def configure(self, xtol, rtol, maxiter, disp, k):
  1005. self.disp = disp
  1006. self.xtol = xtol
  1007. self.rtol = rtol
  1008. self.maxiter = maxiter
  1009. # Silently replace a low value of k with 1
  1010. self.k = max(k, self._K_MIN)
  1011. # Noisily replace a high value of k with self._K_MAX
  1012. if self.k > self._K_MAX:
  1013. msg = f"toms748: Overriding k: ->{self._K_MAX}"
  1014. warnings.warn(msg, RuntimeWarning, stacklevel=3)
  1015. self.k = self._K_MAX
  1016. def _callf(self, x, error=True):
  1017. """Call the user-supplied function, update book-keeping"""
  1018. fx = self.f(x, *self.args)
  1019. self.function_calls += 1
  1020. if not np.isfinite(fx) and error:
  1021. raise ValueError(f"Invalid function value: f({x:f}) -> {fx} ")
  1022. return fx
  1023. def get_result(self, x, flag=_ECONVERGED):
  1024. r"""Package the result and statistics into a tuple."""
  1025. return (x, self.function_calls, self.iterations, flag)
  1026. def _update_bracket(self, c, fc):
  1027. return _update_bracket(self.ab, self.fab, c, fc)
  1028. def start(self, f, a, b, args=()):
  1029. r"""Prepare for the iterations."""
  1030. self.function_calls = 0
  1031. self.iterations = 0
  1032. self.f = f
  1033. self.args = args
  1034. self.ab[:] = [a, b]
  1035. if not np.isfinite(a) or np.imag(a) != 0:
  1036. raise ValueError(f"Invalid x value: {a} ")
  1037. if not np.isfinite(b) or np.imag(b) != 0:
  1038. raise ValueError(f"Invalid x value: {b} ")
  1039. fa = self._callf(a)
  1040. if not np.isfinite(fa) or np.imag(fa) != 0:
  1041. raise ValueError(f"Invalid function value: f({a:f}) -> {fa} ")
  1042. if fa == 0:
  1043. return _ECONVERGED, a
  1044. fb = self._callf(b)
  1045. if not np.isfinite(fb) or np.imag(fb) != 0:
  1046. raise ValueError(f"Invalid function value: f({b:f}) -> {fb} ")
  1047. if fb == 0:
  1048. return _ECONVERGED, b
  1049. if np.sign(fb) * np.sign(fa) > 0:
  1050. raise ValueError("f(a) and f(b) must have different signs, but "
  1051. f"f({a:e})={fa:e}, f({b:e})={fb:e} ")
  1052. self.fab[:] = [fa, fb]
  1053. return _EINPROGRESS, sum(self.ab) / 2.0
  1054. def get_status(self):
  1055. """Determine the current status."""
  1056. a, b = self.ab[:2]
  1057. if np.isclose(a, b, rtol=self.rtol, atol=self.xtol):
  1058. return _ECONVERGED, sum(self.ab) / 2.0
  1059. if self.iterations >= self.maxiter:
  1060. return _ECONVERR, sum(self.ab) / 2.0
  1061. return _EINPROGRESS, sum(self.ab) / 2.0
  1062. def iterate(self):
  1063. """Perform one step in the algorithm.
  1064. Implements Algorithm 4.1(k=1) or 4.2(k=2) in [APS1995]
  1065. """
  1066. self.iterations += 1
  1067. eps = np.finfo(float).eps
  1068. d, fd, e, fe = self.d, self.fd, self.e, self.fe
  1069. ab_width = self.ab[1] - self.ab[0] # Need the start width below
  1070. c = None
  1071. for nsteps in range(2, self.k+2):
  1072. # If the f-values are sufficiently separated, perform an inverse
  1073. # polynomial interpolation step. Otherwise, nsteps repeats of
  1074. # an approximate Newton-Raphson step.
  1075. if _notclose(self.fab + [fd, fe], rtol=0, atol=32*eps):
  1076. c0 = _inverse_poly_zero(self.ab[0], self.ab[1], d, e,
  1077. self.fab[0], self.fab[1], fd, fe)
  1078. if self.ab[0] < c0 < self.ab[1]:
  1079. c = c0
  1080. if c is None:
  1081. c = _newton_quadratic(self.ab, self.fab, d, fd, nsteps)
  1082. fc = self._callf(c)
  1083. if fc == 0:
  1084. return _ECONVERGED, c
  1085. # re-bracket
  1086. e, fe = d, fd
  1087. d, fd = self._update_bracket(c, fc)
  1088. # u is the endpoint with the smallest f-value
  1089. uix = (0 if np.abs(self.fab[0]) < np.abs(self.fab[1]) else 1)
  1090. u, fu = self.ab[uix], self.fab[uix]
  1091. _, A = _compute_divided_differences(self.ab, self.fab,
  1092. forward=(uix == 0), full=False)
  1093. c = u - 2 * fu / A
  1094. if np.abs(c - u) > 0.5 * (self.ab[1] - self.ab[0]):
  1095. c = sum(self.ab) / 2.0
  1096. else:
  1097. if np.isclose(c, u, rtol=eps, atol=0):
  1098. # c didn't change (much).
  1099. # Either because the f-values at the endpoints have vastly
  1100. # differing magnitudes, or because the root is very close to
  1101. # that endpoint
  1102. frs = np.frexp(self.fab)[1]
  1103. if frs[uix] < frs[1 - uix] - 50: # Differ by more than 2**50
  1104. c = (31 * self.ab[uix] + self.ab[1 - uix]) / 32
  1105. else:
  1106. # Make a bigger adjustment, about the
  1107. # size of the requested tolerance.
  1108. mm = (1 if uix == 0 else -1)
  1109. adj = mm * np.abs(c) * self.rtol + mm * self.xtol
  1110. c = u + adj
  1111. if not self.ab[0] < c < self.ab[1]:
  1112. c = sum(self.ab) / 2.0
  1113. fc = self._callf(c)
  1114. if fc == 0:
  1115. return _ECONVERGED, c
  1116. e, fe = d, fd
  1117. d, fd = self._update_bracket(c, fc)
  1118. # If the width of the new interval did not decrease enough, bisect
  1119. if self.ab[1] - self.ab[0] > self._MU * ab_width:
  1120. e, fe = d, fd
  1121. z = sum(self.ab) / 2.0
  1122. fz = self._callf(z)
  1123. if fz == 0:
  1124. return _ECONVERGED, z
  1125. d, fd = self._update_bracket(z, fz)
  1126. # Record d and e for next iteration
  1127. self.d, self.fd = d, fd
  1128. self.e, self.fe = e, fe
  1129. status, xn = self.get_status()
  1130. return status, xn
  1131. def solve(self, f, a, b, args=(),
  1132. xtol=_xtol, rtol=_rtol, k=2, maxiter=_iter, disp=True):
  1133. r"""Solve f(x) = 0 given an interval containing a root."""
  1134. self.configure(xtol=xtol, rtol=rtol, maxiter=maxiter, disp=disp, k=k)
  1135. status, xn = self.start(f, a, b, args)
  1136. if status == _ECONVERGED:
  1137. return self.get_result(xn)
  1138. # The first step only has two x-values.
  1139. c = _secant(self.ab, self.fab)
  1140. if not self.ab[0] < c < self.ab[1]:
  1141. c = sum(self.ab) / 2.0
  1142. fc = self._callf(c)
  1143. if fc == 0:
  1144. return self.get_result(c)
  1145. self.d, self.fd = self._update_bracket(c, fc)
  1146. self.e, self.fe = None, None
  1147. self.iterations += 1
  1148. while True:
  1149. status, xn = self.iterate()
  1150. if status == _ECONVERGED:
  1151. return self.get_result(xn)
  1152. if status == _ECONVERR:
  1153. fmt = "Failed to converge after %d iterations, bracket is %s"
  1154. if disp:
  1155. msg = fmt % (self.iterations + 1, self.ab)
  1156. raise RuntimeError(msg)
  1157. return self.get_result(xn, _ECONVERR)
  1158. def toms748(f, a, b, args=(), k=1,
  1159. xtol=_xtol, rtol=_rtol, maxiter=_iter,
  1160. full_output=False, disp=True):
  1161. """
  1162. Find a root using TOMS Algorithm 748 method.
  1163. Implements the Algorithm 748 method of Alefeld, Potro and Shi to find a
  1164. root of the function `f` on the interval ``[a , b]``, where `f(a)` and
  1165. `f(b)` must have opposite signs.
  1166. It uses a mixture of inverse cubic interpolation and
  1167. "Newton-quadratic" steps. [APS1995].
  1168. Parameters
  1169. ----------
  1170. f : function
  1171. Python function returning a scalar. The function :math:`f`
  1172. must be continuous, and :math:`f(a)` and :math:`f(b)`
  1173. have opposite signs.
  1174. a : scalar,
  1175. lower boundary of the search interval
  1176. b : scalar,
  1177. upper boundary of the search interval
  1178. args : tuple, optional
  1179. containing extra arguments for the function `f`.
  1180. `f` is called by ``f(x, *args)``.
  1181. k : int, optional
  1182. The number of Newton quadratic steps to perform each
  1183. iteration. ``k>=1``.
  1184. xtol : scalar, optional
  1185. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  1186. atol=xtol, rtol=rtol)``, where ``x`` is the exact root. The
  1187. parameter must be positive.
  1188. rtol : scalar, optional
  1189. The computed root ``x0`` will satisfy ``np.isclose(x, x0,
  1190. atol=xtol, rtol=rtol)``, where ``x`` is the exact root.
  1191. maxiter : int, optional
  1192. If convergence is not achieved in `maxiter` iterations, an error is
  1193. raised. Must be >= 0.
  1194. full_output : bool, optional
  1195. If `full_output` is False, the root is returned. If `full_output` is
  1196. True, the return value is ``(x, r)``, where `x` is the root, and `r` is
  1197. a `RootResults` object.
  1198. disp : bool, optional
  1199. If True, raise RuntimeError if the algorithm didn't converge.
  1200. Otherwise, the convergence status is recorded in the `RootResults`
  1201. return object.
  1202. Returns
  1203. -------
  1204. root : float
  1205. Approximate root of `f`
  1206. r : `RootResults` (present if ``full_output = True``)
  1207. Object containing information about the convergence. In particular,
  1208. ``r.converged`` is True if the routine converged.
  1209. See Also
  1210. --------
  1211. brentq, brenth, ridder, bisect, newton
  1212. fsolve : find roots in N dimensions.
  1213. elementwise.find_root : efficient elementwise 1-D root-finder
  1214. Notes
  1215. -----
  1216. `f` must be continuous.
  1217. Algorithm 748 with ``k=2`` is asymptotically the most efficient
  1218. algorithm known for finding roots of a four times continuously
  1219. differentiable function.
  1220. In contrast with Brent's algorithm, which may only decrease the length of
  1221. the enclosing bracket on the last step, Algorithm 748 decreases it each
  1222. iteration with the same asymptotic efficiency as it finds the root.
  1223. For easy statement of efficiency indices, assume that `f` has 4
  1224. continuous deriviatives.
  1225. For ``k=1``, the convergence order is at least 2.7, and with about
  1226. asymptotically 2 function evaluations per iteration, the efficiency
  1227. index is approximately 1.65.
  1228. For ``k=2``, the order is about 4.6 with asymptotically 3 function
  1229. evaluations per iteration, and the efficiency index 1.66.
  1230. For higher values of `k`, the efficiency index approaches
  1231. the kth root of ``(3k-2)``, hence ``k=1`` or ``k=2`` are
  1232. usually appropriate.
  1233. As mentioned in the parameter documentation, the computed root ``x0`` will
  1234. satisfy ``np.isclose(x, x0, atol=xtol, rtol=rtol)``, where ``x`` is the
  1235. exact root. In equation form, this terminating condition is ``abs(x - x0)
  1236. <= xtol + rtol * abs(x0)``.
  1237. The default value ``xtol=2e-12`` may lead to surprising behavior if one
  1238. expects `toms748` to always compute roots with relative error near machine
  1239. precision. Care should be taken to select `xtol` for the use case at hand.
  1240. Setting ``xtol=5e-324``, the smallest subnormal number, will ensure the
  1241. highest level of accuracy. Larger values of `xtol` may be useful for saving
  1242. function evaluations when a root is at or near zero in applications where
  1243. the tiny absolute differences available between floating point numbers near
  1244. zero are not meaningful.
  1245. References
  1246. ----------
  1247. .. [APS1995]
  1248. Alefeld, G. E. and Potra, F. A. and Shi, Yixun,
  1249. *Algorithm 748: Enclosing Zeros of Continuous Functions*,
  1250. ACM Trans. Math. Softw. Volume 221(1995)
  1251. doi = {10.1145/210089.210111}
  1252. Examples
  1253. --------
  1254. >>> def f(x):
  1255. ... return (x**3 - 1) # only one real root at x = 1
  1256. >>> from scipy import optimize
  1257. >>> root, results = optimize.toms748(f, 0, 2, full_output=True)
  1258. >>> root
  1259. 1.0
  1260. >>> results
  1261. converged: True
  1262. flag: converged
  1263. function_calls: 11
  1264. iterations: 5
  1265. root: 1.0
  1266. method: toms748
  1267. """
  1268. if xtol <= 0:
  1269. raise ValueError(f"xtol too small ({xtol:g} <= 0)")
  1270. if rtol < _rtol / 4:
  1271. raise ValueError(f"rtol too small ({rtol:g} < {_rtol/4:g})")
  1272. maxiter = operator.index(maxiter)
  1273. if maxiter < 1:
  1274. raise ValueError("maxiter must be greater than 0")
  1275. if not np.isfinite(a):
  1276. raise ValueError(f"a is not finite {a}")
  1277. if not np.isfinite(b):
  1278. raise ValueError(f"b is not finite {b}")
  1279. if a >= b:
  1280. raise ValueError(f"a and b are not an interval [{a}, {b}]")
  1281. if not k >= 1:
  1282. raise ValueError(f"k too small ({k} < 1)")
  1283. if not isinstance(args, tuple):
  1284. args = (args,)
  1285. f = _wrap_nan_raise(f)
  1286. solver = TOMS748Solver()
  1287. result = solver.solve(f, a, b, args=args, k=k, xtol=xtol, rtol=rtol,
  1288. maxiter=maxiter, disp=disp)
  1289. x, function_calls, iterations, flag = result
  1290. return _results_select(full_output, (x, function_calls, iterations, flag),
  1291. "toms748")