scipy.optimize.brenth#
- scipy.optimize.brenth(f, a, b, args=(), xtol=2e-12, rtol=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 with hyperbolic extrapolation.
A variation on the classic Brent routine to find a zero of the function f between the arguments a and b that uses hyperbolic extrapolation instead of inverse quadratic extrapolation. Bus & Dekker (1975) guarantee convergence for this method, claiming that the upper bound of function evaluations here is 4 or 5 times lesser than that for bisection. f(a) and f(b) cannot have the same signs. Generally, on a par with the brent routine, but not as heavily tested. It is a safe version of the secant method that uses hyperbolic extrapolation. The version here is by Chuck Harris, and implements Algorithm M of [BusAndDekker1975], where further details (convergence properties, additional remarks and such) can be found
- Parameters
- ffunction
Python function returning a number. 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.allclose(x, x0, atol=xtol, rtol=rtol)
, wherex
is the exact root. The parameter must be nonnegative. As withbrentq
, for nice functions the method will often satisfy the above condition withxtol/2
andrtol/2
.- rtolnumber, optional
The computed root
x0
will satisfynp.allclose(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
. As withbrentq
, for nice functions the method will often satisfy the above condition withxtol/2
andrtol/2
.- 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
- x0float
Zero 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
brentq
,brenth
,ridder
,bisect
,newton
1-D root-finding
fixed_point
scalar fixed-point finder
References
- BusAndDekker1975
Bus, J. C. P., Dekker, T. J., “Two Efficient Algorithms with Guaranteed Convergence for Finding a Zero of a Function”, ACM Transactions on Mathematical Software, Vol. 1, Issue 4, Dec. 1975, pp. 330-345. Section 3: “Algorithm M”. DOI:10.1145/355656.355659
Examples
>>> def f(x): ... return (x**2 - 1)
>>> from scipy import optimize
>>> root = optimize.brenth(f, -2, 0) >>> root -1.0
>>> root = optimize.brenth(f, 0, 2) >>> root 1.0