fmin_powell#
- scipy.optimize.fmin_powell(func, x0, args=(), xtol=0.0001, ftol=0.0001, maxiter=None, maxfun=None, full_output=0, disp=1, retall=0, callback=None, direc=None)[source]#
- Minimize a function using modified Powell’s method. - This method only uses function values, not derivatives. - Parameters:
- funccallable f(x,*args)
- Objective function to be minimized. 
- x0ndarray
- Initial guess. 
- argstuple, optional
- Extra arguments passed to func. 
- xtolfloat, optional
- Line-search error tolerance. 
- ftolfloat, optional
- Relative error in - func(xopt)acceptable for convergence.
- maxiterint, optional
- Maximum number of iterations to perform. 
- maxfunint, optional
- Maximum number of function evaluations to make. 
- full_outputbool, optional
- If True, - fopt,- xi,- direc,- iter,- funcalls, and- warnflagare returned.
- dispbool, optional
- If True, print convergence messages. 
- retallbool, optional
- If True, return a list of the solution at each iteration. 
- callbackcallable, optional
- An optional user-supplied function, called after each iteration. Called as - callback(xk), where- xkis the current parameter vector.
- direcndarray, optional
- Initial fitting step and parameter order set as an (N, N) array, where N is the number of fitting parameters in x0. Defaults to step size 1.0 fitting all parameters simultaneously ( - np.eye((N, N))). To prevent initial consideration of values in a step or to change initial step size, set to 0 or desired step size in the Jth position in the Mth block, where J is the position in x0 and M is the desired evaluation step, with steps being evaluated in index order. Step size and ordering will change freely as minimization proceeds.
 
- Returns:
- xoptndarray
- Parameter which minimizes func. 
- foptnumber
- Value of function at minimum: - fopt = func(xopt).
- direcndarray
- Current direction set. 
- iterint
- Number of iterations. 
- funcallsint
- Number of function calls made. 
- warnflagint
- Integer warning flag:
- 1 : Maximum number of function evaluations. 2 : Maximum number of iterations. 3 : NaN result encountered. 4 : The result is out of the provided bounds. 
 
- allvecslist
- List of solutions at each iteration. 
 
 - See also - minimize
- Interface to unconstrained minimization algorithms for multivariate functions. See the ‘Powell’ method in particular. 
 - Notes - Uses a modification of Powell’s method to find the minimum of a function of N variables. Powell’s method is a conjugate direction method. - The algorithm has two loops. The outer loop merely iterates over the inner loop. The inner loop minimizes over each current direction in the direction set. At the end of the inner loop, if certain conditions are met, the direction that gave the largest decrease is dropped and replaced with the difference between the current estimated x and the estimated x from the beginning of the inner-loop. - The technical conditions for replacing the direction of greatest increase amount to checking that - No further gain can be made along the direction of greatest increase from that iteration. 
- The direction of greatest increase accounted for a large sufficient fraction of the decrease in the function value from that iteration of the inner loop. 
 - References - Powell M.J.D. (1964) An efficient method for finding the minimum of a function of several variables without calculating derivatives, Computer Journal, 7 (2):155-162. - Press W., Teukolsky S.A., Vetterling W.T., and Flannery B.P.: Numerical Recipes (any edition), Cambridge University Press - Examples - >>> def f(x): ... return x**2 - >>> from scipy import optimize - >>> minimum = optimize.fmin_powell(f, -1) Optimization terminated successfully. Current function value: 0.000000 Iterations: 2 Function evaluations: 16 >>> minimum array(0.0)