| """Interior-point method for linear programming |
| |
| The *interior-point* method uses the primal-dual path following algorithm |
| outlined in [1]_. This algorithm supports sparse constraint matrices and |
| is typically faster than the simplex methods, especially for large, sparse |
| problems. Note, however, that the solution returned may be slightly less |
| accurate than those of the simplex methods and will not, in general, |
| correspond with a vertex of the polytope defined by the constraints. |
| |
| .. versionadded:: 1.0.0 |
| |
| References |
| ---------- |
| .. [1] 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. |
| """ |
| |
|
|
| import numpy as np |
| import scipy as sp |
| import scipy.sparse as sps |
| from warnings import warn |
| from scipy.linalg import LinAlgError |
| from ._optimize import OptimizeWarning, OptimizeResult, _check_unknown_options |
| from ._linprog_util import _postsolve |
| has_umfpack = True |
| has_cholmod = True |
| try: |
| import sksparse |
| from sksparse.cholmod import cholesky as cholmod |
| from sksparse.cholmod import analyze as cholmod_analyze |
| except ImportError: |
| has_cholmod = False |
| try: |
| import scikits.umfpack |
| except ImportError: |
| has_umfpack = False |
|
|
|
|
| def _get_solver(M, sparse=False, lstsq=False, sym_pos=True, |
| cholesky=True, permc_spec='MMD_AT_PLUS_A'): |
| """ |
| Given solver options, return a handle to the appropriate linear system |
| solver. |
| |
| Parameters |
| ---------- |
| M : 2-D array |
| As defined in [4] Equation 8.31 |
| sparse : bool (default = False) |
| True if the system to be solved is sparse. This is typically set |
| True when the original ``A_ub`` and ``A_eq`` arrays are sparse. |
| lstsq : bool (default = False) |
| True if the system is ill-conditioned and/or (nearly) singular and |
| thus a more robust least-squares solver is desired. This is sometimes |
| needed as the solution is approached. |
| sym_pos : bool (default = True) |
| True if the system matrix is symmetric positive definite |
| Sometimes this needs to be set false as the solution is approached, |
| even when the system should be symmetric positive definite, due to |
| numerical difficulties. |
| cholesky : bool (default = True) |
| True if the system is to be solved by Cholesky, rather than LU, |
| decomposition. This is typically faster unless the problem is very |
| small or prone to numerical difficulties. |
| permc_spec : str (default = 'MMD_AT_PLUS_A') |
| Sparsity preservation strategy used by SuperLU. 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. |
| |
| See SuperLU documentation. |
| |
| Returns |
| ------- |
| solve : function |
| Handle to the appropriate solver function |
| |
| """ |
| try: |
| if sparse: |
| if lstsq: |
| def solve(r, sym_pos=False): |
| return sps.linalg.lsqr(M, r)[0] |
| elif cholesky: |
| try: |
| |
| |
| _get_solver.cholmod_factor.cholesky_inplace(M) |
| except Exception: |
| _get_solver.cholmod_factor = cholmod_analyze(M) |
| _get_solver.cholmod_factor.cholesky_inplace(M) |
| solve = _get_solver.cholmod_factor |
| else: |
| if has_umfpack and sym_pos: |
| solve = sps.linalg.factorized(M) |
| else: |
| solve = sps.linalg.splu(M, permc_spec=permc_spec).solve |
|
|
| else: |
| if lstsq: |
| def solve(r): |
| return sp.linalg.lstsq(M, r)[0] |
| elif cholesky: |
| L = sp.linalg.cho_factor(M) |
|
|
| def solve(r): |
| return sp.linalg.cho_solve(L, r) |
| else: |
| |
| |
| def solve(r, sym_pos=sym_pos): |
| if sym_pos: |
| return sp.linalg.solve(M, r, assume_a="pos") |
| else: |
| return sp.linalg.solve(M, r) |
| |
| |
| |
| |
| except KeyboardInterrupt: |
| raise |
| except Exception: |
| return None |
| return solve |
|
|
|
|
| def _get_delta(A, b, c, x, y, z, tau, kappa, gamma, eta, sparse=False, |
| lstsq=False, sym_pos=True, cholesky=True, pc=True, ip=False, |
| permc_spec='MMD_AT_PLUS_A'): |
| """ |
| Given standard form problem defined by ``A``, ``b``, and ``c``; |
| current variable estimates ``x``, ``y``, ``z``, ``tau``, and ``kappa``; |
| algorithmic parameters ``gamma and ``eta; |
| and options ``sparse``, ``lstsq``, ``sym_pos``, ``cholesky``, ``pc`` |
| (predictor-corrector), and ``ip`` (initial point improvement), |
| get the search direction for increments to the variable estimates. |
| |
| Parameters |
| ---------- |
| As defined in [4], except: |
| sparse : bool |
| True if the system to be solved is sparse. This is typically set |
| True when the original ``A_ub`` and ``A_eq`` arrays are sparse. |
| lstsq : bool |
| True if the system is ill-conditioned and/or (nearly) singular and |
| thus a more robust least-squares solver is desired. This is sometimes |
| needed as the solution is approached. |
| sym_pos : bool |
| True if the system matrix is symmetric positive definite |
| Sometimes this needs to be set false as the solution is approached, |
| even when the system should be symmetric positive definite, due to |
| numerical difficulties. |
| cholesky : bool |
| True if the system is to be solved by Cholesky, rather than LU, |
| decomposition. This is typically faster unless the problem is very |
| small or prone to numerical difficulties. |
| pc : bool |
| True if the predictor-corrector method of Mehrota is to be used. This |
| is almost always (if not always) beneficial. Even though it requires |
| the solution of an additional linear system, the factorization |
| is typically (implicitly) reused so solution is efficient, and the |
| number of algorithm iterations is typically reduced. |
| ip : bool |
| True if the improved initial point suggestion due to [4] section 4.3 |
| is desired. It's unclear whether this is beneficial. |
| permc_spec : str (default = 'MMD_AT_PLUS_A') |
| (Has effect only with ``sparse = True``, ``lstsq = False``, ``sym_pos = |
| True``.) 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``. |
| |
| Returns |
| ------- |
| Search directions as defined in [4] |
| |
| References |
| ---------- |
| .. [4] 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. |
| |
| """ |
| if A.shape[0] == 0: |
| |
| |
| sparse, lstsq, sym_pos, cholesky = False, False, True, False |
| n_x = len(x) |
|
|
| |
| r_P = b * tau - A.dot(x) |
| r_D = c * tau - A.T.dot(y) - z |
| r_G = c.dot(x) - b.transpose().dot(y) + kappa |
| mu = (x.dot(z) + tau * kappa) / (n_x + 1) |
|
|
| |
| Dinv = x / z |
|
|
| if sparse: |
| M = A.dot(sps.diags(Dinv, 0, format="csc").dot(A.T)) |
| else: |
| M = A.dot(Dinv.reshape(-1, 1) * A.T) |
| solve = _get_solver(M, sparse, lstsq, sym_pos, cholesky, permc_spec) |
|
|
| |
| |
| |
| n_corrections = 1 if pc else 0 |
|
|
| i = 0 |
| alpha, d_x, d_z, d_tau, d_kappa = 0, 0, 0, 0, 0 |
| while i <= n_corrections: |
| |
| rhatp = eta(gamma) * r_P |
| rhatd = eta(gamma) * r_D |
| rhatg = eta(gamma) * r_G |
|
|
| |
| rhatxs = gamma * mu - x * z |
| rhattk = gamma * mu - tau * kappa |
|
|
| if i == 1: |
| if ip: |
| |
| rhatxs = ((1 - alpha) * gamma * mu - |
| x * z - alpha**2 * d_x * d_z) |
| rhattk = ((1 - alpha) * gamma * mu - |
| tau * kappa - |
| alpha**2 * d_tau * d_kappa) |
| else: |
| |
| rhatxs -= d_x * d_z |
| rhattk -= d_tau * d_kappa |
|
|
| |
| |
| |
| |
| |
| |
| |
| |
| |
| |
| |
| |
| solved = False |
| while not solved: |
| try: |
| |
| p, q = _sym_solve(Dinv, A, c, b, solve) |
| |
| u, v = _sym_solve(Dinv, A, rhatd - |
| (1 / x) * rhatxs, rhatp, solve) |
| if np.any(np.isnan(p)) or np.any(np.isnan(q)): |
| raise LinAlgError |
| solved = True |
| except (LinAlgError, ValueError, TypeError) as e: |
| |
| |
| |
| if cholesky: |
| cholesky = False |
| warn( |
| "Solving system with option 'cholesky':True " |
| "failed. It is normal for this to happen " |
| "occasionally, especially as the solution is " |
| "approached. However, if you see this frequently, " |
| "consider setting option 'cholesky' to False.", |
| OptimizeWarning, stacklevel=5) |
| elif sym_pos: |
| sym_pos = False |
| warn( |
| "Solving system with option 'sym_pos':True " |
| "failed. It is normal for this to happen " |
| "occasionally, especially as the solution is " |
| "approached. However, if you see this frequently, " |
| "consider setting option 'sym_pos' to False.", |
| OptimizeWarning, stacklevel=5) |
| elif not lstsq: |
| lstsq = True |
| warn( |
| "Solving system with option 'sym_pos':False " |
| "failed. This may happen occasionally, " |
| "especially as the solution is " |
| "approached. However, if you see this frequently, " |
| "your problem may be numerically challenging. " |
| "If you cannot improve the formulation, consider " |
| "setting 'lstsq' to True. Consider also setting " |
| "`presolve` to True, if it is not already.", |
| OptimizeWarning, stacklevel=5) |
| else: |
| raise e |
| solve = _get_solver(M, sparse, lstsq, sym_pos, |
| cholesky, permc_spec) |
| |
| d_tau = ((rhatg + 1 / tau * rhattk - (-c.dot(u) + b.dot(v))) / |
| (1 / tau * kappa + (-c.dot(p) + b.dot(q)))) |
| d_x = u + p * d_tau |
| d_y = v + q * d_tau |
|
|
| |
| d_z = (1 / x) * (rhatxs - z * d_x) |
| d_kappa = 1 / tau * (rhattk - kappa * d_tau) |
|
|
| |
| alpha = _get_step(x, d_x, z, d_z, tau, d_tau, kappa, d_kappa, 1) |
| if ip: |
| gamma = 10 |
| else: |
| beta1 = 0.1 |
| gamma = (1 - alpha)**2 * min(beta1, (1 - alpha)) |
| i += 1 |
|
|
| return d_x, d_y, d_z, d_tau, d_kappa |
|
|
|
|
| def _sym_solve(Dinv, A, r1, r2, solve): |
| """ |
| An implementation of [4] equation 8.31 and 8.32 |
| |
| References |
| ---------- |
| .. [4] 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. |
| |
| """ |
| |
| r = r2 + A.dot(Dinv * r1) |
| v = solve(r) |
| |
| u = Dinv * (A.T.dot(v) - r1) |
| return u, v |
|
|
|
|
| def _get_step(x, d_x, z, d_z, tau, d_tau, kappa, d_kappa, alpha0): |
| """ |
| An implementation of [4] equation 8.21 |
| |
| References |
| ---------- |
| .. [4] 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. |
| |
| """ |
| |
| |
| |
| |
| i_x = d_x < 0 |
| i_z = d_z < 0 |
| alpha_x = alpha0 * np.min(x[i_x] / -d_x[i_x]) if np.any(i_x) else 1 |
| alpha_tau = alpha0 * tau / -d_tau if d_tau < 0 else 1 |
| alpha_z = alpha0 * np.min(z[i_z] / -d_z[i_z]) if np.any(i_z) else 1 |
| alpha_kappa = alpha0 * kappa / -d_kappa if d_kappa < 0 else 1 |
| alpha = np.min([1, alpha_x, alpha_tau, alpha_z, alpha_kappa]) |
| return alpha |
|
|
|
|
| def _get_message(status): |
| """ |
| Given problem status code, return a more detailed message. |
| |
| Parameters |
| ---------- |
| status : int |
| 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 |
| |
| Returns |
| ------- |
| message : str |
| A string descriptor of the exit status of the optimization. |
| |
| """ |
| messages = ( |
| ["Optimization terminated successfully.", |
| "The iteration limit was reached before the algorithm converged.", |
| "The algorithm terminated successfully and determined that the " |
| "problem is infeasible.", |
| "The algorithm terminated successfully and determined that the " |
| "problem is unbounded.", |
| "Numerical difficulties were encountered before the problem " |
| "converged. Please check your problem formulation for errors, " |
| "independence of linear equality constraints, and reasonable " |
| "scaling and matrix condition numbers. If you continue to " |
| "encounter this error, please submit a bug report." |
| ]) |
| return messages[status] |
|
|
|
|
| def _do_step(x, y, z, tau, kappa, d_x, d_y, d_z, d_tau, d_kappa, alpha): |
| """ |
| An implementation of [4] Equation 8.9 |
| |
| References |
| ---------- |
| .. [4] 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. |
| |
| """ |
| x = x + alpha * d_x |
| tau = tau + alpha * d_tau |
| z = z + alpha * d_z |
| kappa = kappa + alpha * d_kappa |
| y = y + alpha * d_y |
| return x, y, z, tau, kappa |
|
|
|
|
| def _get_blind_start(shape): |
| """ |
| Return the starting point from [4] 4.4 |
| |
| References |
| ---------- |
| .. [4] 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. |
| |
| """ |
| m, n = shape |
| x0 = np.ones(n) |
| y0 = np.zeros(m) |
| z0 = np.ones(n) |
| tau0 = 1 |
| kappa0 = 1 |
| return x0, y0, z0, tau0, kappa0 |
|
|
|
|
| def _indicators(A, b, c, c0, x, y, z, tau, kappa): |
| """ |
| Implementation of several equations from [4] used as indicators of |
| the status of optimization. |
| |
| References |
| ---------- |
| .. [4] 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. |
| |
| """ |
|
|
| |
| x0, y0, z0, tau0, kappa0 = _get_blind_start(A.shape) |
|
|
| |
| def r_p(x, tau): |
| return b * tau - A.dot(x) |
|
|
| def r_d(y, z, tau): |
| return c * tau - A.T.dot(y) - z |
|
|
| def r_g(x, y, kappa): |
| return kappa + c.dot(x) - b.dot(y) |
|
|
| |
| def mu(x, tau, z, kappa): |
| return (x.dot(z) + np.dot(tau, kappa)) / (len(x) + 1) |
|
|
| obj = c.dot(x / tau) + c0 |
|
|
| def norm(a): |
| return np.linalg.norm(a) |
|
|
| |
| r_p0 = r_p(x0, tau0) |
| r_d0 = r_d(y0, z0, tau0) |
| r_g0 = r_g(x0, y0, kappa0) |
| mu_0 = mu(x0, tau0, z0, kappa0) |
| rho_A = norm(c.T.dot(x) - b.T.dot(y)) / (tau + norm(b.T.dot(y))) |
| rho_p = norm(r_p(x, tau)) / max(1, norm(r_p0)) |
| rho_d = norm(r_d(y, z, tau)) / max(1, norm(r_d0)) |
| rho_g = norm(r_g(x, y, kappa)) / max(1, norm(r_g0)) |
| rho_mu = mu(x, tau, z, kappa) / mu_0 |
| return rho_p, rho_d, rho_A, rho_g, rho_mu, obj |
|
|
|
|
| def _display_iter(rho_p, rho_d, rho_g, alpha, rho_mu, obj, header=False): |
| """ |
| Print indicators of optimization status to the console. |
| |
| Parameters |
| ---------- |
| rho_p : float |
| The (normalized) primal feasibility, see [4] 4.5 |
| rho_d : float |
| The (normalized) dual feasibility, see [4] 4.5 |
| rho_g : float |
| The (normalized) duality gap, see [4] 4.5 |
| alpha : float |
| The step size, see [4] 4.3 |
| rho_mu : float |
| The (normalized) path parameter, see [4] 4.5 |
| obj : float |
| The objective function value of the current iterate |
| header : bool |
| True if a header is to be printed |
| |
| References |
| ---------- |
| .. [4] 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. |
| |
| """ |
| if header: |
| print("Primal Feasibility ", |
| "Dual Feasibility ", |
| "Duality Gap ", |
| "Step ", |
| "Path Parameter ", |
| "Objective ") |
|
|
| |
| fmt = '{0:<20.13}{1:<20.13}{2:<20.13}{3:<17.13}{4:<20.13}{5:<20.13}' |
| print(fmt.format( |
| float(rho_p), |
| float(rho_d), |
| float(rho_g), |
| alpha if isinstance(alpha, str) else float(alpha), |
| float(rho_mu), |
| float(obj))) |
|
|
|
|
| def _ip_hsd(A, b, c, c0, alpha0, beta, maxiter, disp, tol, sparse, lstsq, |
| sym_pos, cholesky, pc, ip, permc_spec, callback, postsolve_args): |
| r""" |
| Solve a linear programming problem in standard form: |
| |
| Minimize:: |
| |
| c @ x |
| |
| Subject to:: |
| |
| A @ x == b |
| x >= 0 |
| |
| using the interior point method of [4]. |
| |
| Parameters |
| ---------- |
| A : 2-D array |
| 2-D array such that ``A @ x``, gives the values of the equality |
| constraints at ``x``. |
| b : 1-D array |
| 1-D array of values representing the RHS of each equality constraint |
| (row) in ``A`` (for standard form problem). |
| c : 1-D array |
| Coefficients of the linear objective function to be minimized (for |
| standard form problem). |
| c0 : float |
| Constant term in objective function due to fixed (and eliminated) |
| variables. (Purely for display.) |
| alpha0 : float |
| The maximal step size for Mehrota's predictor-corrector search |
| direction; see :math:`\beta_3`of [4] Table 8.1 |
| beta : float |
| The desired reduction of the path parameter :math:`\mu` (see [6]_) |
| maxiter : int |
| The maximum number of iterations of the algorithm. |
| disp : bool |
| Set to ``True`` if indicators of optimization status are to be printed |
| to the console each iteration. |
| tol : float |
| Termination tolerance; see [4]_ Section 4.5. |
| sparse : bool |
| Set to ``True`` if the problem is to be treated as sparse. However, |
| the inputs ``A_eq`` and ``A_ub`` should nonetheless be provided as |
| (dense) arrays rather than sparse matrices. |
| lstsq : bool |
| Set to ``True`` if the problem is expected to be very poorly |
| conditioned. This should always be left as ``False`` unless severe |
| numerical difficulties are frequently encountered, and a better option |
| would be to improve the formulation of the problem. |
| sym_pos : bool |
| Leave ``True`` if the problem is expected to yield a well conditioned |
| symmetric positive definite normal equation matrix (almost always). |
| cholesky : bool |
| Set to ``True`` if the normal equations are to be solved by explicit |
| Cholesky decomposition followed by explicit forward/backward |
| substitution. This is typically faster for moderate, dense problems |
| that are numerically well-behaved. |
| pc : bool |
| Leave ``True`` if the predictor-corrector method of Mehrota is to be |
| used. This is almost always (if not always) beneficial. |
| ip : bool |
| Set to ``True`` if the improved initial point suggestion due to [4]_ |
| Section 4.3 is desired. It's unclear whether this is beneficial. |
| permc_spec : str (default = 'MMD_AT_PLUS_A') |
| (Has effect only with ``sparse = True``, ``lstsq = False``, ``sym_pos = |
| True``.) 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``. |
| callback : callable, optional |
| If a callback function is provided, it will be called within each |
| iteration of the algorithm. The callback function must accept a single |
| `scipy.optimize.OptimizeResult` consisting of the following fields: |
| |
| x : 1-D array |
| Current solution vector |
| fun : float |
| Current value of the objective function |
| success : bool |
| True only when an algorithm has completed successfully, |
| so this is always False as the callback function is called |
| only while the algorithm is still iterating. |
| slack : 1-D array |
| The values of the slack variables. Each slack variable |
| corresponds to an inequality constraint. If the slack is zero, |
| the corresponding constraint is active. |
| con : 1-D array |
| The (nominally zero) residuals of the equality constraints, |
| that is, ``b - A_eq @ x`` |
| phase : int |
| The phase of the algorithm being executed. This is always |
| 1 for the interior-point method because it has only one phase. |
| status : int |
| For revised simplex, this is always 0 because if a different |
| status is detected, the algorithm terminates. |
| nit : int |
| The number of iterations performed. |
| message : str |
| A string descriptor of the exit status of the optimization. |
| postsolve_args : tuple |
| Data needed by _postsolve to convert the solution to the standard-form |
| problem into the solution to the original problem. |
| |
| Returns |
| ------- |
| x_hat : float |
| Solution vector (for standard form problem). |
| status : int |
| 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 |
| |
| message : str |
| A string descriptor of the exit status of the optimization. |
| iteration : int |
| The number of iterations taken to solve the problem |
| |
| References |
| ---------- |
| .. [4] 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] 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 |
| |
| """ |
|
|
| iteration = 0 |
|
|
| |
| x, y, z, tau, kappa = _get_blind_start(A.shape) |
|
|
| |
| ip = ip if pc else False |
|
|
| |
| rho_p, rho_d, rho_A, rho_g, rho_mu, obj = _indicators( |
| A, b, c, c0, x, y, z, tau, kappa) |
| go = rho_p > tol or rho_d > tol or rho_A > tol |
|
|
| if disp: |
| _display_iter(rho_p, rho_d, rho_g, "-", rho_mu, obj, header=True) |
| if callback is not None: |
| x_o, fun, slack, con = _postsolve(x/tau, postsolve_args) |
| res = OptimizeResult({'x': x_o, 'fun': fun, 'slack': slack, |
| 'con': con, 'nit': iteration, 'phase': 1, |
| 'complete': False, 'status': 0, |
| 'message': "", 'success': False}) |
| callback(res) |
|
|
| status = 0 |
| message = "Optimization terminated successfully." |
|
|
| if sparse: |
| A = sps.csc_matrix(A) |
|
|
| while go: |
|
|
| iteration += 1 |
|
|
| if ip: |
| |
| gamma = 1 |
|
|
| def eta(g): |
| return 1 |
| else: |
| |
| |
| |
| gamma = 0 if pc else beta * np.mean(z * x) |
| |
|
|
| def eta(g=gamma): |
| return 1 - g |
|
|
| try: |
| |
| d_x, d_y, d_z, d_tau, d_kappa = _get_delta( |
| A, b, c, x, y, z, tau, kappa, gamma, eta, |
| sparse, lstsq, sym_pos, cholesky, pc, ip, permc_spec) |
|
|
| if ip: |
| |
| |
| |
| alpha = 1.0 |
| x, y, z, tau, kappa = _do_step( |
| x, y, z, tau, kappa, d_x, d_y, |
| d_z, d_tau, d_kappa, alpha) |
| x[x < 1] = 1 |
| z[z < 1] = 1 |
| tau = max(1, tau) |
| kappa = max(1, kappa) |
| ip = False |
| else: |
| |
| alpha = _get_step(x, d_x, z, d_z, tau, |
| d_tau, kappa, d_kappa, alpha0) |
| |
| x, y, z, tau, kappa = _do_step( |
| x, y, z, tau, kappa, d_x, d_y, d_z, d_tau, d_kappa, alpha) |
|
|
| except (LinAlgError, FloatingPointError, |
| ValueError, ZeroDivisionError): |
| |
| |
| |
| status = 4 |
| message = _get_message(status) |
| break |
|
|
| |
| rho_p, rho_d, rho_A, rho_g, rho_mu, obj = _indicators( |
| A, b, c, c0, x, y, z, tau, kappa) |
| go = rho_p > tol or rho_d > tol or rho_A > tol |
|
|
| if disp: |
| _display_iter(rho_p, rho_d, rho_g, alpha, rho_mu, obj) |
| if callback is not None: |
| x_o, fun, slack, con = _postsolve(x/tau, postsolve_args) |
| res = OptimizeResult({'x': x_o, 'fun': fun, 'slack': slack, |
| 'con': con, 'nit': iteration, 'phase': 1, |
| 'complete': False, 'status': 0, |
| 'message': "", 'success': False}) |
| callback(res) |
|
|
| |
| inf1 = (rho_p < tol and rho_d < tol and rho_g < tol and tau < tol * |
| max(1, kappa)) |
| inf2 = rho_mu < tol and tau < tol * min(1, kappa) |
| if inf1 or inf2: |
| |
| if b.transpose().dot(y) > tol: |
| status = 2 |
| else: |
| status = 3 |
| message = _get_message(status) |
| break |
| elif iteration >= maxiter: |
| status = 1 |
| message = _get_message(status) |
| break |
|
|
| x_hat = x / tau |
| |
| return x_hat, status, message, iteration |
|
|
|
|
| def _linprog_ip(c, c0, A, b, callback, postsolve_args, maxiter=1000, tol=1e-8, |
| disp=False, alpha0=.99995, beta=0.1, sparse=False, lstsq=False, |
| sym_pos=True, cholesky=None, pc=True, ip=False, |
| permc_spec='MMD_AT_PLUS_A', **unknown_options): |
| r""" |
| 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 |
| |
| User-facing documentation is in _linprog_doc.py. |
| |
| Parameters |
| ---------- |
| c : 1-D array |
| Coefficients of the linear objective function to be minimized. |
| c0 : float |
| Constant term in objective function due to fixed (and eliminated) |
| variables. (Purely for display.) |
| A : 2-D array |
| 2-D array such that ``A @ x``, gives the values of the equality |
| constraints at ``x``. |
| b : 1-D array |
| 1-D array of values representing the right hand side of each equality |
| constraint (row) in ``A``. |
| callback : callable, optional |
| Callback function to be executed once per iteration. |
| postsolve_args : tuple |
| Data needed by _postsolve to convert the solution to the standard-form |
| problem into the solution to the original problem. |
| |
| Options |
| ------- |
| maxiter : int (default = 1000) |
| The maximum number of iterations of the algorithm. |
| tol : float (default = 1e-8) |
| Termination tolerance to be used for all termination criteria; |
| see [4]_ Section 4.5. |
| disp : bool (default = False) |
| Set to ``True`` if indicators of optimization status are to be printed |
| to the console each iteration. |
| alpha0 : float (default = 0.99995) |
| The maximal step size for Mehrota's predictor-corrector search |
| direction; see :math:`\beta_{3}` of [4]_ Table 8.1. |
| beta : float (default = 0.1) |
| The desired reduction of the path parameter :math:`\mu` (see [6]_) |
| when Mehrota's predictor-corrector is not in use (uncommon). |
| sparse : bool (default = False) |
| Set to ``True`` if the problem is to be treated as sparse after |
| presolve. If either ``A_eq`` or ``A_ub`` is 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 ``True`` |
| or providing ``A_eq`` and ``A_ub`` as sparse matrices. |
| lstsq : bool (default = False) |
| Set to ``True`` if the problem is expected to be very poorly |
| conditioned. This should always be left ``False`` unless severe |
| numerical difficulties are encountered. Leave this at the default |
| unless you receive a warning message suggesting otherwise. |
| sym_pos : bool (default = True) |
| Leave ``True`` if 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. |
| cholesky : bool (default = True) |
| Set to ``True`` if 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. |
| pc : bool (default = True) |
| Leave ``True`` if the predictor-corrector method of Mehrota is to be |
| used. This is almost always (if not always) beneficial. |
| ip : bool (default = False) |
| Set to ``True`` if the improved initial point suggestion due to [4]_ |
| Section 4.3 is desired. Whether this is beneficial or not |
| depends on the problem. |
| permc_spec : str (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_options : dict |
| Optional arguments not used by this particular solver. If |
| `unknown_options` is non-empty a warning is issued listing all |
| unused options. |
| |
| Returns |
| ------- |
| x : 1-D array |
| Solution vector. |
| status : int |
| 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 |
| |
| message : str |
| A string descriptor of the exit status of the optimization. |
| iteration : int |
| The number of iterations taken to solve the problem. |
| |
| 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: |
| |
| 1. ``scipy.linalg.cho_factor`` |
| |
| 2. ``scipy.linalg.solve`` with option ``sym_pos=True`` |
| |
| 3. ``scipy.linalg.solve`` with option ``sym_pos=False`` |
| |
| 4. ``scipy.linalg.lstsq`` |
| |
| For sparse problems: |
| |
| 1. ``sksparse.cholmod.cholesky`` (if scikit-sparse and SuiteSparse are installed) |
| |
| 2. ``scipy.sparse.linalg.factorized`` |
| (if scikit-umfpack and SuiteSparse are installed) |
| |
| 3. ``scipy.sparse.linalg.splu`` (which uses SuperLU distributed with SciPy) |
| |
| 4. ``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=False`` skips to solver 2, |
| ``sym_pos=False`` skips to solver 3, and ``lstsq=True`` skips |
| to solver 4 for both sparse and dense problems. |
| |
| Potential improvements for combating 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 ``tol`` option, |
| 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 ``linprog`` |
| module 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 ``linprog`` module 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 = 0`` and ``ub = None`` unless 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. |
| |
| ``linprog`` module 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] 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] 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. |
| |
| """ |
|
|
| _check_unknown_options(unknown_options) |
|
|
| |
| if (cholesky or cholesky is None) and sparse and not has_cholmod: |
| if cholesky: |
| warn("Sparse cholesky is only available with scikit-sparse. " |
| "Setting `cholesky = False`", |
| OptimizeWarning, stacklevel=3) |
| cholesky = False |
|
|
| if sparse and lstsq: |
| warn("Option combination 'sparse':True and 'lstsq':True " |
| "is not recommended.", |
| OptimizeWarning, stacklevel=3) |
|
|
| if lstsq and cholesky: |
| warn("Invalid option combination 'lstsq':True " |
| "and 'cholesky':True; option 'cholesky' has no effect when " |
| "'lstsq' is set True.", |
| OptimizeWarning, stacklevel=3) |
|
|
| valid_permc_spec = ('NATURAL', 'MMD_ATA', 'MMD_AT_PLUS_A', 'COLAMD') |
| if permc_spec.upper() not in valid_permc_spec: |
| warn("Invalid permc_spec option: '" + str(permc_spec) + "'. " |
| "Acceptable values are 'NATURAL', 'MMD_ATA', 'MMD_AT_PLUS_A', " |
| "and 'COLAMD'. Reverting to default.", |
| OptimizeWarning, stacklevel=3) |
| permc_spec = 'MMD_AT_PLUS_A' |
|
|
| |
| if not sym_pos and cholesky: |
| raise ValueError( |
| "Invalid option combination 'sym_pos':False " |
| "and 'cholesky':True: Cholesky decomposition is only possible " |
| "for symmetric positive definite matrices.") |
|
|
| cholesky = cholesky or (cholesky is None and sym_pos and not lstsq) |
|
|
| x, status, message, iteration = _ip_hsd(A, b, c, c0, alpha0, beta, |
| maxiter, disp, tol, sparse, |
| lstsq, sym_pos, cholesky, |
| pc, ip, permc_spec, callback, |
| postsolve_args) |
|
|
| return x, status, message, iteration |
|
|