bundles / scipy latest / scipy / optimize / _lbfgsb_py / fmin_l_bfgs_b
function
scipy.optimize._lbfgsb_py:fmin_l_bfgs_b
Signature
def fmin_l_bfgs_b ( func , x0 , fprime = None , args = () , approx_grad = 0 , bounds = None , m = 10 , factr = 10000000.0 , pgtol = 1e-05 , epsilon = 1e-08 , iprint = <object object at 0x0000> , maxfun = 15000 , maxiter = 15000 , disp = <object object at 0x0000> , callback = None , maxls = 20 ) Summary
Minimize a function func using the L-BFGS-B algorithm.
Parameters
func: callable f(x,*args)Function to minimize.
x0: ndarrayInitial guess.
fprime: callable fprime(x,*args), optionalThe gradient of
func. If None, thenfuncreturns the function value and the gradient (f, g = func(x, *args)), unlessapprox_gradis True in which casefuncreturns onlyf.args: sequence, optionalArguments to pass to
funcandfprime.approx_grad: bool, optionalWhether to approximate the gradient numerically (in which case
funcreturns only the function value).bounds: list, optional(min, max)pairs for each element inx, defining the bounds on that parameter. Use None or +-inf for one ofminormaxwhen there is no bound in that direction.m: int, optionalThe maximum number of variable metric corrections used to define the limited memory matrix. (The limited memory BFGS method does not store the full hessian but uses this many terms in an approximation to it.)
factr: float, optionalThe iteration stops when
(f^k - f^{k+1})/max{|f^k|,|f^{k+1}|,1} <= factr * eps, whereepsis the machine precision, which is automatically generated by the code. Typical values forfactrare: 1e12 for low accuracy; 1e7 for moderate accuracy; 10.0 for extremely high accuracy. See Notes for relationship toftol, which is exposed (instead offactr) by the scipy.optimize.minimize interface to L-BFGS-B.pgtol: float, optionalThe iteration will stop when
max{|proj g_i | i = 1, ..., n} <= pgtolwhereproj g_iis the i-th component of the projected gradient.epsilon: float, optionalStep size used when
approx_gradis True, for numerically calculating the gradientiprint: int, optionalDeprecated option that previously controlled the text printed on the screen during the problem solution. Now the code does not emit any output and this keyword has no function.
disp: int, optionalDeprecated option that previously controlled the text printed on the screen during the problem solution. Now the code does not emit any output and this keyword has no function.
maxfun: int, optionalMaximum number of function evaluations. Note that this function may violate the limit because of evaluating gradients by numerical differentiation.
maxiter: int, optionalMaximum number of iterations.
callback: callable, optionalCalled after each iteration, as
callback(xk), wherexkis the current parameter vector.maxls: int, optionalMaximum number of line search steps (per iteration). Default is 20.
Returns
x: array_likeEstimated position of the minimum.
f: floatValue of
funcat the minimum.d: dictInformation dictionary.
d['warnflag'] is
0 if converged,
1 if too many function evaluations or too many iterations,
2 if stopped for another reason, given in d['task']
d['grad'] is the gradient at the minimum (should be 0 ish)
d['funcalls'] is the number of function calls made.
d['nit'] is the number of iterations.
Notes
SciPy uses a C-translated and modified version of the Fortran code, L-BFGS-B v3.0 (released April 25, 2011, BSD-3 licensed). Original Fortran version was written by Ciyou Zhu, Richard Byrd, Jorge Nocedal and, Jose Luis Morales.
Examples
Solve a linear regression problem via `fmin_l_bfgs_b`. To do this, first we define an objective function ``f(m, b) = (y - y_model)**2``, where `y` describes the observations and `y_model` the prediction of the linear model as ``y_model = m*x + b``. The bounds for the parameters, ``m`` and ``b``, are arbitrarily chosen as ``(0,5)`` and ``(5,10)`` for this example.import numpy as np from scipy.optimize import fmin_l_bfgs_b X = np.arange(0, 10, 1) M = 2 B = 3 Y = M * X + B def func(parameters, *args): x = args[0] y = args[1] m, b = parameters y_model = m*x + b error = sum(np.power((y - y_model), 2)) return error✓
initial_values = np.array([0.0, 1.0])
✓x_opt, f_opt, info = fmin_l_bfgs_b(func, x0=initial_values, args=(X, Y), approx_grad=True)✓
x_opt, f_opt
✗bounds = [(0, 5), (5, 10)] x_opt, f_op, info = fmin_l_bfgs_b(func, x0=initial_values, args=(X, Y), approx_grad=True, bounds=bounds)✓
x_opt, f_opt
✗See also
- minimize
Interface to minimization algorithms for multivariate functions. See the 'L-BFGS-B'
methodin particular. Note that theftoloption is made available via that interface, whilefactris provided via this interface, wherefactris the factor multiplying the default machine floating-point precision to arrive atftol:ftol = factr * numpy.finfo(float).eps.
Aliases
-
scipy.optimize.fmin_l_bfgs_b