scipy 1.8.0 Pypi GitHub Homepage
Other Docs
NotesParametersReturnsBackRef
minimize(fun, x0, args=(), method=None, jac=None, hess=None, hessp=None, bounds=None, constraints=(), tol=None, callback=None, options=None)

Notes

This section describes the available solvers that can be selected by the 'method' parameter. The default method is BFGS.

Unconstrained minimization

Method CG <optimize.minimize-cg> uses a nonlinear conjugate gradient algorithm by Polak and Ribiere, a variant of the Fletcher-Reeves method described in pp.120-122. Only the first derivatives are used.

Method BFGS <optimize.minimize-bfgs> uses the quasi-Newton method of Broyden, Fletcher, Goldfarb, and Shanno (BFGS) pp. 136. It uses the first derivatives only. BFGS has proven good performance even for non-smooth optimizations. This method also returns an approximation of the Hessian inverse, stored as :None:None:`hess_inv` in the OptimizeResult object.

Method Newton-CG <optimize.minimize-newtoncg> uses a Newton-CG algorithm pp. 168 (also known as the truncated Newton method). It uses a CG method to the compute the search direction. See also TNC method for a box-constrained minimization with a similar algorithm. Suitable for large-scale problems.

Method dogleg <optimize.minimize-dogleg> uses the dog-leg trust-region algorithm for unconstrained minimization. This algorithm requires the gradient and Hessian; furthermore the Hessian is required to be positive definite.

Method trust-ncg <optimize.minimize-trustncg> uses the Newton conjugate gradient trust-region algorithm for unconstrained minimization. This algorithm requires the gradient and either the Hessian or a function that computes the product of the Hessian with a given vector. Suitable for large-scale problems.

Method trust-krylov <optimize.minimize-trustkrylov> uses the Newton GLTR trust-region algorithm , for unconstrained minimization. This algorithm requires the gradient and either the Hessian or a function that computes the product of the Hessian with a given vector. Suitable for large-scale problems. On indefinite problems it requires usually less iterations than the :None:None:`trust-ncg` method and is recommended for medium and large-scale problems.

Method trust-exact <optimize.minimize-trustexact> is a trust-region method for unconstrained minimization in which quadratic subproblems are solved almost exactly . This algorithm requires the gradient and the Hessian (which is not required to be positive definite). It is, in many situations, the Newton method to converge in fewer iterations and the most recommended for small and medium-size problems.

Bound-Constrained minimization

Method Nelder-Mead <optimize.minimize-neldermead> uses the Simplex algorithm , . This algorithm is robust in many applications. However, if numerical computation of derivative can be trusted, other algorithms using the first and/or second derivatives information might be preferred for their better performance in general.

Method L-BFGS-B <optimize.minimize-lbfgsb> uses the L-BFGS-B algorithm , for bound constrained minimization.

Method Powell <optimize.minimize-powell> is a modification of Powell's method , which is a conjugate direction method. It performs sequential one-dimensional minimizations along each vector of the directions set (:None:None:`direc` field in :None:None:`options` and :None:None:`info`), which is updated at each iteration of the main minimization loop. The function need not be differentiable, and no derivatives are taken. If bounds are not provided, then an unbounded line search will be used. If bounds are provided and the initial guess is within the bounds, then every function evaluation throughout the minimization procedure will be within the bounds. If bounds are provided, the initial guess is outside the bounds, and :None:None:`direc` is full rank (default has full rank), then some function evaluations during the first iteration may be outside the bounds, but every function evaluation after the first iteration will be within the bounds. If :None:None:`direc` is not full rank, then some parameters may not be optimized and the solution is not guaranteed to be within the bounds.

Method TNC <optimize.minimize-tnc> uses a truncated Newton algorithm , to minimize a function with variables subject to bounds. This algorithm uses gradient information; it is also called Newton Conjugate-Gradient. It differs from the Newton-CG method described above as it wraps a C implementation and allows each variable to be given upper and lower bounds.

Constrained Minimization

Method COBYLA <optimize.minimize-cobyla> uses the Constrained Optimization BY Linear Approximation (COBYLA) method , , . The algorithm is based on linear approximations to the objective function and each constraint. The method wraps a FORTRAN implementation of the algorithm. The constraints functions 'fun' may return either a single number or an array or list of numbers.

