Find a root of a vector function.
New in version 0.11.0.
Parameters : | fun : callable
x0 : ndarray
args : tuple, optional
method : str, optional
jac : bool or callable, optional
tol : float, optional
callback : function, optional
options : dict, optional
|
---|---|
Returns : | sol : Result
|
Notes
This section describes the available solvers that can be selected by the ‘method’ parameter. The default method is hybr.
Method hybr uses a modification of the Powell hybrid method as implemented in MINPACK [R88].
Method lm solves the system of nonlinear equations in a least squares sense using a modification of the Levenberg-Marquardt algorithm as implemented in MINPACK [R88].
Methods broyden1, broyden2, anderson, linearmixing, diagbroyden, excitingmixing, krylov are inexact Newton methods, with backtracking or full line searches [R89]. Each method corresponds to a particular Jacobian approximations. See nonlin for details.
Warning
The algorithms implemented for methods diagbroyden, linearmixing and excitingmixing may be useful for specific problems, but whether they will work may depend strongly on the problem.
References
[R88] | (1, 2, 3) More, Jorge J., Burton S. Garbow, and Kenneth E. Hillstrom. 1980. User Guide for MINPACK-1. |
[R89] | (1, 2) C. T. Kelley. 1995. Iterative Methods for Linear and Nonlinear Equations. Society for Industrial and Applied Mathematics. <http://www.siam.org/books/kelley/> |
Examples
The following functions define a system of nonlinear equations and its jacobian.
>>> def fun(x):
... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0,
... 0.5 * (x[1] - x[0])**3 + x[1]]
>>> def jac(x):
... return np.array([[1 + 1.5 * (x[0] - x[1])**2,
... -1.5 * (x[0] - x[1])**2],
... [-1.5 * (x[1] - x[0])**2,
... 1 + 1.5 * (x[1] - x[0])**2]])
A solution can be obtained as follows.
>>> from scipy import optimize
>>> sol = optimize.root(fun, [0, 0], jac=jac, method='hybr')
>>> sol.x
array([ 0.8411639, 0.1588361])