linprog(method=’interior-point’)¶
- 
scipy.optimize.linprog(c, method='interior-point', callback=None, options={'c0': None, 'A': None, 'b': None, 'postsolve_args': None, 'maxiter': 1000, 'tol': 1e-08, 'disp': False, 'alpha0': 0.99995, 'beta': 0.1, 'sparse': False, 'lstsq': False, 'sym_pos': True, 'cholesky': None, 'pc': True, 'ip': False, 'permc_spec': 'MMD_AT_PLUS_A'}, x0=None)
- Minimize a linear objective function subject to linear equality and non-negativity constraints using the interior point method of [4]. Linear programming is intended to solve problems of the following form: - Minimize: - c @ x - Subject to: - A @ x == b x >= 0 - Parameters
- c1D array
- Coefficients of the linear objective function to be minimized. 
- c0float
- Constant term in objective function due to fixed (and eliminated) variables. (Purely for display.) 
- A2D array
- 2D array such that - A @ x, gives the values of the equality constraints at- x.
- b1D array
- 1D array of values representing the right hand side of each equality constraint (row) in - A.
- callbackcallable, optional
- Callback function to be executed once per iteration. 
- postsolve_argstuple
- Data needed by _postsolve to convert the solution to the standard-form problem into the solution to the original problem. 
 
- Returns
- x1D array
- Solution vector. 
- statusint
- An integer representing the exit status of the optimization: - 0 : Optimization terminated successfully 1 : Iteration limit reached 2 : Problem appears to be infeasible 3 : Problem appears to be unbounded 4 : Serious numerical difficulties encountered 
- messagestr
- A string descriptor of the exit status of the optimization. 
- iterationint
- The number of iterations taken to solve the problem. 
 
 - See also - For documentation for the rest of the parameters, see - scipy.optimize.linprog- Options
- maxiterint (default = 1000)
- The maximum number of iterations of the algorithm. 
- tolfloat (default = 1e-8)
- Termination tolerance to be used for all termination criteria; see [4] Section 4.5. 
- dispbool (default = False)
- Set to - Trueif indicators of optimization status are to be printed to the console each iteration.
- alpha0float (default = 0.99995)
- The maximal step size for Mehrota’s predictor-corrector search direction; see \(\beta_{3}\) of [4] Table 8.1. 
- betafloat (default = 0.1)
- The desired reduction of the path parameter \(\mu\) (see [6]) when Mehrota’s predictor-corrector is not in use (uncommon). 
- sparsebool (default = False)
- Set to - Trueif the problem is to be treated as sparse after presolve. If either- A_eqor- A_ubis a sparse matrix, this option will automatically be set- True, and the problem will be treated as sparse even during presolve. If your constraint matrices contain mostly zeros and the problem is not very small (less than about 100 constraints or variables), consider setting- Trueor providing- A_eqand- A_ubas sparse matrices.
- lstsqbool (default = False)
- Set to - Trueif the problem is expected to be very poorly conditioned. This should always be left- Falseunless severe numerical difficulties are encountered. Leave this at the default unless you receive a warning message suggesting otherwise.
- sym_posbool (default = True)
- Leave - Trueif the problem is expected to yield a well conditioned symmetric positive definite normal equation matrix (almost always). Leave this at the default unless you receive a warning message suggesting otherwise.
- choleskybool (default = True)
- Set to - Trueif the normal equations are to be solved by explicit Cholesky decomposition followed by explicit forward/backward substitution. This is typically faster for problems that are numerically well-behaved.
- pcbool (default = True)
- Leave - Trueif the predictor-corrector method of Mehrota is to be used. This is almost always (if not always) beneficial.
- ipbool (default = False)
- Set to - Trueif the improved initial point suggestion due to [4] Section 4.3 is desired. Whether this is beneficial or not depends on the problem.
- permc_specstr (default = ‘MMD_AT_PLUS_A’)
- (Has effect only with - sparse = True,- lstsq = False,- sym_pos = True, and no SuiteSparse.) A matrix is factorized in each iteration of the algorithm. This option specifies how to permute the columns of the matrix for sparsity preservation. Acceptable values are:- NATURAL: natural ordering.
- MMD_ATA: minimum degree ordering on the structure of A^T A.
- MMD_AT_PLUS_A: minimum degree ordering on the structure of A^T+A.
- COLAMD: approximate minimum degree column ordering.
 - This option can impact the convergence of the interior point algorithm; test different values to determine which performs best for your problem. For more information, refer to - scipy.sparse.linalg.splu.
