_linprog_ip.py 41 KB

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