brentq#
- scipy.optimize.brentq(f, a, b, args=(), xtol=2e-12, rtol=np.float64(8.881784197001252e-16), maxiter=100, full_output=False, disp=True)[source]#
Find a root of a function in a bracketing interval using Brent’s method.
Uses the classic Brent’s method to find a root of the function f on the sign changing interval [a , b]. Generally considered the best of the rootfinding routines here. It is a safe version of the secant method that uses inverse quadratic extrapolation. Brent’s method combines root bracketing, interval bisection, and inverse quadratic interpolation. It is sometimes known as the van Wijngaarden-Dekker-Brent method. Brent (1973) claims convergence is guaranteed for functions computable within [a,b].
[Brent1973] provides the classic description of the algorithm. Another description can be found in a recent edition of Numerical Recipes, including [PressEtal1992]. A third description is at http://mathworld.wolfram.com/BrentsMethod.html. It should be easy to understand the algorithm just by reading our code. Our code diverges a bit from standard presentations: we choose a different formula for the extrapolation step.
- Parameters:
- ffunction
Python function returning a number. The function \(f\) must be continuous, and \(f(a)\) and \(f(b)\) must have opposite signs.
- ascalar
One end of the bracketing interval \([a, b]\).
- bscalar
The other end of the bracketing interval \([a, b]\).
- xtolnumber, optional
The computed root
x0
will satisfynp.isclose(x, x0, atol=xtol, rtol=rtol)
, wherex
is the exact root. The parameter must be positive. For nice functions, Brent’s method will often satisfy the above condition withxtol/2
andrtol/2
. [Brent1973]- rtolnumber, optional
The computed root
x0
will satisfynp.isclose(x, x0, atol=xtol, rtol=rtol)
, wherex
is the exact root. The parameter cannot be smaller than its default value of4*np.finfo(float).eps
. For nice functions, Brent’s method will often satisfy the above condition withxtol/2
andrtol/2
. [Brent1973]- maxiterint, optional
If convergence is not achieved in maxiter iterations, an error is raised. Must be >= 0.
- argstuple, optional
Containing extra arguments for the function f. f is called by
apply(f, (x)+args)
.- full_outputbool, optional
If full_output is False, the root is returned. If full_output is True, the return value is
(x, r)
, where x is the root, and r is aRootResults
object.- dispbool, optional
If True, raise RuntimeError if the algorithm didn’t converge. Otherwise, the convergence status is recorded in any
RootResults
return object.
- Returns:
- rootfloat
Root of f between a and b.
- r
RootResults
(present iffull_output = True
) Object containing information about the convergence. In particular,
r.converged
is True if the routine converged.
See also
fmin
,fmin_powell
,fmin_cg
,fmin_bfgs
,fmin_ncg
multivariate local optimizers
leastsq
nonlinear least squares minimizer
fmin_l_bfgs_b
,fmin_tnc
,fmin_cobyla
constrained multivariate optimizers
basinhopping
,differential_evolution
,brute
global optimizers
fminbound
,brent
,golden
,bracket
local scalar minimizers
fsolve
N-D root-finding
brenth
,ridder
,bisect
,newton
1-D root-finding
fixed_point
scalar fixed-point finder
elementwise.find_root
efficient elementwise 1-D root-finder
Notes
f must be continuous. f(a) and f(b) must have opposite signs.
As mentioned in the parameter documentation, the computed root
x0
will satisfynp.isclose(x, x0, atol=xtol, rtol=rtol)
, wherex
is the exact root. In equation form, this terminating condition isabs(x - x0) <= xtol + rtol * abs(x0)
.The default value
xtol=2e-12
may lead to surprising behavior if one expectsbrentq
to always compute roots with relative error near machine precision. Care should be taken to select xtol for the use case at hand. Settingxtol=5e-324
, the smallest subnormal number, will ensure the highest level of accuracy. Larger values of xtol may be useful for saving function evaluations when a root is at or near zero in applications where the tiny absolute differences available between floating point numbers near zero are not meaningful.References
[Brent1973] (1,2,3)Brent, R. P., Algorithms for Minimization Without Derivatives. Englewood Cliffs, NJ: Prentice-Hall, 1973. Ch. 3-4.
[PressEtal1992]Press, W. H.; Flannery, B. P.; Teukolsky, S. A.; and Vetterling, W. T. Numerical Recipes in FORTRAN: The Art of Scientific Computing, 2nd ed. Cambridge, England: Cambridge University Press, pp. 352-355, 1992. Section 9.3: “Van Wijngaarden-Dekker-Brent Method.”
Examples
>>> def f(x): ... return (x**2 - 1)
>>> from scipy import optimize
>>> root = optimize.brentq(f, -2, 0) >>> root -1.0
>>> root = optimize.brentq(f, 0, 2) >>> root 1.0