1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069 |
- """
- An interior-point method for linear programming.
- """
- # Author: Matt Haberland
- from __future__ import print_function, division, absolute_import
- 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, _check_unknown_options
- def _get_solver(sparse=False, lstsq=False, sym_pos=True, cholesky=True):
- """
- Given solver options, return a handle to the appropriate linear system
- solver.
- Parameters
- ----------
- 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.
- Returns
- -------
- solve : function
- Handle to the appropriate solver function
- """
- if sparse:
- if lstsq or not(sym_pos):
- def solve(M, r, sym_pos=False):
- return sps.linalg.lsqr(M, r)[0]
- else:
- # this is not currently used; it is replaced by splu solve
- # TODO: expose use of this as an option
- def solve(M, r):
- return sps.linalg.spsolve(M, r, permc_spec="MMD_AT_PLUS_A")
- else:
- if lstsq: # sometimes necessary as solution is approached
- def solve(M, r):
- return sp.linalg.lstsq(M, r)[0]
- elif cholesky:
- solve = sp.linalg.cho_solve
- else:
- # this seems to cache the matrix factorization, so solving
- # with multiple right hand sides is much faster
- def solve(M, r, sym_pos=sym_pos):
- return sp.linalg.solve(M, r, sym_pos=sym_pos)
- 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:
- # If there are no constraints, some solvers fail (understandably)
- # rather than returning empty solution. This gets the job done.
- sparse, lstsq, sym_pos, cholesky = False, False, True, False
- solve = _get_solver(sparse, lstsq, sym_pos, cholesky)
- n_x = len(x)
- # [4] Equation 8.8
- 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)
- # Assemble M from [4] Equation 8.31
- Dinv = x / z
- splu = False
- if sparse and not lstsq:
- # sparse requires Dinv to be diag matrix
- M = A.dot(sps.diags(Dinv, 0, format="csc").dot(A.T))
- try:
- # TODO: should use linalg.factorized instead, but I don't have
- # umfpack and therefore cannot test its performance
- solve = sps.linalg.splu(M, permc_spec=permc_spec).solve
- splu = True
- except Exception:
- lstsq = True
- solve = _get_solver(sparse, lstsq, sym_pos, cholesky)
- else:
- # dense does not; use broadcasting
- M = A.dot(Dinv.reshape(-1, 1) * A.T)
- # For some small problems, calling sp.linalg.solve w/ sym_pos = True
- # may be faster. I am pretty certain it caches the factorization for
- # multiple uses and checks the incoming matrix to see if it's the same as
- # the one it already factorized. (I can't explain the speed otherwise.)
- if cholesky:
- try:
- L = sp.linalg.cho_factor(M)
- except Exception:
- cholesky = False
- solve = _get_solver(sparse, lstsq, sym_pos, cholesky)
- # pc: "predictor-corrector" [4] Section 4.1
- # In development this option could be turned off
- # but it always seems to improve performance substantially
- 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:
- # Reference [4] Eq. 8.6
- rhatp = eta(gamma) * r_P
- rhatd = eta(gamma) * r_D
- rhatg = np.array(eta(gamma) * r_G).reshape((1,))
- # Reference [4] Eq. 8.7
- rhatxs = gamma * mu - x * z
- rhattk = np.array(gamma * mu - tau * kappa).reshape((1,))
- if i == 1:
- if ip: # if the correction is to get "initial point"
- # Reference [4] Eq. 8.23
- rhatxs = ((1 - alpha) * gamma * mu -
- x * z - alpha**2 * d_x * d_z)
- rhattk = np.array(
- (1 -
- alpha) *
- gamma *
- mu -
- tau *
- kappa -
- alpha**2 *
- d_tau *
- d_kappa).reshape(
- (1,
- ))
- else: # if the correction is for "predictor-corrector"
- # Reference [4] Eq. 8.13
- rhatxs -= d_x * d_z
- rhattk -= d_tau * d_kappa
- # sometimes numerical difficulties arise as the solution is approached
- # this loop tries to solve the equations using a sequence of functions
- # for solve. For dense systems, the order is:
- # 1. scipy.linalg.cho_factor/scipy.linalg.cho_solve,
- # 2. scipy.linalg.solve w/ sym_pos = True,
- # 3. scipy.linalg.solve w/ sym_pos = False, and if all else fails
- # 4. scipy.linalg.lstsq
- # For sparse systems, the order is:
- # 1. scipy.sparse.linalg.splu
- # 2. scipy.sparse.linalg.lsqr
- # TODO: if umfpack is installed, use factorized instead of splu.
- # Can't do that now because factorized doesn't pass permc_spec
- # to splu if umfpack isn't installed. Also, umfpack not tested.
- solved = False
- while(not solved):
- try:
- solve_this = L if cholesky else M
- # [4] Equation 8.28
- p, q = _sym_solve(Dinv, solve_this, A, c, b, solve, splu)
- # [4] Equation 8.29
- u, v = _sym_solve(Dinv, solve_this, A, rhatd -
- (1 / x) * rhatxs, rhatp, solve, splu)
- if np.any(np.isnan(p)) or np.any(np.isnan(q)):
- raise LinAlgError
- solved = True
- except (LinAlgError, ValueError) as e:
- # Usually this doesn't happen. If it does, it happens when
- # there are redundant constraints or when approaching the
- # solution. If so, change solver.
- cholesky = False
- if not lstsq:
- if sym_pos:
- 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)
- sym_pos = False
- else:
- 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)
- lstsq = True
- else:
- raise e
- solve = _get_solver(sparse, lstsq, sym_pos)
- # [4] Results after 8.29
- 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
- # [4] Relations between after 8.25 and 8.26
- d_z = (1 / x) * (rhatxs - z * d_x)
- d_kappa = 1 / tau * (rhattk - kappa * d_tau)
- # [4] 8.12 and "Let alpha be the maximal possible step..." before 8.23
- alpha = _get_step(x, d_x, z, d_z, tau, d_tau, kappa, d_kappa, 1)
- if ip: # initial point - see [4] 4.4
- gamma = 10
- else: # predictor-corrector, [4] definition after 8.12
- beta1 = 0.1 # [4] pg. 220 (Table 8.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, M, A, r1, r2, solve, splu=False):
- """
- 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.
- """
- # [4] 8.31
- r = r2 + A.dot(Dinv * r1)
- if splu:
- v = solve(r)
- else:
- v = solve(M, r)
- # [4] 8.32
- 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.
- """
- # [4] 4.3 Equation 8.21, ignoring 8.20 requirement
- # same step is taken in primal and dual spaces
- # alpha0 is basically beta3 from [4] Table 8.1, but instead of beta3
- # the value 1 is used in Mehrota corrector and initial point correction
- 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.
- """
- # residuals for termination are relative to initial values
- x0, y0, z0, tau0, kappa0 = _get_blind_start(A.shape)
- # See [4], Section 4 - The Homogeneous Algorithm, Equation 8.8
- 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)
- # np.dot unpacks if they are arrays of size one
- 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)
- # See [4], Section 4.5 - The Stopping Criteria
- 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 ")
- # no clue why this works
- fmt = '{0:<20.13}{1:<20.13}{2:<20.13}{3:<17.13}{4:<20.13}{5:<20.13}'
- print(fmt.format(
- rho_p,
- rho_d,
- rho_g,
- alpha,
- rho_mu,
- obj))
- def _ip_hsd(A, b, c, c0, alpha0, beta, maxiter, disp, tol,
- sparse, lstsq, sym_pos, cholesky, pc, ip, permc_spec):
- 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 : 2D array
- 2D array such that ``A @ x``, gives the values of the equality
- constraints at ``x``.
- b : 1D array
- 1D array of values representing the RHS of each equality constraint
- (row) in ``A`` (for standard form problem).
- c : 1D 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``.
- 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
- # default initial point
- x, y, z, tau, kappa = _get_blind_start(A.shape)
- # first iteration is special improvement of initial point
- ip = ip if pc else False
- # [4] 4.5
- 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 # we might get lucky : )
- if disp:
- _display_iter(rho_p, rho_d, rho_g, "-", rho_mu, obj, header=True)
- status = 0
- message = "Optimization terminated successfully."
- if sparse:
- A = sps.csc_matrix(A)
- A.T = A.transpose() # A.T is defined for sparse matrices but is slow
- # Redefine it to avoid calculating again
- # This is fine as long as A doesn't change
- while go:
- iteration += 1
- if ip: # initial point
- # [4] Section 4.4
- gamma = 1
- def eta(g):
- return 1
- else:
- # gamma = 0 in predictor step according to [4] 4.1
- # if predictor/corrector is off, use mean of complementarity [6]
- # 5.1 / [4] Below Figure 10-4
- gamma = 0 if pc else beta * np.mean(z * x)
- # [4] Section 4.1
- def eta(g=gamma):
- return 1 - g
- try:
- # Solve [4] 8.6 and 8.7/8.13/8.23
- 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: # initial point
- # [4] 4.4
- # Formula after 8.23 takes a full step regardless if this will
- # take it negative
- 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 # done with initial point
- else:
- # [4] Section 4.3
- alpha = _get_step(x, d_x, z, d_z, tau,
- d_tau, kappa, d_kappa, alpha0)
- # [4] Equation 8.9
- 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):
- # this can happen when sparse solver is used and presolve
- # is turned off. Also observed ValueError in AppVeyor Python 3.6
- # Win32 build (PR #8676). I've never seen it otherwise.
- status = 4
- message = _get_message(status)
- break
- # [4] 4.5
- 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, float(rho_mu), obj)
- # [4] 4.5
- 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:
- # [4] Lemma 8.4 / Theorem 8.3
- if b.transpose().dot(y) > tol:
- status = 2
- else: # elif c.T.dot(x) < tol: ? Probably not necessary.
- status = 3
- message = _get_message(status)
- break
- elif iteration >= maxiter:
- status = 1
- message = _get_message(status)
- break
- x_hat = x / tau
- # [4] Statement after Theorem 8.2
- return x_hat, status, message, iteration
- def _linprog_ip(
- c,
- c0=0,
- A=None,
- b=None,
- callback=None,
- alpha0=.99995,
- beta=0.1,
- maxiter=1000,
- disp=False,
- tol=1e-8,
- 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
- Parameters
- ----------
- c : 1D 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 : 2D array
- 2D array such that ``A @ x``, gives the values of the equality
- constraints at ``x``.
- b : 1D array
- 1D array of values representing the right hand side of each equality
- constraint (row) in ``A``.
- Options
- -------
- maxiter : int (default = 1000)
- The maximum number of iterations of the algorithm.
- disp : bool (default = False)
- Set to ``True`` if indicators of optimization status are to be printed
- to the console each iteration.
- tol : float (default = 1e-8)
- Termination tolerance to be used for all termination criteria;
- see [4]_ Section 4.5.
- 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 moderate, dense 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``.) 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
- -------
- x : 1D 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]_ and [4]_.
- 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 the default ``cholesky=True``, this is accomplished using
- ``scipy.linalg.cho_factor`` followed by forward/backward substitutions
- via ``scipy.linalg.cho_solve``. With ``cholesky=False`` and
- ``sym_pos=True``, Cholesky decomposition is performed instead by
- ``scipy.linalg.solve``. Based on speed tests, this also appears to retain
- the Cholesky decomposition of the matrix for later use, which is beneficial
- as the same system is solved four times with different right hand sides
- in each iteration of the algorithm.
- In problems with redundancy (e.g. if presolve is turned off with option
- ``presolve=False``) or if the matrices become ill-conditioned (e.g. as the
- solution is approached and some decision variables approach zero),
- Cholesky decomposition can fail. Should this occur, successively more
- robust solvers (``scipy.linalg.solve`` with ``sym_pos=False`` then
- ``scipy.linalg.lstsq``) are tried, at the cost of computational efficiency.
- These solvers can be used from the outset by setting the options
- ``sym_pos=False`` and ``lstsq=True``, respectively.
- Note that with the option ``sparse=True``, the normal equations are solved
- using ``scipy.sparse.linalg.spsolve``. Unfortunately, this uses the more
- expensive LU decomposition from the outset, but for large, sparse problems,
- the use of sparse linear algebra techniques improves the solve speed
- despite the use of LU rather than Cholesky decomposition. A simple
- improvement would be to use the sparse Cholesky decomposition of
- ``CHOLMOD`` via ``scikit-sparse`` when available.
- Other 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 ``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 callback is not None:
- raise NotImplementedError("method 'interior-point' does not support "
- "callback functions.")
- # These should be warnings, not errors
- if sparse and lstsq:
- warn("Invalid option combination 'sparse':True "
- "and 'lstsq':True; Sparse least squares is not recommended.",
- OptimizeWarning)
- if sparse and not sym_pos:
- warn("Invalid option combination 'sparse':True "
- "and 'sym_pos':False; the effect is the same as sparse least "
- "squares, which is not recommended.",
- OptimizeWarning)
- if sparse and cholesky:
- # Cholesky decomposition is not available for sparse problems
- warn("Invalid option combination 'sparse':True "
- "and 'cholesky':True; sparse Colesky decomposition is not "
- "available.",
- OptimizeWarning)
- if lstsq and cholesky:
- warn("Invalid option combination 'lstsq':True "
- "and 'cholesky':True; option 'cholesky' has no effect when "
- "'lstsq' is set True.",
- OptimizeWarning)
- 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)
- permc_spec = 'MMD_AT_PLUS_A'
- # This can be an error
- 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 is None and sym_pos and not sparse 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)
- return x, status, message, iteration
|