minimize(method=’trust-constr’)¶
-
scipy.optimize.
minimize
(fun, x0, args=(), method='trust-constr', hess=None, hessp=None, bounds=None, constraints=(), tol=None, callback=None, options={'grad': None, 'xtol': 1e-08, 'gtol': 1e-08, 'barrier_tol': 1e-08, 'sparse_jacobian': None, 'maxiter': 1000, 'verbose': 0, 'finite_diff_rel_step': None, 'initial_constr_penalty': 1.0, 'initial_tr_radius': 1.0, 'initial_barrier_parameter': 0.1, 'initial_barrier_tolerance': 0.1, 'factorization_method': None, 'disp': False}) Minimize a scalar function subject to constraints.
- Parameters
- gtolfloat, optional
Tolerance for termination by the norm of the Lagrangian gradient. The algorithm will terminate when both the infinity norm (i.e. max abs value) of the Lagrangian gradient and the constraint violation are smaller than
gtol
. Default is 1e-8.- xtolfloat, optional
Tolerance for termination by the change of the independent variable. The algorithm will terminate when
tr_radius < xtol
, wheretr_radius
is the radius of the trust region used in the algorithm. Default is 1e-8.- barrier_tolfloat, optional
Threshold on the barrier parameter for the algorithm termination. When inequality constraints are present the algorithm will terminate only when the barrier parameter is less than barrier_tol. Default is 1e-8.
- sparse_jacobian{bool, None}, optional
Determines how to represent Jacobians of the constraints. If bool, then Jacobians of all the constraints will be converted to the corresponding format. If None (default), then Jacobians won’t be converted, but the algorithm can proceed only if they all have the same format.
- initial_tr_radius: float, optional
Initial trust radius. The trust radius gives the maximum distance between solution points in consecutive iterations. It reflects the trust the algorithm puts in the local approximation of the optimization problem. For an accurate local approximation the trust-region should be large and for an approximation valid only close to the current point it should be a small one. The trust radius is automatically updated throughout the optimization process, with
initial_tr_radius
being its initial value. Default is 1 (recommended in [1], p. 19).- initial_constr_penaltyfloat, optional
Initial constraints penalty parameter. The penalty parameter is used for balancing the requirements of decreasing the objective function and satisfying the constraints. It is used for defining the merit function:
merit_function(x) = fun(x) + constr_penalty * constr_norm_l2(x)
, whereconstr_norm_l2(x)
is the l2 norm of a vector containing all the constraints. The merit function is used for accepting or rejecting trial points andconstr_penalty
weights the two conflicting goals of reducing objective function and constraints. The penalty is automatically updated throughout the optimization process, withinitial_constr_penalty
being its initial value. Default is 1 (recommended in [1], p 19).- initial_barrier_parameter, initial_barrier_tolerance: float, optional
Initial barrier parameter and initial tolerance for the barrier subproblem. Both are used only when inequality constraints are present. For dealing with optimization problems
min_x f(x)
subject to inequality constraintsc(x) <= 0
the algorithm introduces slack variables, solving the problemmin_(x,s) f(x) + barrier_parameter*sum(ln(s))
subject to the equality constraintsc(x) + s = 0
instead of the original problem. This subproblem is solved for increasing values ofbarrier_parameter
and with decreasing tolerances for the termination, starting withinitial_barrier_parameter
for the barrier parameter andinitial_barrier_tolerance
for the barrier subproblem barrier. Default is 0.1 for both values (recommended in [1] p. 19).- factorization_methodstring or None, optional
Method to factorize the Jacobian of the constraints. Use None (default) for the auto selection or one of:
‘NormalEquation’ (requires scikit-sparse)
‘AugmentedSystem’
‘QRFactorization’
‘SVDFactorization’
The methods ‘NormalEquation’ and ‘AugmentedSystem’ can be used only with sparse constraints. The projections required by the algorithm will be computed using, respectively, the the normal equation and the augmented system approaches explained in [1]. ‘NormalEquation’ computes the Cholesky factorization of
A A.T
and ‘AugmentedSystem’ performs the LU factorization of an augmented system. They usually provide similar results. ‘AugmentedSystem’ is used by default for sparse matrices.The methods ‘QRFactorization’ and ‘SVDFactorization’ can be used only with dense constraints. They compute the required projections using, respectively, QR and SVD factorizations. The ‘SVDFactorization’ method can cope with Jacobian matrices with deficient row rank and will be used whenever other factorization methods fail (which may imply the conversion of sparse matrices to a dense format when required). By default ‘QRFactorization’ is used for dense matrices.
- finite_diff_rel_stepNone or array_like, optional
Relative step size for the finite difference approximation.
- maxiterint, optional
Maximum number of algorithm iterations. Default is 1000.
- verbose{0, 1, 2}, optional
Level of algorithm’s verbosity:
0 (default) : work silently.
1 : display a termination report.
2 : display progress during iterations.
3 : display progress during iterations (more complete report).
- dispbool, optional
If True (default) then verbose will be set to 1 if it was 0.
- Returns
- `OptimizeResult` with the fields documented below. Note the following:
All values corresponding to the constraints are ordered as they were passed to the solver. And values corresponding to bounds constraints are put after other constraints.
All numbers of function, Jacobian or Hessian evaluations correspond to numbers of actual Python function calls. It means, for example, that if a Jacobian is estimated by finite differences then the number of Jacobian evaluations will be zero and the number of function evaluations will be incremented by all calls during the finite difference estimation.
- xndarray, shape (n,)
Solution found.
- optimalityfloat
Infinity norm of the Lagrangian gradient at the solution.
- constr_violationfloat
Maximum constraint violation at the solution.
- funfloat
Objective function at the solution.
- gradndarray, shape (n,)
Gradient of the objective function at the solution.
- lagrangian_gradndarray, shape (n,)
Gradient of the Lagrangian function at the solution.
- nitint
Total number of iterations.
- nfevinteger
Number of the objective function evaluations.
- ngevinteger
Number of the objective function gradient evaluations.
- nhevinteger
Number of the objective function Hessian evaluations.
- cg_niterint
Total number of the conjugate gradient method iterations.
- method{‘equality_constrained_sqp’, ‘tr_interior_point’}
Optimization method used.
- constrlist of ndarray
List of constraint values at the solution.
- jaclist of {ndarray, sparse matrix}
List of the Jacobian matrices of the constraints at the solution.
- vlist of ndarray
List of the Lagrange multipliers for the constraints at the solution. For an inequality constraint a positive multiplier means that the upper bound is active, a negative multiplier means that the lower bound is active and if a multiplier is zero it means the constraint is not active.
- constr_nfevlist of int
Number of constraint evaluations for each of the constraints.
- constr_njevlist of int
Number of Jacobian matrix evaluations for each of the constraints.
- constr_nhevlist of int
Number of Hessian evaluations for each of the constraints.
- tr_radiusfloat
Radius of the trust region at the last iteration.
- constr_penaltyfloat
Penalty parameter at the last iteration, see initial_constr_penalty.
- barrier_tolerancefloat
Tolerance for the barrier subproblem at the last iteration. Only for problems with inequality constraints.
- barrier_parameterfloat
Barrier parameter at the last iteration. Only for problems with inequality constraints.
- execution_timefloat
Total execution time.
- messagestr
Termination message.
- status{0, 1, 2, 3}
Termination status:
0 : The maximum number of function evaluations is exceeded.
1 : gtol termination condition is satisfied.
2 : xtol termination condition is satisfied.
3 : callback function requested termination.
- cg_stop_condint
Reason for CG subproblem termination at the last iteration:
0 : CG subproblem not evaluated.
1 : Iteration limit was reached.
2 : Reached the trust-region boundary.
3 : Negative curvature detected.
4 : Tolerance was satisfied.
References