_linprog_ip.py 45 KB

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