Method SLSQP <optimize.minimize-slsqp> uses Sequential Least SQuares Programming to minimize a function of several variables with any combination of bounds, equality and inequality constraints. The method wraps the SLSQP Optimization subroutine originally implemented by Dieter Kraft . Note that the wrapper handles infinite values in bounds by converting them into large floating values.

Method trust-constr <optimize.minimize-trustconstr> is a trust-region algorithm for constrained optimization. It swiches between two implementations depending on the problem definition. It is the most versatile constrained minimization algorithm implemented in SciPy and the most appropriate for large-scale problems. For equality constrained problems it is an implementation of Byrd-Omojokun Trust-Region SQP method described in and in , p. 549. When inequality constraints are imposed as well, it swiches to the trust-region interior point method described in . This interior point algorithm, in turn, solves inequality constraints by introducing slack variables and solving a sequence of equality-constrained barrier problems for progressively smaller values of the barrier parameter. The previously described equality constrained SQP method is used to solve the subproblems with increasing levels of accuracy as the iterate gets closer to a solution.

Finite-Difference Options

For Method trust-constr <optimize.minimize-trustconstr> the gradient and the Hessian may be approximated using three finite-difference schemes: {'2-point', '3-point', 'cs'}. The scheme 'cs' is, potentially, the most accurate but it requires the function to correctly handle complex inputs and to be differentiable in the complex plane. The scheme '3-point' is more accurate than '2-point' but requires twice as many operations. If the gradient is estimated via finite-differences the Hessian must be estimated using one of the quasi-Newton strategies.

Method specific options for the :None:None:`hess` keyword

+--------------+------+----------+-------------------------+-----+ | method/Hess | None | callable | '2-point/'3-point'/'cs' | HUS | +==============+======+==========+=========================+=====+ | Newton-CG | x | (n, n) | x | x | | | | LO | | | +--------------+------+----------+-------------------------+-----+ | dogleg | | (n, n) | | | +--------------+------+----------+-------------------------+-----+ | trust-ncg | | (n, n) | x | x | +--------------+------+----------+-------------------------+-----+ | trust-krylov | | (n, n) | x | x | +--------------+------+----------+-------------------------+-----+ | trust-exact | | (n, n) | | | +--------------+------+----------+-------------------------+-----+ | trust-constr | x | (n, n) | x | x | | | | LO | | | | | | sp | | | +--------------+------+----------+-------------------------+-----+

where LO=LinearOperator, sp=Sparse matrix, HUS=HessianUpdateStrategy

Custom minimizers

It may be useful to pass a custom minimization method, for example when using a frontend to this method such as scipy.optimize.basinhopping or a different library. You can simply pass a callable as the method parameter.

The callable is called as method(fun, x0, args, **kwargs, **options) where kwargs corresponds to any other parameters passed to minimize (such as :None:None:`callback`, :None:None:`hess`, etc.), except the :None:None:`options` dict, which has its contents also passed as method parameters pair by pair. Also, if :None:None:`jac` has been passed as a bool type, :None:None:`jac` and :None:None:`fun` are mangled so that :None:None:`fun` returns just the function values and :None:None:`jac` is converted to a function returning the Jacobian. The method shall return an OptimizeResult object.

The provided method callable must be able to accept (and possibly ignore) arbitrary parameters; the set of parameters accepted by minimize may expand in future versions and then these parameters will be passed to the method. You can find an example in the scipy.optimize tutorial.

versionadded

Parameters

fun : callable

The objective function to be minimized.

fun(x, *args) -> float

where x is an 1-D array with shape (n,) and args is a tuple of the fixed parameters needed to completely specify the function.

x0 : ndarray, shape (n,)

Initial guess. Array of real elements of size (n,), where n is the number of independent variables.

args : tuple, optional

Extra arguments passed to the objective function and its derivatives (:None:None:`fun`, :None:None:`jac` and :None:None:`hess` functions).

method : str or callable, optional

