SciPy

scipy.optimize.root_scalar

scipy.optimize.root_scalar(f, args=(), method=None, bracket=None, fprime=None, fprime2=None, x0=None, x1=None, xtol=None, rtol=None, maxiter=None, options=None)[source]

Find a root of a scalar function.

Parameters
fcallable

A function to find a root of.

argstuple, optional

Extra arguments passed to the objective function and its derivative(s).

methodstr, optional

Type of solver. Should be one of

bracket: A sequence of 2 floats, optional

An interval bracketing a root. f(x, *args) must have different signs at the two endpoints.

x0float, optional

Initial guess.

x1float, optional

A second guess.

fprimebool or callable, optional

If fprime is a boolean and is True, f is assumed to return the value of the objective function and of the derivative. fprime can also be a callable returning the derivative of f. In this case, it must accept the same arguments as f.

fprime2bool or callable, optional

If fprime2 is a boolean and is True, f is assumed to return the value of the objective function and of the first and second derivatives. fprime2 can also be a callable returning the second derivative of f. In this case, it must accept the same arguments as f.

xtolfloat, optional

Tolerance (absolute) for termination.

rtolfloat, optional

Tolerance (relative) for termination.

maxiterint, optional

Maximum number of iterations.

optionsdict, optional

A dictionary of solver options. E.g., k, see show_options() for details.

Returns
solRootResults

The solution represented as a RootResults object. Important attributes are: root the solution , converged a boolean flag indicating if the algorithm exited successfully and flag which describes the cause of the termination. See RootResults for a description of other attributes.

See also

show_options

Additional options accepted by the solvers

root

Find a root of a vector function.

Notes

This section describes the available solvers that can be selected by the ‘method’ parameter.

The default is to use the best method available for the situation presented. If a bracket is provided, it may use one of the bracketing methods. If a derivative and an initial value are specified, it may select one of the derivative-based methods. If no method is judged applicable, it will raise an Exception.

Examples

Find the root of a simple cubic

>>> from scipy import optimize
>>> def f(x):
...     return (x**3 - 1)  # only one real root at x = 1
>>> def fprime(x):
...     return 3*x**2

The brentq method takes as input a bracket

>>> sol = optimize.root_scalar(f, bracket=[0, 3], method='brentq')
>>> sol.root, sol.iterations, sol.function_calls
(1.0, 10, 11)

The newton method takes as input a single point and uses the derivative(s)

>>> sol = optimize.root_scalar(f, x0=0.2, fprime=fprime, method='newton')
>>> sol.root, sol.iterations, sol.function_calls
(1.0, 11, 22)

The function can provide the value and derivative(s) in a single call.

>>> def f_p_pp(x):
...     return (x**3 - 1), 3*x**2, 6*x
>>> sol = optimize.root_scalar(f_p_pp, x0=0.2, fprime=True, method='newton')
>>> sol.root, sol.iterations, sol.function_calls
(1.0, 11, 11)
>>> sol = optimize.root_scalar(f_p_pp, x0=0.2, fprime=True, fprime2=True, method='halley')
>>> sol.root, sol.iterations, sol.function_calls
(1.0, 7, 8)

Previous topic

scipy.optimize.curve_fit

Next topic

scipy.optimize.brentq