_linprog_ip.py 46 KB

1234567891011121314151617181920212223242526272829303132333435363738394041424344454647484950515253545556575859606162636465666768697071727374757677787980818283848586878889909192939495969798991001011021031041051061071081091101111121131141151161171181191201211221231241251261271281291301311321331341351361371381391401411421431441451461471481491501511521531541551561571581591601611621631641651661671681691701711721731741751761771781791801811821831841851861871881891901911921931941951961971981992002012022032042052062072082092102112122132142152162172182192202212222232242252262272282292302312322332342352362372382392402412422432442452462472482492502512522532542552562572582592602612622632642652662672682692702712722732742752762772782792802812822832842852862872882892902912922932942952962972982993003013023033043053063073083093103113123133143153163173183193203213223233243253263273283293303313323333343353363373383393403413423433443453463473483493503513523533543553563573583593603613623633643653663673683693703713723733743753763773783793803813823833843853863873883893903913923933943953963973983994004014024034044054064074084094104114124134144154164174184194204214224234244254264274284294304314324334344354364374384394404414424434444454464474484494504514524534544554564574584594604614624634644654664674684694704714724734744754764774784794804814824834844854864874884894904914924934944954964974984995005015025035045055065075085095105115125135145155165175185195205215225235245255265275285295305315325335345355365375385395405415425435445455465475485495505515525535545555565575585595605615625635645655665675685695705715725735745755765775785795805815825835845855865875885895905915925935945955965975985996006016026036046056066076086096106116126136146156166176186196206216226236246256266276286296306316326336346356366376386396406416426436446456466476486496506516526536546556566576586596606616626636646656666676686696706716726736746756766776786796806816826836846856866876886896906916926936946956966976986997007017027037047057067077087097107117127137147157167177187197207217227237247257267277287297307317327337347357367377387397407417427437447457467477487497507517527537547557567577587597607617627637647657667677687697707717727737747757767777787797807817827837847857867877887897907917927937947957967977987998008018028038048058068078088098108118128138148158168178188198208218228238248258268278288298308318328338348358368378388398408418428438448458468478488498508518528538548558568578588598608618628638648658668678688698708718728738748758768778788798808818828838848858868878888898908918928938948958968978988999009019029039049059069079089099109119129139149159169179189199209219229239249259269279289299309319329339349359369379389399409419429439449459469479489499509519529539549559569579589599609619629639649659669679689699709719729739749759769779789799809819829839849859869879889899909919929939949959969979989991000100110021003100410051006100710081009101010111012101310141015101610171018101910201021102210231024102510261027102810291030103110321033103410351036103710381039104010411042104310441045104610471048104910501051105210531054105510561057105810591060106110621063106410651066106710681069107010711072107310741075107610771078107910801081108210831084108510861087108810891090109110921093109410951096109710981099110011011102110311041105110611071108110911101111111211131114111511161117111811191120112111221123112411251126112711281129113011311132113311341135113611371138113911401141
  1. """Interior-point method for linear programming
  2. The *interior-point* method uses the primal-dual path following algorithm
  3. outlined in [1]_. This algorithm supports sparse constraint matrices and
  4. is typically faster than the simplex methods, especially for large, sparse
  5. problems. Note, however, that the solution returned may be slightly less
  6. accurate than those of the simplex methods and will not, in general,
  7. correspond with a vertex of the polytope defined by the constraints.
  8. .. versionadded:: 1.0.0
  9. References
  10. ----------
  11. .. [1] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  12. optimizer for linear programming: an implementation of the
  13. homogeneous algorithm." High performance optimization. Springer US,
  14. 2000. 197-232.
  15. """
  16. # Author: Matt Haberland
  17. import numpy as np
  18. import scipy as sp
  19. import scipy.sparse as sps
  20. from warnings import warn
  21. from scipy.linalg import LinAlgError
  22. from ._optimize import OptimizeWarning, OptimizeResult, _check_unknown_options
  23. from ._linprog_util import _postsolve
  24. has_umfpack = True
  25. has_cholmod = True
  26. try:
  27. import sksparse # noqa: F401
  28. from sksparse.cholmod import cholesky as cholmod # noqa: F401
  29. from sksparse.cholmod import analyze as cholmod_analyze
  30. from sksparse.cholmod import CholmodTypeConversionWarning
  31. from warnings import catch_warnings
  32. except ImportError:
  33. has_cholmod = False
  34. try:
  35. import scikits.umfpack # test whether to use factorized # noqa: F401
  36. except ImportError:
  37. has_umfpack = False
  38. def _get_solver(M, sparse=False, lstsq=False, sym_pos=True,
  39. cholesky=True, permc_spec='MMD_AT_PLUS_A'):
  40. """
  41. Given solver options, return a handle to the appropriate linear system
  42. solver.
  43. Parameters
  44. ----------
  45. M : 2-D array
  46. As defined in [4] Equation 8.31
  47. sparse : bool (default = False)
  48. True if the system to be solved is sparse. This is typically set
  49. True when the original ``A_ub`` and ``A_eq`` arrays are sparse.
  50. lstsq : bool (default = False)
  51. True if the system is ill-conditioned and/or (nearly) singular and
  52. thus a more robust least-squares solver is desired. This is sometimes
  53. needed as the solution is approached.
  54. sym_pos : bool (default = True)
  55. True if the system matrix is symmetric positive definite
  56. Sometimes this needs to be set false as the solution is approached,
  57. even when the system should be symmetric positive definite, due to
  58. numerical difficulties.
  59. cholesky : bool (default = True)
  60. True if the system is to be solved by Cholesky, rather than LU,
  61. decomposition. This is typically faster unless the problem is very
  62. small or prone to numerical difficulties.
  63. permc_spec : str (default = 'MMD_AT_PLUS_A')
  64. Sparsity preservation strategy used by SuperLU. Acceptable values are:
  65. - ``NATURAL``: natural ordering.
  66. - ``MMD_ATA``: minimum degree ordering on the structure of A^T A.
  67. - ``MMD_AT_PLUS_A``: minimum degree ordering on the structure of A^T+A.
  68. - ``COLAMD``: approximate minimum degree column ordering.
  69. See SuperLU documentation.
  70. Returns
  71. -------
  72. solve : function
  73. Handle to the appropriate solver function
  74. """
  75. try:
  76. if sparse:
  77. if lstsq:
  78. def solve(r, sym_pos=False):
  79. return sps.linalg.lsqr(M, r)[0]
  80. elif cholesky:
  81. # TODO: revert this suppress_warning once the warning bug fix in
  82. # sksparse is merged/released
  83. # Suppress spurious warning bug from sksparse with csc_array gh-22089
  84. # try:
  85. # # Will raise an exception in the first call,
  86. # # or when the matrix changes due to a new problem
  87. # _get_solver.cholmod_factor.cholesky_inplace(M)
  88. # except Exception:
  89. # _get_solver.cholmod_factor = cholmod_analyze(M)
  90. # _get_solver.cholmod_factor.cholesky_inplace(M)
  91. with catch_warnings(
  92. action='ignore', category=CholmodTypeConversionWarning
  93. ):
  94. try:
  95. # Will raise an exception in the first call,
  96. # or when the matrix changes due to a new problem
  97. _get_solver.cholmod_factor.cholesky_inplace(M)
  98. except Exception:
  99. _get_solver.cholmod_factor = cholmod_analyze(M)
  100. _get_solver.cholmod_factor.cholesky_inplace(M)
  101. solve = _get_solver.cholmod_factor
  102. else:
  103. if has_umfpack and sym_pos:
  104. solve = sps.linalg.factorized(M)
  105. else: # factorized doesn't pass permc_spec
  106. solve = sps.linalg.splu(M, permc_spec=permc_spec).solve
  107. else:
  108. if lstsq: # sometimes necessary as solution is approached
  109. def solve(r):
  110. return sp.linalg.lstsq(M, r)[0]
  111. elif cholesky:
  112. L = sp.linalg.cho_factor(M)
  113. def solve(r):
  114. return sp.linalg.cho_solve(L, r)
  115. else:
  116. # this seems to cache the matrix factorization, so solving
  117. # with multiple right hand sides is much faster
  118. def solve(r, sym_pos=sym_pos):
  119. if sym_pos:
  120. return sp.linalg.solve(M, r, assume_a="pos")
  121. else:
  122. return sp.linalg.solve(M, r)
  123. # There are many things that can go wrong here, and it's hard to say
  124. # what all of them are. It doesn't really matter: if the matrix can't be
  125. # factorized, return None. get_solver will be called again with different
  126. # inputs, and a new routine will try to factorize the matrix.
  127. except KeyboardInterrupt:
  128. raise
  129. except Exception:
  130. return None
  131. return solve
  132. def _get_delta(A, b, c, x, y, z, tau, kappa, gamma, eta, sparse=False,
  133. lstsq=False, sym_pos=True, cholesky=True, pc=True, ip=False,
  134. permc_spec='MMD_AT_PLUS_A'):
  135. """
  136. Given standard form problem defined by ``A``, ``b``, and ``c``;
  137. current variable estimates ``x``, ``y``, ``z``, ``tau``, and ``kappa``;
  138. algorithmic parameters ``gamma and ``eta;
  139. and options ``sparse``, ``lstsq``, ``sym_pos``, ``cholesky``, ``pc``
  140. (predictor-corrector), and ``ip`` (initial point improvement),
  141. get the search direction for increments to the variable estimates.
  142. Parameters
  143. ----------
  144. As defined in [4], except:
  145. sparse : bool
  146. True if the system to be solved is sparse. This is typically set
  147. True when the original ``A_ub`` and ``A_eq`` arrays are sparse.
  148. lstsq : bool
  149. True if the system is ill-conditioned and/or (nearly) singular and
  150. thus a more robust least-squares solver is desired. This is sometimes
  151. needed as the solution is approached.
  152. sym_pos : bool
  153. True if the system matrix is symmetric positive definite
  154. Sometimes this needs to be set false as the solution is approached,
  155. even when the system should be symmetric positive definite, due to
  156. numerical difficulties.
  157. cholesky : bool
  158. True if the system is to be solved by Cholesky, rather than LU,
  159. decomposition. This is typically faster unless the problem is very
  160. small or prone to numerical difficulties.
  161. pc : bool
  162. True if the predictor-corrector method of Mehrota is to be used. This
  163. is almost always (if not always) beneficial. Even though it requires
  164. the solution of an additional linear system, the factorization
  165. is typically (implicitly) reused so solution is efficient, and the
  166. number of algorithm iterations is typically reduced.
  167. ip : bool
  168. True if the improved initial point suggestion due to [4] section 4.3
  169. is desired. It's unclear whether this is beneficial.
  170. permc_spec : str (default = 'MMD_AT_PLUS_A')
  171. (Has effect only with ``sparse = True``, ``lstsq = False``, ``sym_pos =
  172. True``.) A matrix is factorized in each iteration of the algorithm.
  173. This option specifies how to permute the columns of the matrix for
  174. sparsity preservation. Acceptable values are:
  175. - ``NATURAL``: natural ordering.
  176. - ``MMD_ATA``: minimum degree ordering on the structure of A^T A.
  177. - ``MMD_AT_PLUS_A``: minimum degree ordering on the structure of A^T+A.
  178. - ``COLAMD``: approximate minimum degree column ordering.
  179. This option can impact the convergence of the
  180. interior point algorithm; test different values to determine which
  181. performs best for your problem. For more information, refer to
  182. ``scipy.sparse.linalg.splu``.
  183. Returns
  184. -------
  185. Search directions as defined in [4]
  186. References
  187. ----------
  188. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  189. optimizer for linear programming: an implementation of the
  190. homogeneous algorithm." High performance optimization. Springer US,
  191. 2000. 197-232.
  192. """
  193. if A.shape[0] == 0:
  194. # If there are no constraints, some solvers fail (understandably)
  195. # rather than returning empty solution. This gets the job done.
  196. lstsq, sym_pos, cholesky = False, True, False
  197. n_x = len(x)
  198. # [4] Equation 8.8
  199. r_P = b * tau - A.dot(x)
  200. r_D = c * tau - A.T.dot(y) - z
  201. r_G = c.dot(x) - b.transpose().dot(y) + kappa
  202. mu = (x.dot(z) + tau * kappa) / (n_x + 1)
  203. # Assemble M from [4] Equation 8.31
  204. Dinv = x / z
  205. if sparse:
  206. M = A.dot(sps.diags_array(Dinv, format="csc").dot(A.T))
  207. else:
  208. M = A.dot(Dinv.reshape(-1, 1) * A.T)
  209. solve = _get_solver(M, sparse, lstsq, sym_pos, cholesky, permc_spec)
  210. # pc: "predictor-corrector" [4] Section 4.1
  211. # In development this option could be turned off
  212. # but it always seems to improve performance substantially
  213. n_corrections = 1 if pc else 0
  214. i = 0
  215. alpha, d_x, d_z, d_tau, d_kappa = 0, 0, 0, 0, 0
  216. while i <= n_corrections:
  217. # Reference [4] Eq. 8.6
  218. rhatp = eta(gamma) * r_P
  219. rhatd = eta(gamma) * r_D
  220. rhatg = eta(gamma) * r_G
  221. # Reference [4] Eq. 8.7
  222. rhatxs = gamma * mu - x * z
  223. rhattk = gamma * mu - tau * kappa
  224. if i == 1:
  225. if ip: # if the correction is to get "initial point"
  226. # Reference [4] Eq. 8.23
  227. rhatxs = ((1 - alpha) * gamma * mu -
  228. x * z - alpha**2 * d_x * d_z)
  229. rhattk = ((1 - alpha) * gamma * mu -
  230. tau * kappa -
  231. alpha**2 * d_tau * d_kappa)
  232. else: # if the correction is for "predictor-corrector"
  233. # Reference [4] Eq. 8.13
  234. rhatxs -= d_x * d_z
  235. rhattk -= d_tau * d_kappa
  236. # sometimes numerical difficulties arise as the solution is approached
  237. # this loop tries to solve the equations using a sequence of functions
  238. # for solve. For dense systems, the order is:
  239. # 1. scipy.linalg.cho_factor/scipy.linalg.cho_solve,
  240. # 2. scipy.linalg.solve w/ sym_pos = True,
  241. # 3. scipy.linalg.solve w/ sym_pos = False, and if all else fails
  242. # 4. scipy.linalg.lstsq
  243. # For sparse systems, the order is:
  244. # 1. sksparse.cholmod.cholesky (if available)
  245. # 2. scipy.sparse.linalg.factorized (if umfpack available)
  246. # 3. scipy.sparse.linalg.splu
  247. # 4. scipy.sparse.linalg.lsqr
  248. solved = False
  249. while not solved:
  250. try:
  251. # [4] Equation 8.28
  252. p, q = _sym_solve(Dinv, A, c, b, solve)
  253. # [4] Equation 8.29
  254. u, v = _sym_solve(Dinv, A, rhatd -
  255. (1 / x) * rhatxs, rhatp, solve)
  256. if np.any(np.isnan(p)) or np.any(np.isnan(q)):
  257. raise LinAlgError
  258. solved = True
  259. except (LinAlgError, ValueError, TypeError) as e:
  260. # Usually this doesn't happen. If it does, it happens when
  261. # there are redundant constraints or when approaching the
  262. # solution. If so, change solver.
  263. if cholesky:
  264. cholesky = False
  265. warn(
  266. "Solving system with option 'cholesky':True "
  267. "failed. It is normal for this to happen "
  268. "occasionally, especially as the solution is "
  269. "approached. However, if you see this frequently, "
  270. "consider setting option 'cholesky' to False.",
  271. OptimizeWarning, stacklevel=5)
  272. elif sym_pos:
  273. sym_pos = False
  274. warn(
  275. "Solving system with option 'sym_pos':True "
  276. "failed. It is normal for this to happen "
  277. "occasionally, especially as the solution is "
  278. "approached. However, if you see this frequently, "
  279. "consider setting option 'sym_pos' to False.",
  280. OptimizeWarning, stacklevel=5)
  281. elif not lstsq:
  282. lstsq = True
  283. warn(
  284. "Solving system with option 'sym_pos':False "
  285. "failed. This may happen occasionally, "
  286. "especially as the solution is "
  287. "approached. However, if you see this frequently, "
  288. "your problem may be numerically challenging. "
  289. "If you cannot improve the formulation, consider "
  290. "setting 'lstsq' to True. Consider also setting "
  291. "`presolve` to True, if it is not already.",
  292. OptimizeWarning, stacklevel=5)
  293. else:
  294. raise e
  295. solve = _get_solver(M, sparse, lstsq, sym_pos,
  296. cholesky, permc_spec)
  297. # [4] Results after 8.29
  298. d_tau = ((rhatg + 1 / tau * rhattk - (-c.dot(u) + b.dot(v))) /
  299. (1 / tau * kappa + (-c.dot(p) + b.dot(q))))
  300. d_x = u + p * d_tau
  301. d_y = v + q * d_tau
  302. # [4] Relations between after 8.25 and 8.26
  303. d_z = (1 / x) * (rhatxs - z * d_x)
  304. d_kappa = 1 / tau * (rhattk - kappa * d_tau)
  305. # [4] 8.12 and "Let alpha be the maximal possible step..." before 8.23
  306. alpha = _get_step(x, d_x, z, d_z, tau, d_tau, kappa, d_kappa, 1)
  307. if ip: # initial point - see [4] 4.4
  308. gamma = 10
  309. else: # predictor-corrector, [4] definition after 8.12
  310. beta1 = 0.1 # [4] pg. 220 (Table 8.1)
  311. gamma = (1 - alpha)**2 * min(beta1, (1 - alpha))
  312. i += 1
  313. return d_x, d_y, d_z, d_tau, d_kappa
  314. def _sym_solve(Dinv, A, r1, r2, solve):
  315. """
  316. An implementation of [4] equation 8.31 and 8.32
  317. References
  318. ----------
  319. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  320. optimizer for linear programming: an implementation of the
  321. homogeneous algorithm." High performance optimization. Springer US,
  322. 2000. 197-232.
  323. """
  324. # [4] 8.31
  325. r = r2 + A.dot(Dinv * r1)
  326. v = solve(r)
  327. # [4] 8.32
  328. u = Dinv * (A.T.dot(v) - r1)
  329. return u, v
  330. def _get_step(x, d_x, z, d_z, tau, d_tau, kappa, d_kappa, alpha0):
  331. """
  332. An implementation of [4] equation 8.21
  333. References
  334. ----------
  335. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  336. optimizer for linear programming: an implementation of the
  337. homogeneous algorithm." High performance optimization. Springer US,
  338. 2000. 197-232.
  339. """
  340. # [4] 4.3 Equation 8.21, ignoring 8.20 requirement
  341. # same step is taken in primal and dual spaces
  342. # alpha0 is basically beta3 from [4] Table 8.1, but instead of beta3
  343. # the value 1 is used in Mehrota corrector and initial point correction
  344. i_x = d_x < 0
  345. i_z = d_z < 0
  346. alpha_x = alpha0 * np.min(x[i_x] / -d_x[i_x]) if np.any(i_x) else 1
  347. alpha_tau = alpha0 * tau / -d_tau if d_tau < 0 else 1
  348. alpha_z = alpha0 * np.min(z[i_z] / -d_z[i_z]) if np.any(i_z) else 1
  349. alpha_kappa = alpha0 * kappa / -d_kappa if d_kappa < 0 else 1
  350. alpha = np.min([1, alpha_x, alpha_tau, alpha_z, alpha_kappa])
  351. return alpha
  352. def _get_message(status):
  353. """
  354. Given problem status code, return a more detailed message.
  355. Parameters
  356. ----------
  357. status : int
  358. An integer representing the exit status of the optimization::
  359. 0 : Optimization terminated successfully
  360. 1 : Iteration limit reached
  361. 2 : Problem appears to be infeasible
  362. 3 : Problem appears to be unbounded
  363. 4 : Serious numerical difficulties encountered
  364. Returns
  365. -------
  366. message : str
  367. A string descriptor of the exit status of the optimization.
  368. """
  369. messages = (
  370. ["Optimization terminated successfully.",
  371. "The iteration limit was reached before the algorithm converged.",
  372. "The algorithm terminated successfully and determined that the "
  373. "problem is infeasible.",
  374. "The algorithm terminated successfully and determined that the "
  375. "problem is unbounded.",
  376. "Numerical difficulties were encountered before the problem "
  377. "converged. Please check your problem formulation for errors, "
  378. "independence of linear equality constraints, and reasonable "
  379. "scaling and matrix condition numbers. If you continue to "
  380. "encounter this error, please submit a bug report."
  381. ])
  382. return messages[status]
  383. def _do_step(x, y, z, tau, kappa, d_x, d_y, d_z, d_tau, d_kappa, alpha):
  384. """
  385. An implementation of [4] Equation 8.9
  386. References
  387. ----------
  388. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  389. optimizer for linear programming: an implementation of the
  390. homogeneous algorithm." High performance optimization. Springer US,
  391. 2000. 197-232.
  392. """
  393. x = x + alpha * d_x
  394. tau = tau + alpha * d_tau
  395. z = z + alpha * d_z
  396. kappa = kappa + alpha * d_kappa
  397. y = y + alpha * d_y
  398. return x, y, z, tau, kappa
  399. def _get_blind_start(shape):
  400. """
  401. Return the starting point from [4] 4.4
  402. References
  403. ----------
  404. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  405. optimizer for linear programming: an implementation of the
  406. homogeneous algorithm." High performance optimization. Springer US,
  407. 2000. 197-232.
  408. """
  409. m, n = shape
  410. x0 = np.ones(n)
  411. y0 = np.zeros(m)
  412. z0 = np.ones(n)
  413. tau0 = 1
  414. kappa0 = 1
  415. return x0, y0, z0, tau0, kappa0
  416. def _indicators(A, b, c, c0, x, y, z, tau, kappa):
  417. """
  418. Implementation of several equations from [4] used as indicators of
  419. the status of optimization.
  420. References
  421. ----------
  422. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  423. optimizer for linear programming: an implementation of the
  424. homogeneous algorithm." High performance optimization. Springer US,
  425. 2000. 197-232.
  426. """
  427. # residuals for termination are relative to initial values
  428. x0, y0, z0, tau0, kappa0 = _get_blind_start(A.shape)
  429. # See [4], Section 4 - The Homogeneous Algorithm, Equation 8.8
  430. def r_p(x, tau):
  431. return b * tau - A.dot(x)
  432. def r_d(y, z, tau):
  433. return c * tau - A.T.dot(y) - z
  434. def r_g(x, y, kappa):
  435. return kappa + c.dot(x) - b.dot(y)
  436. # np.dot unpacks if they are arrays of size one
  437. def mu(x, tau, z, kappa):
  438. return (x.dot(z) + np.dot(tau, kappa)) / (len(x) + 1)
  439. obj = c.dot(x / tau) + c0
  440. def norm(a):
  441. return np.linalg.norm(a)
  442. # See [4], Section 4.5 - The Stopping Criteria
  443. r_p0 = r_p(x0, tau0)
  444. r_d0 = r_d(y0, z0, tau0)
  445. r_g0 = r_g(x0, y0, kappa0)
  446. mu_0 = mu(x0, tau0, z0, kappa0)
  447. rho_A = norm(c.T.dot(x) - b.T.dot(y)) / (tau + norm(b.T.dot(y)))
  448. rho_p = norm(r_p(x, tau)) / max(1, norm(r_p0))
  449. rho_d = norm(r_d(y, z, tau)) / max(1, norm(r_d0))
  450. rho_g = norm(r_g(x, y, kappa)) / max(1, norm(r_g0))
  451. rho_mu = mu(x, tau, z, kappa) / mu_0
  452. return rho_p, rho_d, rho_A, rho_g, rho_mu, obj
  453. def _display_iter(rho_p, rho_d, rho_g, alpha, rho_mu, obj, header=False):
  454. """
  455. Print indicators of optimization status to the console.
  456. Parameters
  457. ----------
  458. rho_p : float
  459. The (normalized) primal feasibility, see [4] 4.5
  460. rho_d : float
  461. The (normalized) dual feasibility, see [4] 4.5
  462. rho_g : float
  463. The (normalized) duality gap, see [4] 4.5
  464. alpha : float
  465. The step size, see [4] 4.3
  466. rho_mu : float
  467. The (normalized) path parameter, see [4] 4.5
  468. obj : float
  469. The objective function value of the current iterate
  470. header : bool
  471. True if a header is to be printed
  472. References
  473. ----------
  474. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  475. optimizer for linear programming: an implementation of the
  476. homogeneous algorithm." High performance optimization. Springer US,
  477. 2000. 197-232.
  478. """
  479. if header:
  480. print("Primal Feasibility ",
  481. "Dual Feasibility ",
  482. "Duality Gap ",
  483. "Step ",
  484. "Path Parameter ",
  485. "Objective ")
  486. # no clue why this works
  487. fmt = '{0:<20.13}{1:<20.13}{2:<20.13}{3:<17.13}{4:<20.13}{5:<20.13}'
  488. print(fmt.format(
  489. float(rho_p),
  490. float(rho_d),
  491. float(rho_g),
  492. alpha if isinstance(alpha, str) else float(alpha),
  493. float(rho_mu),
  494. float(obj)))
  495. def _ip_hsd(A, b, c, c0, alpha0, beta, maxiter, disp, tol, sparse, lstsq,
  496. sym_pos, cholesky, pc, ip, permc_spec, callback, postsolve_args):
  497. r"""
  498. Solve a linear programming problem in standard form:
  499. Minimize::
  500. c @ x
  501. Subject to::
  502. A @ x == b
  503. x >= 0
  504. using the interior point method of [4].
  505. Parameters
  506. ----------
  507. A : 2-D array
  508. 2-D array such that ``A @ x``, gives the values of the equality
  509. constraints at ``x``.
  510. b : 1-D array
  511. 1-D array of values representing the RHS of each equality constraint
  512. (row) in ``A`` (for standard form problem).
  513. c : 1-D array
  514. Coefficients of the linear objective function to be minimized (for
  515. standard form problem).
  516. c0 : float
  517. Constant term in objective function due to fixed (and eliminated)
  518. variables. (Purely for display.)
  519. alpha0 : float
  520. The maximal step size for Mehrota's predictor-corrector search
  521. direction; see :math:`\beta_3`of [4] Table 8.1
  522. beta : float
  523. The desired reduction of the path parameter :math:`\mu` (see [6]_)
  524. maxiter : int
  525. The maximum number of iterations of the algorithm.
  526. disp : bool
  527. Set to ``True`` if indicators of optimization status are to be printed
  528. to the console each iteration.
  529. tol : float
  530. Termination tolerance; see [4]_ Section 4.5.
  531. sparse : bool
  532. Set to ``True`` if the problem is to be treated as sparse. However,
  533. the inputs ``A_eq`` and ``A_ub`` should nonetheless be provided as
  534. (dense) arrays rather than sparse matrices.
  535. lstsq : bool
  536. Set to ``True`` if the problem is expected to be very poorly
  537. conditioned. This should always be left as ``False`` unless severe
  538. numerical difficulties are frequently encountered, and a better option
  539. would be to improve the formulation of the problem.
  540. sym_pos : bool
  541. Leave ``True`` if the problem is expected to yield a well conditioned
  542. symmetric positive definite normal equation matrix (almost always).
  543. cholesky : bool
  544. Set to ``True`` if the normal equations are to be solved by explicit
  545. Cholesky decomposition followed by explicit forward/backward
  546. substitution. This is typically faster for moderate, dense problems
  547. that are numerically well-behaved.
  548. pc : bool
  549. Leave ``True`` if the predictor-corrector method of Mehrota is to be
  550. used. This is almost always (if not always) beneficial.
  551. ip : bool
  552. Set to ``True`` if the improved initial point suggestion due to [4]_
  553. Section 4.3 is desired. It's unclear whether this is beneficial.
  554. permc_spec : str (default = 'MMD_AT_PLUS_A')
  555. (Has effect only with ``sparse = True``, ``lstsq = False``, ``sym_pos =
  556. True``.) A matrix is factorized in each iteration of the algorithm.
  557. This option specifies how to permute the columns of the matrix for
  558. sparsity preservation. Acceptable values are:
  559. - ``NATURAL``: natural ordering.
  560. - ``MMD_ATA``: minimum degree ordering on the structure of A^T A.
  561. - ``MMD_AT_PLUS_A``: minimum degree ordering on the structure of A^T+A.
  562. - ``COLAMD``: approximate minimum degree column ordering.
  563. This option can impact the convergence of the
  564. interior point algorithm; test different values to determine which
  565. performs best for your problem. For more information, refer to
  566. ``scipy.sparse.linalg.splu``.
  567. callback : callable, optional
  568. If a callback function is provided, it will be called within each
  569. iteration of the algorithm. The callback function must accept a single
  570. `scipy.optimize.OptimizeResult` consisting of the following fields:
  571. x : 1-D array
  572. Current solution vector
  573. fun : float
  574. Current value of the objective function
  575. success : bool
  576. True only when an algorithm has completed successfully,
  577. so this is always False as the callback function is called
  578. only while the algorithm is still iterating.
  579. slack : 1-D array
  580. The values of the slack variables. Each slack variable
  581. corresponds to an inequality constraint. If the slack is zero,
  582. the corresponding constraint is active.
  583. con : 1-D array
  584. The (nominally zero) residuals of the equality constraints,
  585. that is, ``b - A_eq @ x``
  586. phase : int
  587. The phase of the algorithm being executed. This is always
  588. 1 for the interior-point method because it has only one phase.
  589. status : int
  590. For revised simplex, this is always 0 because if a different
  591. status is detected, the algorithm terminates.
  592. nit : int
  593. The number of iterations performed.
  594. message : str
  595. A string descriptor of the exit status of the optimization.
  596. postsolve_args : tuple
  597. Data needed by _postsolve to convert the solution to the standard-form
  598. problem into the solution to the original problem.
  599. Returns
  600. -------
  601. x_hat : float
  602. Solution vector (for standard form problem).
  603. status : int
  604. An integer representing the exit status of the optimization::
  605. 0 : Optimization terminated successfully
  606. 1 : Iteration limit reached
  607. 2 : Problem appears to be infeasible
  608. 3 : Problem appears to be unbounded
  609. 4 : Serious numerical difficulties encountered
  610. message : str
  611. A string descriptor of the exit status of the optimization.
  612. iteration : int
  613. The number of iterations taken to solve the problem
  614. References
  615. ----------
  616. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  617. optimizer for linear programming: an implementation of the
  618. homogeneous algorithm." High performance optimization. Springer US,
  619. 2000. 197-232.
  620. .. [6] Freund, Robert M. "Primal-Dual Interior-Point Methods for Linear
  621. Programming based on Newton's Method." Unpublished Course Notes,
  622. March 2004. Available 2/25/2017 at:
  623. https://ocw.mit.edu/courses/sloan-school-of-management/15-084j-nonlinear-programming-spring-2004/lecture-notes/lec14_int_pt_mthd.pdf
  624. """
  625. iteration = 0
  626. # default initial point
  627. x, y, z, tau, kappa = _get_blind_start(A.shape)
  628. # first iteration is special improvement of initial point
  629. ip = ip if pc else False
  630. # [4] 4.5
  631. rho_p, rho_d, rho_A, rho_g, rho_mu, obj = _indicators(
  632. A, b, c, c0, x, y, z, tau, kappa)
  633. go = rho_p > tol or rho_d > tol or rho_A > tol # we might get lucky : )
  634. if disp:
  635. _display_iter(rho_p, rho_d, rho_g, "-", rho_mu, obj, header=True)
  636. if callback is not None:
  637. x_o, fun, slack, con = _postsolve(x/tau, postsolve_args)
  638. res = OptimizeResult({'x': x_o, 'fun': fun, 'slack': slack,
  639. 'con': con, 'nit': iteration, 'phase': 1,
  640. 'complete': False, 'status': 0,
  641. 'message': "", 'success': False})
  642. callback(res)
  643. status = 0
  644. message = "Optimization terminated successfully."
  645. if sparse:
  646. A = sps.csc_array(A)
  647. while go:
  648. iteration += 1
  649. if ip: # initial point
  650. # [4] Section 4.4
  651. gamma = 1
  652. def eta(g):
  653. return 1
  654. else:
  655. # gamma = 0 in predictor step according to [4] 4.1
  656. # if predictor/corrector is off, use mean of complementarity [6]
  657. # 5.1 / [4] Below Figure 10-4
  658. gamma = 0 if pc else beta * np.mean(z * x)
  659. # [4] Section 4.1
  660. def eta(g=gamma):
  661. return 1 - g
  662. try:
  663. # Solve [4] 8.6 and 8.7/8.13/8.23
  664. d_x, d_y, d_z, d_tau, d_kappa = _get_delta(
  665. A, b, c, x, y, z, tau, kappa, gamma, eta,
  666. sparse, lstsq, sym_pos, cholesky, pc, ip, permc_spec)
  667. if ip: # initial point
  668. # [4] 4.4
  669. # Formula after 8.23 takes a full step regardless if this will
  670. # take it negative
  671. alpha = 1.0
  672. x, y, z, tau, kappa = _do_step(
  673. x, y, z, tau, kappa, d_x, d_y,
  674. d_z, d_tau, d_kappa, alpha)
  675. x[x < 1] = 1
  676. z[z < 1] = 1
  677. tau = max(1, tau)
  678. kappa = max(1, kappa)
  679. ip = False # done with initial point
  680. else:
  681. # [4] Section 4.3
  682. alpha = _get_step(x, d_x, z, d_z, tau,
  683. d_tau, kappa, d_kappa, alpha0)
  684. # [4] Equation 8.9
  685. x, y, z, tau, kappa = _do_step(
  686. x, y, z, tau, kappa, d_x, d_y, d_z, d_tau, d_kappa, alpha)
  687. except (LinAlgError, FloatingPointError,
  688. ValueError, ZeroDivisionError):
  689. # this can happen when sparse solver is used and presolve
  690. # is turned off. Also observed ValueError in AppVeyor Python 3.6
  691. # Win32 build (PR #8676). I've never seen it otherwise.
  692. status = 4
  693. message = _get_message(status)
  694. break
  695. # [4] 4.5
  696. rho_p, rho_d, rho_A, rho_g, rho_mu, obj = _indicators(
  697. A, b, c, c0, x, y, z, tau, kappa)
  698. go = rho_p > tol or rho_d > tol or rho_A > tol
  699. if disp:
  700. _display_iter(rho_p, rho_d, rho_g, alpha, rho_mu, obj)
  701. if callback is not None:
  702. x_o, fun, slack, con = _postsolve(x/tau, postsolve_args)
  703. res = OptimizeResult({'x': x_o, 'fun': fun, 'slack': slack,
  704. 'con': con, 'nit': iteration, 'phase': 1,
  705. 'complete': False, 'status': 0,
  706. 'message': "", 'success': False})
  707. callback(res)
  708. # [4] 4.5
  709. inf1 = (rho_p < tol and rho_d < tol and rho_g < tol and tau < tol *
  710. max(1, kappa))
  711. inf2 = rho_mu < tol and tau < tol * min(1, kappa)
  712. if inf1 or inf2:
  713. # [4] Lemma 8.4 / Theorem 8.3
  714. if b.transpose().dot(y) > tol:
  715. status = 2
  716. else: # elif c.T.dot(x) < tol: ? Probably not necessary.
  717. status = 3
  718. message = _get_message(status)
  719. break
  720. elif iteration >= maxiter:
  721. status = 1
  722. message = _get_message(status)
  723. break
  724. x_hat = x / tau
  725. # [4] Statement after Theorem 8.2
  726. return x_hat, status, message, iteration
  727. def _linprog_ip(c, c0, A, b, callback, postsolve_args, maxiter=1000, tol=1e-8,
  728. disp=False, alpha0=.99995, beta=0.1, sparse=False, lstsq=False,
  729. sym_pos=True, cholesky=None, pc=True, ip=False,
  730. permc_spec='MMD_AT_PLUS_A', **unknown_options):
  731. r"""
  732. Minimize a linear objective function subject to linear
  733. equality and non-negativity constraints using the interior point method
  734. of [4]_. Linear programming is intended to solve problems
  735. of the following form:
  736. Minimize::
  737. c @ x
  738. Subject to::
  739. A @ x == b
  740. x >= 0
  741. User-facing documentation is in _linprog_doc.py.
  742. Parameters
  743. ----------
  744. c : 1-D array
  745. Coefficients of the linear objective function to be minimized.
  746. c0 : float
  747. Constant term in objective function due to fixed (and eliminated)
  748. variables. (Purely for display.)
  749. A : 2-D array
  750. 2-D array such that ``A @ x``, gives the values of the equality
  751. constraints at ``x``.
  752. b : 1-D array
  753. 1-D array of values representing the right hand side of each equality
  754. constraint (row) in ``A``.
  755. callback : callable, optional
  756. Callback function to be executed once per iteration.
  757. postsolve_args : tuple
  758. Data needed by _postsolve to convert the solution to the standard-form
  759. problem into the solution to the original problem.
  760. Options
  761. -------
  762. maxiter : int (default = 1000)
  763. The maximum number of iterations of the algorithm.
  764. tol : float (default = 1e-8)
  765. Termination tolerance to be used for all termination criteria;
  766. see [4]_ Section 4.5.
  767. disp : bool (default = False)
  768. Set to ``True`` if indicators of optimization status are to be printed
  769. to the console each iteration.
  770. alpha0 : float (default = 0.99995)
  771. The maximal step size for Mehrota's predictor-corrector search
  772. direction; see :math:`\beta_{3}` of [4]_ Table 8.1.
  773. beta : float (default = 0.1)
  774. The desired reduction of the path parameter :math:`\mu` (see [6]_)
  775. when Mehrota's predictor-corrector is not in use (uncommon).
  776. sparse : bool (default = False)
  777. Set to ``True`` if the problem is to be treated as sparse after
  778. presolve. If either ``A_eq`` or ``A_ub`` is a sparse matrix,
  779. this option will automatically be set ``True``, and the problem
  780. will be treated as sparse even during presolve. If your constraint
  781. matrices contain mostly zeros and the problem is not very small (less
  782. than about 100 constraints or variables), consider setting ``True``
  783. or providing ``A_eq`` and ``A_ub`` as sparse matrices.
  784. lstsq : bool (default = False)
  785. Set to ``True`` if the problem is expected to be very poorly
  786. conditioned. This should always be left ``False`` unless severe
  787. numerical difficulties are encountered. Leave this at the default
  788. unless you receive a warning message suggesting otherwise.
  789. sym_pos : bool (default = True)
  790. Leave ``True`` if the problem is expected to yield a well conditioned
  791. symmetric positive definite normal equation matrix
  792. (almost always). Leave this at the default unless you receive
  793. a warning message suggesting otherwise.
  794. cholesky : bool (default = True)
  795. Set to ``True`` if the normal equations are to be solved by explicit
  796. Cholesky decomposition followed by explicit forward/backward
  797. substitution. This is typically faster for problems
  798. that are numerically well-behaved.
  799. pc : bool (default = True)
  800. Leave ``True`` if the predictor-corrector method of Mehrota is to be
  801. used. This is almost always (if not always) beneficial.
  802. ip : bool (default = False)
  803. Set to ``True`` if the improved initial point suggestion due to [4]_
  804. Section 4.3 is desired. Whether this is beneficial or not
  805. depends on the problem.
  806. permc_spec : str (default = 'MMD_AT_PLUS_A')
  807. (Has effect only with ``sparse = True``, ``lstsq = False``, ``sym_pos =
  808. True``, and no SuiteSparse.)
  809. A matrix is factorized in each iteration of the algorithm.
  810. This option specifies how to permute the columns of the matrix for
  811. sparsity preservation. Acceptable values are:
  812. - ``NATURAL``: natural ordering.
  813. - ``MMD_ATA``: minimum degree ordering on the structure of A^T A.
  814. - ``MMD_AT_PLUS_A``: minimum degree ordering on the structure of A^T+A.
  815. - ``COLAMD``: approximate minimum degree column ordering.
  816. This option can impact the convergence of the
  817. interior point algorithm; test different values to determine which
  818. performs best for your problem. For more information, refer to
  819. ``scipy.sparse.linalg.splu``.
  820. unknown_options : dict
  821. Optional arguments not used by this particular solver. If
  822. `unknown_options` is non-empty a warning is issued listing all
  823. unused options.
  824. Returns
  825. -------
  826. x : 1-D array
  827. Solution vector.
  828. status : int
  829. An integer representing the exit status of the optimization::
  830. 0 : Optimization terminated successfully
  831. 1 : Iteration limit reached
  832. 2 : Problem appears to be infeasible
  833. 3 : Problem appears to be unbounded
  834. 4 : Serious numerical difficulties encountered
  835. message : str
  836. A string descriptor of the exit status of the optimization.
  837. iteration : int
  838. The number of iterations taken to solve the problem.
  839. Notes
  840. -----
  841. This method implements the algorithm outlined in [4]_ with ideas from [8]_
  842. and a structure inspired by the simpler methods of [6]_.
  843. The primal-dual path following method begins with initial 'guesses' of
  844. the primal and dual variables of the standard form problem and iteratively
  845. attempts to solve the (nonlinear) Karush-Kuhn-Tucker conditions for the
  846. problem with a gradually reduced logarithmic barrier term added to the
  847. objective. This particular implementation uses a homogeneous self-dual
  848. formulation, which provides certificates of infeasibility or unboundedness
  849. where applicable.
  850. The default initial point for the primal and dual variables is that
  851. defined in [4]_ Section 4.4 Equation 8.22. Optionally (by setting initial
  852. point option ``ip=True``), an alternate (potentially improved) starting
  853. point can be calculated according to the additional recommendations of
  854. [4]_ Section 4.4.
  855. A search direction is calculated using the predictor-corrector method
  856. (single correction) proposed by Mehrota and detailed in [4]_ Section 4.1.
  857. (A potential improvement would be to implement the method of multiple
  858. corrections described in [4]_ Section 4.2.) In practice, this is
  859. accomplished by solving the normal equations, [4]_ Section 5.1 Equations
  860. 8.31 and 8.32, derived from the Newton equations [4]_ Section 5 Equations
  861. 8.25 (compare to [4]_ Section 4 Equations 8.6-8.8). The advantage of
  862. solving the normal equations rather than 8.25 directly is that the
  863. matrices involved are symmetric positive definite, so Cholesky
  864. decomposition can be used rather than the more expensive LU factorization.
  865. With default options, the solver used to perform the factorization depends
  866. on third-party software availability and the conditioning of the problem.
  867. For dense problems, solvers are tried in the following order:
  868. 1. ``scipy.linalg.cho_factor``
  869. 2. ``scipy.linalg.solve`` with option ``sym_pos=True``
  870. 3. ``scipy.linalg.solve`` with option ``sym_pos=False``
  871. 4. ``scipy.linalg.lstsq``
  872. For sparse problems:
  873. 1. ``sksparse.cholmod.cholesky`` (if scikit-sparse and SuiteSparse are installed)
  874. 2. ``scipy.sparse.linalg.factorized``
  875. (if scikit-umfpack and SuiteSparse are installed)
  876. 3. ``scipy.sparse.linalg.splu`` (which uses SuperLU distributed with SciPy)
  877. 4. ``scipy.sparse.linalg.lsqr``
  878. If the solver fails for any reason, successively more robust (but slower)
  879. solvers are attempted in the order indicated. Attempting, failing, and
  880. re-starting factorization can be time consuming, so if the problem is
  881. numerically challenging, options can be set to bypass solvers that are
  882. failing. Setting ``cholesky=False`` skips to solver 2,
  883. ``sym_pos=False`` skips to solver 3, and ``lstsq=True`` skips
  884. to solver 4 for both sparse and dense problems.
  885. Potential improvements for combating issues associated with dense
  886. columns in otherwise sparse problems are outlined in [4]_ Section 5.3 and
  887. [10]_ Section 4.1-4.2; the latter also discusses the alleviation of
  888. accuracy issues associated with the substitution approach to free
  889. variables.
  890. After calculating the search direction, the maximum possible step size
  891. that does not activate the non-negativity constraints is calculated, and
  892. the smaller of this step size and unity is applied (as in [4]_ Section
  893. 4.1.) [4]_ Section 4.3 suggests improvements for choosing the step size.
  894. The new point is tested according to the termination conditions of [4]_
  895. Section 4.5. The same tolerance, which can be set using the ``tol`` option,
  896. is used for all checks. (A potential improvement would be to expose
  897. the different tolerances to be set independently.) If optimality,
  898. unboundedness, or infeasibility is detected, the solve procedure
  899. terminates; otherwise it repeats.
  900. The expected problem formulation differs between the top level ``linprog``
  901. module and the method specific solvers. The method specific solvers expect a
  902. problem in standard form:
  903. Minimize::
  904. c @ x
  905. Subject to::
  906. A @ x == b
  907. x >= 0
  908. Whereas the top level ``linprog`` module expects a problem of form:
  909. Minimize::
  910. c @ x
  911. Subject to::
  912. A_ub @ x <= b_ub
  913. A_eq @ x == b_eq
  914. lb <= x <= ub
  915. where ``lb = 0`` and ``ub = None`` unless set in ``bounds``.
  916. The original problem contains equality, upper-bound and variable constraints
  917. whereas the method specific solver requires equality constraints and
  918. variable non-negativity.
  919. ``linprog`` module converts the original problem to standard form by
  920. converting the simple bounds to upper bound constraints, introducing
  921. non-negative slack variables for inequality constraints, and expressing
  922. unbounded variables as the difference between two non-negative variables.
  923. References
  924. ----------
  925. .. [4] Andersen, Erling D., and Knud D. Andersen. "The MOSEK interior point
  926. optimizer for linear programming: an implementation of the
  927. homogeneous algorithm." High performance optimization. Springer US,
  928. 2000. 197-232.
  929. .. [6] Freund, Robert M. "Primal-Dual Interior-Point Methods for Linear
  930. Programming based on Newton's Method." Unpublished Course Notes,
  931. March 2004. Available 2/25/2017 at
  932. https://ocw.mit.edu/courses/sloan-school-of-management/15-084j-nonlinear-programming-spring-2004/lecture-notes/lec14_int_pt_mthd.pdf
  933. .. [8] Andersen, Erling D., and Knud D. Andersen. "Presolving in linear
  934. programming." Mathematical Programming 71.2 (1995): 221-245.
  935. .. [9] Bertsimas, Dimitris, and J. Tsitsiklis. "Introduction to linear
  936. programming." Athena Scientific 1 (1997): 997.
  937. .. [10] Andersen, Erling D., et al. Implementation of interior point methods
  938. for large scale linear programming. HEC/Universite de Geneve, 1996.
  939. """
  940. _check_unknown_options(unknown_options)
  941. # These should be warnings, not errors
  942. if (cholesky or cholesky is None) and sparse and not has_cholmod:
  943. if cholesky:
  944. warn("Sparse cholesky is only available with scikit-sparse. "
  945. "Setting `cholesky = False`",
  946. OptimizeWarning, stacklevel=3)
  947. cholesky = False
  948. if sparse and lstsq:
  949. warn("Option combination 'sparse':True and 'lstsq':True "
  950. "is not recommended.",
  951. OptimizeWarning, stacklevel=3)
  952. if lstsq and cholesky:
  953. warn("Invalid option combination 'lstsq':True "
  954. "and 'cholesky':True; option 'cholesky' has no effect when "
  955. "'lstsq' is set True.",
  956. OptimizeWarning, stacklevel=3)
  957. valid_permc_spec = ('NATURAL', 'MMD_ATA', 'MMD_AT_PLUS_A', 'COLAMD')
  958. if permc_spec.upper() not in valid_permc_spec:
  959. warn("Invalid permc_spec option: '" + str(permc_spec) + "'. "
  960. "Acceptable values are 'NATURAL', 'MMD_ATA', 'MMD_AT_PLUS_A', "
  961. "and 'COLAMD'. Reverting to default.",
  962. OptimizeWarning, stacklevel=3)
  963. permc_spec = 'MMD_AT_PLUS_A'
  964. # This can be an error
  965. if not sym_pos and cholesky:
  966. raise ValueError(
  967. "Invalid option combination 'sym_pos':False "
  968. "and 'cholesky':True: Cholesky decomposition is only possible "
  969. "for symmetric positive definite matrices.")
  970. cholesky = cholesky or (cholesky is None and sym_pos and not lstsq)
  971. x, status, message, iteration = _ip_hsd(A, b, c, c0, alpha0, beta,
  972. maxiter, disp, tol, sparse,
  973. lstsq, sym_pos, cholesky,
  974. pc, ip, permc_spec, callback,
  975. postsolve_args)
  976. return x, status, message, iteration