Type of solver. Should be one of

  • 'Nelder-Mead' (see here) <optimize.minimize-neldermead>

  • 'Powell' (see here) <optimize.minimize-powell>

  • 'CG' (see here) <optimize.minimize-cg>

  • 'BFGS' (see here) <optimize.minimize-bfgs>

  • 'Newton-CG' (see here) <optimize.minimize-newtoncg>

  • 'L-BFGS-B' (see here) <optimize.minimize-lbfgsb>

  • 'TNC' (see here) <optimize.minimize-tnc>

  • 'COBYLA' (see here) <optimize.minimize-cobyla>

  • 'SLSQP' (see here) <optimize.minimize-slsqp>

  • 'trust-constr' (see here) <optimize.minimize-trustconstr>

  • 'dogleg' (see here) <optimize.minimize-dogleg>

  • 'trust-ncg' (see here) <optimize.minimize-trustncg>

  • 'trust-exact' (see here) <optimize.minimize-trustexact>

  • 'trust-krylov' (see here) <optimize.minimize-trustkrylov>

  • custom - a callable object (added in version 0.14.0), see below for description.

If not given, chosen to be one of BFGS , L-BFGS-B , SLSQP , depending on whether or not the problem has constraints or bounds.

jac : {callable, '2-point', '3-point', 'cs', bool}, optional

Method for computing the gradient vector. Only for CG, BFGS, Newton-CG, L-BFGS-B, TNC, SLSQP, dogleg, trust-ncg, trust-krylov, trust-exact and trust-constr. If it is a callable, it should be a function that returns the gradient vector:

jac(x, *args) -> array_like, shape (n,)

where x is an array with shape (n,) and args is a tuple with the fixed parameters. If :None:None:`jac` is a Boolean and is True, :None:None:`fun` is assumed to return a tuple (f, g) containing the objective function and the gradient. Methods 'Newton-CG', 'trust-ncg', 'dogleg', 'trust-exact', and 'trust-krylov' require that either a callable be supplied, or that :None:None:`fun` return the objective and gradient. If None or False, the gradient will be estimated using 2-point finite difference estimation with an absolute step size. Alternatively, the keywords {'2-point', '3-point', 'cs'} can be used to select a finite difference scheme for numerical estimation of the gradient with a relative step size. These finite difference schemes obey any specified :None:None:`bounds`.

hess : {callable, '2-point', '3-point', 'cs', HessianUpdateStrategy}, optional

Method for computing the Hessian matrix. Only for Newton-CG, dogleg, trust-ncg, trust-krylov, trust-exact and trust-constr. If it is callable, it should return the Hessian matrix:

hess(x, *args) -> {LinearOperator, spmatrix, array}, (n, n)

where x is a (n,) ndarray and args is a tuple with the fixed parameters. The keywords {'2-point', '3-point', 'cs'} can also be used to select a finite difference scheme for numerical estimation of the hessian. Alternatively, objects implementing the HessianUpdateStrategy interface can be used to approximate the Hessian. Available quasi-Newton methods implementing this interface are:

Not all of the options are available for each of the methods; for availability refer to the notes.

hessp : callable, optional

Hessian of objective function times an arbitrary vector p. Only for Newton-CG, trust-ncg, trust-krylov, trust-constr. Only one of hessp or :None:None:`hess` needs to be given. If :None:None:`hess` is provided, then hessp will be ignored. hessp must compute the Hessian times an arbitrary vector:

hessp(x, p, *args) -> ndarray shape (n,)

where x is a (n,) ndarray, p is an arbitrary vector with dimension (n,) and args is a tuple with the fixed parameters.

bounds : sequence or `Bounds`, optional

Bounds on variables for Nelder-Mead, L-BFGS-B, TNC, SLSQP, Powell, and trust-constr methods. There are two ways to specify the bounds:

  1. Instance of Bounds class.

  2. Sequence of (min, max) pairs for each element in x. None is used to specify no bound.

constraints : {Constraint, dict} or List of {Constraint, dict}, optional

Constraints definition. Only for COBYLA, SLSQP and trust-constr.

Constraints for 'trust-constr' are defined as a single object or a list of objects specifying constraints to the optimization problem. Available constraints are:

Constraints for COBYLA, SLSQP are defined as a list of dictionaries. Each dictionary with fields:

type

type

fun

fun

jac

jac

args

args

Equality constraint means that the constraint function result is to be zero whereas inequality means that it is to be non-negative. Note that COBYLA only supports inequality constraints.

tol : float, optional

Tolerance for termination. When :None:None:`tol` is specified, the selected minimization algorithm sets some relevant solver-specific tolerance(s) equal to :None:None:`tol`. For detailed control, use solver-specific options.

options : dict, optional

A dictionary of solver options. All methods accept the following generic options:

maxiter

maxiter

disp

disp