- unknown_optionsdict
- Optional arguments not used by this particular solver. If unknown_options is non-empty a warning is issued listing all unused options. 
 
 - Notes - This method implements the algorithm outlined in [4] with ideas from [8] and a structure inspired by the simpler methods of [6]. - The primal-dual path following method begins with initial ‘guesses’ of the primal and dual variables of the standard form problem and iteratively attempts to solve the (nonlinear) Karush-Kuhn-Tucker conditions for the problem with a gradually reduced logarithmic barrier term added to the objective. This particular implementation uses a homogeneous self-dual formulation, which provides certificates of infeasibility or unboundedness where applicable. - The default initial point for the primal and dual variables is that defined in [4] Section 4.4 Equation 8.22. Optionally (by setting initial point option - ip=True), an alternate (potentially improved) starting point can be calculated according to the additional recommendations of [4] Section 4.4.- A search direction is calculated using the predictor-corrector method (single correction) proposed by Mehrota and detailed in [4] Section 4.1. (A potential improvement would be to implement the method of multiple corrections described in [4] Section 4.2.) In practice, this is accomplished by solving the normal equations, [4] Section 5.1 Equations 8.31 and 8.32, derived from the Newton equations [4] Section 5 Equations 8.25 (compare to [4] Section 4 Equations 8.6-8.8). The advantage of solving the normal equations rather than 8.25 directly is that the matrices involved are symmetric positive definite, so Cholesky decomposition can be used rather than the more expensive LU factorization. - With default options, the solver used to perform the factorization depends on third-party software availability and the conditioning of the problem. - For dense problems, solvers are tried in the following order: - scipy.linalg.cho_factor
- scipy.linalg.solvewith option- sym_pos=True
- scipy.linalg.solvewith option- sym_pos=False
- scipy.linalg.lstsq
 - For sparse problems: - sksparse.cholmod.cholesky(if scikit-sparse and SuiteSparse are installed)
- scipy.sparse.linalg.factorized(if scikit-umfpack and SuiteSparse are installed)
- scipy.sparse.linalg.splu(which uses SuperLU distributed with SciPy)
- scipy.sparse.linalg.lsqr
 - If the solver fails for any reason, successively more robust (but slower) solvers are attempted in the order indicated. Attempting, failing, and re-starting factorization can be time consuming, so if the problem is numerically challenging, options can be set to bypass solvers that are failing. Setting - cholesky=Falseskips to solver 2,- sym_pos=Falseskips to solver 3, and- lstsq=Trueskips to solver 4 for both sparse and dense problems.- Potential improvements for combatting issues associated with dense columns in otherwise sparse problems are outlined in [4] Section 5.3 and [10] Section 4.1-4.2; the latter also discusses the alleviation of accuracy issues associated with the substitution approach to free variables. - After calculating the search direction, the maximum possible step size that does not activate the non-negativity constraints is calculated, and the smaller of this step size and unity is applied (as in [4] Section 4.1.) [4] Section 4.3 suggests improvements for choosing the step size. - The new point is tested according to the termination conditions of [4] Section 4.5. The same tolerance, which can be set using the - toloption, is used for all checks. (A potential improvement would be to expose the different tolerances to be set independently.) If optimality, unboundedness, or infeasibility is detected, the solve procedure terminates; otherwise it repeats.- The expected problem formulation differs between the top level - linprogmodule and the method specific solvers. The method specific solvers expect a problem in standard form:- Minimize: - c @ x - Subject to: - A @ x == b x >= 0 - Whereas the top level - linprogmodule expects a problem of form:- Minimize: - c @ x - Subject to: - A_ub @ x <= b_ub A_eq @ x == b_eq lb <= x <= ub - where - lb = 0and- ub = Noneunless set in- bounds.- The original problem contains equality, upper-bound and variable constraints whereas the method specific solver requires equality constraints and variable non-negativity. - linprogmodule converts the original problem to standard form by converting the simple bounds to upper bound constraints, introducing non-negative slack variables for inequality constraints, and expressing unbounded variables as the difference between two non-negative variables.- References - 4(1,2,3,4,5,6,7,8,9,10,11,12,13,14,15,16)
- Andersen, Erling D., and Knud D. Andersen. “The MOSEK interior point optimizer for linear programming: an implementation of the homogeneous algorithm.” High performance optimization. Springer US, 2000. 197-232. 
- 6(1,2)
- Freund, Robert M. “Primal-Dual Interior-Point Methods for Linear Programming based on Newton’s Method.” Unpublished Course Notes, March 2004. Available 2/25/2017 at https://ocw.mit.edu/courses/sloan-school-of-management/15-084j-nonlinear-programming-spring-2004/lecture-notes/lec14_int_pt_mthd.pdf 
- 8
- Andersen, Erling D., and Knud D. Andersen. “Presolving in linear programming.” Mathematical Programming 71.2 (1995): 221-245. 
- 9
- Bertsimas, Dimitris, and J. Tsitsiklis. “Introduction to linear programming.” Athena Scientific 1 (1997): 997. 
- 10
- Andersen, Erling D., et al. Implementation of interior point methods for large scale linear programming. HEC/Universite de Geneve, 1996. 
 
