scipy.optimize.leastsq#
- scipy.optimize.leastsq(func, x0, args=(), Dfun=None, full_output=0, col_deriv=0, ftol=1.49012e-08, xtol=1.49012e-08, gtol=0.0, maxfev=0, epsfcn=None, factor=100, diag=None)[source]#
- Minimize the sum of squares of a set of equations. - x = arg min(sum(func(y)**2,axis=0)) y - Parameters
- funccallable
- Should take at least one (possibly length - Nvector) argument and returns- Mfloating point numbers. It must not return NaNs or fitting might fail.- Mmust be greater than or equal to- N.
- x0ndarray
- The starting estimate for the minimization. 
- argstuple, optional
- Any extra arguments to func are placed in this tuple. 
- Dfuncallable, optional
- A function or method to compute the Jacobian of func with derivatives across the rows. If this is None, the Jacobian will be estimated. 
- full_outputbool, optional
- non-zero to return all optional outputs. 
- col_derivbool, optional
- non-zero to specify that the Jacobian function computes derivatives down the columns (faster, because there is no transpose operation). 
- ftolfloat, optional
- Relative error desired in the sum of squares. 
- xtolfloat, optional
- Relative error desired in the approximate solution. 
- gtolfloat, optional
- Orthogonality desired between the function vector and the columns of the Jacobian. 
- maxfevint, optional
- The maximum number of calls to the function. If Dfun is provided, then the default maxfev is 100*(N+1) where N is the number of elements in x0, otherwise the default maxfev is 200*(N+1). 
- epsfcnfloat, optional
- A variable used in determining a suitable step length for the forward- difference approximation of the Jacobian (for Dfun=None). Normally the actual step length will be sqrt(epsfcn)*x If epsfcn is less than the machine precision, it is assumed that the relative errors are of the order of the machine precision. 
- factorfloat, optional
- A parameter determining the initial step bound ( - factor * || diag * x||). Should be in interval- (0.1, 100).
- diagsequence, optional
- N positive entries that serve as a scale factors for the variables. 
 
- Returns
- xndarray
- The solution (or the result of the last iteration for an unsuccessful call). 
- cov_xndarray
- The inverse of the Hessian. fjac and ipvt are used to construct an estimate of the Hessian. A value of None indicates a singular matrix, which means the curvature in parameters x is numerically flat. To obtain the covariance matrix of the parameters x, cov_x must be multiplied by the variance of the residuals – see curve_fit. 
- infodictdict
- a dictionary of optional outputs with the keys: - nfev
- The number of function calls 
- fvec
- The function evaluated at the output 
- fjac
- A permutation of the R matrix of a QR factorization of the final approximate Jacobian matrix, stored column wise. Together with ipvt, the covariance of the estimate can be approximated. 
- ipvt
- An integer array of length N which defines a permutation matrix, p, such that fjac*p = q*r, where r is upper triangular with diagonal elements of nonincreasing magnitude. Column j of p is column ipvt(j) of the identity matrix. 
- qtf
- The vector (transpose(q) * fvec). 
 
- mesgstr
- A string message giving information about the cause of failure. 
- ierint
- An integer flag. If it is equal to 1, 2, 3 or 4, the solution was found. Otherwise, the solution was not found. In either case, the optional output variable ‘mesg’ gives more information. 
 
 - See also - least_squares
- Newer interface to solve nonlinear least-squares problems with bounds on the variables. See - method=='lm'in particular.
 - Notes - “leastsq” is a wrapper around MINPACK’s lmdif and lmder algorithms. - cov_x is a Jacobian approximation to the Hessian of the least squares objective function. This approximation assumes that the objective function is based on the difference between some observed target data (ydata) and a (non-linear) function of the parameters f(xdata, params) - func(params) = ydata - f(xdata, params) - so that the objective function is - min sum((ydata - f(xdata, params))**2, axis=0) params - The solution, x, is always a 1-D array, regardless of the shape of x0, or whether x0 is a scalar. - Examples - >>> from scipy.optimize import leastsq >>> def func(x): ... return 2*(x-3)**2+1 >>> leastsq(func, 0) (array([2.99999999]), 1)