For method-specific options, see show_options() .

callback : callable, optional

Called after each iteration. For 'trust-constr' it is a callable with the signature:

callback(xk, OptimizeResult state) -> bool

where xk is the current parameter vector. and state is an OptimizeResult object, with the same fields as the ones from the return. If callback returns True the algorithm execution is terminated. For all the other methods, the signature is:

callback(xk)

where xk is the current parameter vector.

Returns

res : OptimizeResult

The optimization result represented as a OptimizeResult object. Important attributes are: x the solution array, success a Boolean flag indicating if the optimizer exited successfully and message which describes the cause of the termination. See OptimizeResult for a description of other attributes.

Minimization of scalar function of one or more variables.

See Also

minimize_scalar

Interface to minimization algorithms for scalar univariate functions

show_options

Additional options accepted by the solvers

Examples

Let us consider the problem of minimizing the Rosenbrock function. This function (and its respective derivatives) is implemented in rosen (resp. rosen_der , rosen_hess ) in the scipy.optimize .

>>> from scipy.optimize import minimize, rosen, rosen_der

A simple application of the Nelder-Mead method is:

>>> x0 = [1.3, 0.7, 0.8, 1.9, 1.2]
... res = minimize(rosen, x0, method='Nelder-Mead', tol=1e-6)
... res.x array([ 1., 1., 1., 1., 1.])

Now using the BFGS algorithm, using the first derivative and a few options:

>>> res = minimize(rosen, x0, method='BFGS', jac=rosen_der,
...  options={'gtol': 1e-6, 'disp': True}) Optimization terminated successfully. Current function value: 0.000000 Iterations: 26 Function evaluations: 31 Gradient evaluations: 31
>>> res.x
array([ 1.,  1.,  1.,  1.,  1.])
>>> print(res.message)
Optimization terminated successfully.
>>> res.hess_inv
array([[ 0.00749589,  0.01255155,  0.02396251,  0.04750988,  0.09495377],  # may vary
       [ 0.01255155,  0.02510441,  0.04794055,  0.09502834,  0.18996269],
       [ 0.02396251,  0.04794055,  0.09631614,  0.19092151,  0.38165151],
       [ 0.04750988,  0.09502834,  0.19092151,  0.38341252,  0.7664427 ],
       [ 0.09495377,  0.18996269,  0.38165151,  0.7664427,   1.53713523]])

Next, consider a minimization problem with several constraints (namely Example 16.4 from ). The objective function is:

>>> fun = lambda x: (x[0] - 1)**2 + (x[1] - 2.5)**2

There are three constraints defined as:

>>> cons = ({'type': 'ineq', 'fun': lambda x:  x[0] - 2 * x[1] + 2},
...  {'type': 'ineq', 'fun': lambda x: -x[0] - 2 * x[1] + 6},
...  {'type': 'ineq', 'fun': lambda x: -x[0] + 2 * x[1] + 2})

And variables must be positive, hence the following bounds:

>>> bnds = ((0, None), (0, None))

The optimization problem is solved using the SLSQP method as:

>>> res = minimize(fun, (2, 0), method='SLSQP', bounds=bnds,
...  constraints=cons)

It should converge to the theoretical solution (1.4 ,1.7).

See :

Back References

The following pages refer to to this document either explicitly or contain code examples using this.

scipy.optimize._optimize.fmin_cg scipy.optimize._differentialevolution.DifferentialEvolutionSolver scipy.optimize._optimize.show_options scipy.optimize._minimize.minimize scipy.optimize._differentialevolution.differential_evolution papyri scipy.optimize._minimize.minimize_scalar scipy.optimize._lbfgsb_py.fmin_l_bfgs_b scipy.optimize._lbfgsb_py._minimize_lbfgsb

Local connectivity graph

Hover to see nodes names; edges to Self not shown, Caped at 50 nodes.

Using a canvas is more power efficient and can get hundred of nodes ; but does not allow hyperlinks; , arrows or text (beyond on hover)

SVG is more flexible but power hungry; and does not scale well to 50 + nodes.

All aboves nodes referred to, (or are referred from) current nodes; Edges from Self to other have been omitted (or all nodes would be connected to the central node "self" which is not useful). Nodes are colored by the library they belong to, and scaled with the number of references pointing them


GitHub : /scipy/optimize/_minimize.py#45
type: <class 'function'>
Commit: