scipy.optimize.diagbroyden#
- scipy.optimize.diagbroyden(F, xin, iter=None, alpha=None, verbose=False, maxiter=None, f_tol=None, f_rtol=None, x_tol=None, x_rtol=None, tol_norm=None, line_search='armijo', callback=None, **kw)#
Find a root of a function, using diagonal Broyden Jacobian approximation.
The Jacobian approximation is derived from previous iterations, by retaining only the diagonal of Broyden matrices.
Warning
This algorithm may be useful for specific problems, but whether it will work may depend strongly on the problem.
- Parameters
- Ffunction(x) -> f
Function whose root to find; should take and return an array-like object.
- xinarray_like
Initial guess for the solution
- alphafloat, optional
Initial guess for the Jacobian is (-1/alpha).
- iterint, optional
Number of iterations to make. If omitted (default), make as many as required to meet tolerances.
- verbosebool, optional
Print status to stdout on every iteration.
- maxiterint, optional
Maximum number of iterations to make. If more are needed to meet convergence, NoConvergence is raised.
- f_tolfloat, optional
Absolute tolerance (in max-norm) for the residual. If omitted, default is 6e-6.
- f_rtolfloat, optional
Relative tolerance for the residual. If omitted, not used.
- x_tolfloat, optional
Absolute minimum step size, as determined from the Jacobian approximation. If the step size is smaller than this, optimization is terminated as successful. If omitted, not used.
- x_rtolfloat, optional
Relative minimum step size. If omitted, not used.
- tol_normfunction(vector) -> scalar, optional
Norm to use in convergence check. Default is the maximum norm.
- line_search{None, ‘armijo’ (default), ‘wolfe’}, optional
Which type of a line search to use to determine the step size in the direction given by the Jacobian approximation. Defaults to ‘armijo’.
- callbackfunction, optional
Optional callback function. It is called on every iteration as
callback(x, f)
where x is the current solution and f the corresponding residual.
- Returns
- solndarray
An array (of similar array type as x0) containing the final solution.
- Raises
- NoConvergence
When a solution was not found.
See also
root
Interface to root finding algorithms for multivariate functions. See
method=='diagbroyden'
in particular.
Examples
The following functions define a system of nonlinear equations
>>> def fun(x): ... return [x[0] + 0.5 * (x[0] - x[1])**3 - 1.0, ... 0.5 * (x[1] - x[0])**3 + x[1]]
A solution can be obtained as follows.
>>> from scipy import optimize >>> sol = optimize.diagbroyden(fun, [0, 0]) >>> sol array([0.84116403, 0.15883384])