_remove_redundancy.py 18 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522
  1. """
  2. Routines for removing redundant (linearly dependent) equations from linear
  3. programming equality constraints.
  4. """
  5. # Author: Matt Haberland
  6. import numpy as np
  7. from scipy.linalg import svd
  8. from scipy.linalg.interpolative import interp_decomp
  9. import scipy
  10. from scipy.linalg.blas import dtrsm
  11. def _row_count(A):
  12. """
  13. Counts the number of nonzeros in each row of input array A.
  14. Nonzeros are defined as any element with absolute value greater than
  15. tol = 1e-13. This value should probably be an input to the function.
  16. Parameters
  17. ----------
  18. A : 2-D array
  19. An array representing a matrix
  20. Returns
  21. -------
  22. rowcount : 1-D array
  23. Number of nonzeros in each row of A
  24. """
  25. tol = 1e-13
  26. return np.array((abs(A) > tol).sum(axis=1)).flatten()
  27. def _get_densest(A, eligibleRows):
  28. """
  29. Returns the index of the densest row of A. Ignores rows that are not
  30. eligible for consideration.
  31. Parameters
  32. ----------
  33. A : 2-D array
  34. An array representing a matrix
  35. eligibleRows : 1-D logical array
  36. Values indicate whether the corresponding row of A is eligible
  37. to be considered
  38. Returns
  39. -------
  40. i_densest : int
  41. Index of the densest row in A eligible for consideration
  42. """
  43. rowCounts = _row_count(A)
  44. return np.argmax(rowCounts * eligibleRows)
  45. def _remove_zero_rows(A, b):
  46. """
  47. Eliminates trivial equations from system of equations defined by Ax = b
  48. and identifies trivial infeasibilities
  49. Parameters
  50. ----------
  51. A : 2-D array
  52. An array representing the left-hand side of a system of equations
  53. b : 1-D array
  54. An array representing the right-hand side of a system of equations
  55. Returns
  56. -------
  57. A : 2-D array
  58. An array representing the left-hand side of a system of equations
  59. b : 1-D array
  60. An array representing the right-hand side of a system of equations
  61. status: int
  62. An integer indicating the status of the removal operation
  63. 0: No infeasibility identified
  64. 2: Trivially infeasible
  65. message : str
  66. A string descriptor of the exit status of the optimization.
  67. """
  68. status = 0
  69. message = ""
  70. i_zero = _row_count(A) == 0
  71. A = A[np.logical_not(i_zero), :]
  72. if not np.allclose(b[i_zero], 0):
  73. status = 2
  74. message = "There is a zero row in A_eq with a nonzero corresponding " \
  75. "entry in b_eq. The problem is infeasible."
  76. b = b[np.logical_not(i_zero)]
  77. return A, b, status, message
  78. def bg_update_dense(plu, perm_r, v, j):
  79. LU, p = plu
  80. vperm = v[perm_r]
  81. u = dtrsm(1, LU, vperm, lower=1, diag=1)
  82. LU[:j+1, j] = u[:j+1]
  83. l = u[j+1:]
  84. piv = LU[j, j]
  85. LU[j+1:, j] += (l/piv)
  86. return LU, p
  87. def _remove_redundancy_pivot_dense(A, rhs, true_rank=None):
  88. """
  89. Eliminates redundant equations from system of equations defined by Ax = b
  90. and identifies infeasibilities.
  91. Parameters
  92. ----------
  93. A : 2-D sparse matrix
  94. An matrix representing the left-hand side of a system of equations
  95. rhs : 1-D array
  96. An array representing the right-hand side of a system of equations
  97. Returns
  98. -------
  99. A : 2-D sparse matrix
  100. A matrix representing the left-hand side of a system of equations
  101. rhs : 1-D array
  102. An array representing the right-hand side of a system of equations
  103. status: int
  104. An integer indicating the status of the system
  105. 0: No infeasibility identified
  106. 2: Trivially infeasible
  107. message : str
  108. A string descriptor of the exit status of the optimization.
  109. References
  110. ----------
  111. .. [2] Andersen, Erling D. "Finding all linearly dependent rows in
  112. large-scale linear programming." Optimization Methods and Software
  113. 6.3 (1995): 219-227.
  114. """
  115. tolapiv = 1e-8
  116. tolprimal = 1e-8
  117. status = 0
  118. message = ""
  119. inconsistent = ("There is a linear combination of rows of A_eq that "
  120. "results in zero, suggesting a redundant constraint. "
  121. "However the same linear combination of b_eq is "
  122. "nonzero, suggesting that the constraints conflict "
  123. "and the problem is infeasible.")
  124. A, rhs, status, message = _remove_zero_rows(A, rhs)
  125. if status != 0:
  126. return A, rhs, status, message
  127. m, n = A.shape
  128. v = list(range(m)) # Artificial column indices.
  129. b = list(v) # Basis column indices.
  130. # This is better as a list than a set because column order of basis matrix
  131. # needs to be consistent.
  132. d = [] # Indices of dependent rows
  133. perm_r = None
  134. A_orig = A
  135. A = np.zeros((m, m + n), order='F')
  136. np.fill_diagonal(A, 1)
  137. A[:, m:] = A_orig
  138. e = np.zeros(m)
  139. js_candidates = np.arange(m, m+n, dtype=int) # candidate columns for basis
  140. # manual masking was faster than masked array
  141. js_mask = np.ones(js_candidates.shape, dtype=bool)
  142. # Implements basic algorithm from [2]
  143. # Uses some of the suggested improvements (removing zero rows and
  144. # Bartels-Golub update idea).
  145. # Removing column singletons would be easy, but it is not as important
  146. # because the procedure is performed only on the equality constraint
  147. # matrix from the original problem - not on the canonical form matrix,
  148. # which would have many more column singletons due to slack variables
  149. # from the inequality constraints.
  150. # The thoughts on "crashing" the initial basis are only really useful if
  151. # the matrix is sparse.
  152. lu = np.eye(m, order='F'), np.arange(m) # initial LU is trivial
  153. perm_r = lu[1]
  154. for i in v:
  155. e[i] = 1
  156. if i > 0:
  157. e[i-1] = 0
  158. try: # fails for i==0 and any time it gets ill-conditioned
  159. j = b[i-1]
  160. lu = bg_update_dense(lu, perm_r, A[:, j], i-1)
  161. except Exception:
  162. lu = scipy.linalg.lu_factor(A[:, b])
  163. LU, p = lu
  164. perm_r = list(range(m))
  165. for i1, i2 in enumerate(p):
  166. perm_r[i1], perm_r[i2] = perm_r[i2], perm_r[i1]
  167. pi = scipy.linalg.lu_solve(lu, e, trans=1)
  168. js = js_candidates[js_mask]
  169. batch = 50
  170. # This is a tiny bit faster than looping over columns indivually,
  171. # like for j in js: if abs(A[:,j].transpose().dot(pi)) > tolapiv:
  172. for j_index in range(0, len(js), batch):
  173. j_indices = js[j_index: min(j_index+batch, len(js))]
  174. c = abs(A[:, j_indices].transpose().dot(pi))
  175. if (c > tolapiv).any():
  176. j = js[j_index + np.argmax(c)] # very independent column
  177. b[i] = j
  178. js_mask[j-m] = False
  179. break
  180. else:
  181. bibar = pi.T.dot(rhs.reshape(-1, 1))
  182. bnorm = np.linalg.norm(rhs)
  183. if abs(bibar)/(1+bnorm) > tolprimal: # inconsistent
  184. status = 2
  185. message = inconsistent
  186. return A_orig, rhs, status, message
  187. else: # dependent
  188. d.append(i)
  189. if true_rank is not None and len(d) == m - true_rank:
  190. break # found all redundancies
  191. keep = set(range(m))
  192. keep = list(keep - set(d))
  193. return A_orig[keep, :], rhs[keep], status, message
  194. def _remove_redundancy_pivot_sparse(A, rhs):
  195. """
  196. Eliminates redundant equations from system of equations defined by Ax = b
  197. and identifies infeasibilities.
  198. Parameters
  199. ----------
  200. A : 2-D sparse matrix
  201. An matrix representing the left-hand side of a system of equations
  202. rhs : 1-D array
  203. An array representing the right-hand side of a system of equations
  204. Returns
  205. -------
  206. A : 2-D sparse matrix
  207. A matrix representing the left-hand side of a system of equations
  208. rhs : 1-D array
  209. An array representing the right-hand side of a system of equations
  210. status: int
  211. An integer indicating the status of the system
  212. 0: No infeasibility identified
  213. 2: Trivially infeasible
  214. message : str
  215. A string descriptor of the exit status of the optimization.
  216. References
  217. ----------
  218. .. [2] Andersen, Erling D. "Finding all linearly dependent rows in
  219. large-scale linear programming." Optimization Methods and Software
  220. 6.3 (1995): 219-227.
  221. """
  222. tolapiv = 1e-8
  223. tolprimal = 1e-8
  224. status = 0
  225. message = ""
  226. inconsistent = ("There is a linear combination of rows of A_eq that "
  227. "results in zero, suggesting a redundant constraint. "
  228. "However the same linear combination of b_eq is "
  229. "nonzero, suggesting that the constraints conflict "
  230. "and the problem is infeasible.")
  231. A, rhs, status, message = _remove_zero_rows(A, rhs)
  232. if status != 0:
  233. return A, rhs, status, message
  234. m, n = A.shape
  235. v = list(range(m)) # Artificial column indices.
  236. b = list(v) # Basis column indices.
  237. # This is better as a list than a set because column order of basis matrix
  238. # needs to be consistent.
  239. k = set(range(m, m+n)) # Structural column indices.
  240. d = [] # Indices of dependent rows
  241. A_orig = A
  242. A = scipy.sparse.hstack((scipy.sparse.eye(m), A)).tocsc()
  243. e = np.zeros(m)
  244. # Implements basic algorithm from [2]
  245. # Uses only one of the suggested improvements (removing zero rows).
  246. # Removing column singletons would be easy, but it is not as important
  247. # because the procedure is performed only on the equality constraint
  248. # matrix from the original problem - not on the canonical form matrix,
  249. # which would have many more column singletons due to slack variables
  250. # from the inequality constraints.
  251. # The thoughts on "crashing" the initial basis sound useful, but the
  252. # description of the procedure seems to assume a lot of familiarity with
  253. # the subject; it is not very explicit. I already went through enough
  254. # trouble getting the basic algorithm working, so I was not interested in
  255. # trying to decipher this, too. (Overall, the paper is fraught with
  256. # mistakes and ambiguities - which is strange, because the rest of
  257. # Andersen's papers are quite good.)
  258. # I tried and tried and tried to improve performance using the
  259. # Bartels-Golub update. It works, but it's only practical if the LU
  260. # factorization can be specialized as described, and that is not possible
  261. # until the SciPy SuperLU interface permits control over column
  262. # permutation - see issue #7700.
  263. for i in v:
  264. B = A[:, b]
  265. e[i] = 1
  266. if i > 0:
  267. e[i-1] = 0
  268. pi = scipy.sparse.linalg.spsolve(B.transpose(), e).reshape(-1, 1)
  269. js = list(k-set(b)) # not efficient, but this is not the time sink...
  270. # Due to overhead, it tends to be faster (for problems tested) to
  271. # compute the full matrix-vector product rather than individual
  272. # vector-vector products (with the chance of terminating as soon
  273. # as any are nonzero). For very large matrices, it might be worth
  274. # it to compute, say, 100 or 1000 at a time and stop when a nonzero
  275. # is found.
  276. c = (np.abs(A[:, js].transpose().dot(pi)) > tolapiv).nonzero()[0]
  277. if len(c) > 0: # independent
  278. j = js[c[0]]
  279. # in a previous commit, the previous line was changed to choose
  280. # index j corresponding with the maximum dot product.
  281. # While this avoided issues with almost
  282. # singular matrices, it slowed the routine in most NETLIB tests.
  283. # I think this is because these columns were denser than the
  284. # first column with nonzero dot product (c[0]).
  285. # It would be nice to have a heuristic that balances sparsity with
  286. # high dot product, but I don't think it's worth the time to
  287. # develop one right now. Bartels-Golub update is a much higher
  288. # priority.
  289. b[i] = j # replace artificial column
  290. else:
  291. bibar = pi.T.dot(rhs.reshape(-1, 1))
  292. bnorm = np.linalg.norm(rhs)
  293. if abs(bibar)/(1 + bnorm) > tolprimal:
  294. status = 2
  295. message = inconsistent
  296. return A_orig, rhs, status, message
  297. else: # dependent
  298. d.append(i)
  299. keep = set(range(m))
  300. keep = list(keep - set(d))
  301. return A_orig[keep, :], rhs[keep], status, message
  302. def _remove_redundancy_svd(A, b):
  303. """
  304. Eliminates redundant equations from system of equations defined by Ax = b
  305. and identifies infeasibilities.
  306. Parameters
  307. ----------
  308. A : 2-D array
  309. An array representing the left-hand side of a system of equations
  310. b : 1-D array
  311. An array representing the right-hand side of a system of equations
  312. Returns
  313. -------
  314. A : 2-D array
  315. An array representing the left-hand side of a system of equations
  316. b : 1-D array
  317. An array representing the right-hand side of a system of equations
  318. status: int
  319. An integer indicating the status of the system
  320. 0: No infeasibility identified
  321. 2: Trivially infeasible
  322. message : str
  323. A string descriptor of the exit status of the optimization.
  324. References
  325. ----------
  326. .. [2] Andersen, Erling D. "Finding all linearly dependent rows in
  327. large-scale linear programming." Optimization Methods and Software
  328. 6.3 (1995): 219-227.
  329. """
  330. A, b, status, message = _remove_zero_rows(A, b)
  331. if status != 0:
  332. return A, b, status, message
  333. U, s, Vh = svd(A)
  334. eps = np.finfo(float).eps
  335. tol = s.max() * max(A.shape) * eps
  336. m, n = A.shape
  337. s_min = s[-1] if m <= n else 0
  338. # this algorithm is faster than that of [2] when the nullspace is small
  339. # but it could probably be improvement by randomized algorithms and with
  340. # a sparse implementation.
  341. # it relies on repeated singular value decomposition to find linearly
  342. # dependent rows (as identified by columns of U that correspond with zero
  343. # singular values). Unfortunately, only one row can be removed per
  344. # decomposition (I tried otherwise; doing so can cause problems.)
  345. # It would be nice if we could do truncated SVD like sp.sparse.linalg.svds
  346. # but that function is unreliable at finding singular values near zero.
  347. # Finding max eigenvalue L of A A^T, then largest eigenvalue (and
  348. # associated eigenvector) of -A A^T + L I (I is identity) via power
  349. # iteration would also work in theory, but is only efficient if the
  350. # smallest nonzero eigenvalue of A A^T is close to the largest nonzero
  351. # eigenvalue.
  352. while abs(s_min) < tol:
  353. v = U[:, -1] # TODO: return these so user can eliminate from problem?
  354. # rows need to be represented in significant amount
  355. eligibleRows = np.abs(v) > tol * 10e6
  356. if not np.any(eligibleRows) or np.any(np.abs(v.dot(A)) > tol):
  357. status = 4
  358. message = ("Due to numerical issues, redundant equality "
  359. "constraints could not be removed automatically. "
  360. "Try providing your constraint matrices as sparse "
  361. "matrices to activate sparse presolve, try turning "
  362. "off redundancy removal, or try turning off presolve "
  363. "altogether.")
  364. break
  365. if np.any(np.abs(v.dot(b)) > tol * 100): # factor of 100 to fix 10038 and 10349
  366. status = 2
  367. message = ("There is a linear combination of rows of A_eq that "
  368. "results in zero, suggesting a redundant constraint. "
  369. "However the same linear combination of b_eq is "
  370. "nonzero, suggesting that the constraints conflict "
  371. "and the problem is infeasible.")
  372. break
  373. i_remove = _get_densest(A, eligibleRows)
  374. A = np.delete(A, i_remove, axis=0)
  375. b = np.delete(b, i_remove)
  376. U, s, Vh = svd(A)
  377. m, n = A.shape
  378. s_min = s[-1] if m <= n else 0
  379. return A, b, status, message
  380. def _remove_redundancy_id(A, rhs, rank=None, randomized=True):
  381. """Eliminates redundant equations from a system of equations.
  382. Eliminates redundant equations from system of equations defined by Ax = b
  383. and identifies infeasibilities.
  384. Parameters
  385. ----------
  386. A : 2-D array
  387. An array representing the left-hand side of a system of equations
  388. rhs : 1-D array
  389. An array representing the right-hand side of a system of equations
  390. rank : int, optional
  391. The rank of A
  392. randomized: bool, optional
  393. True for randomized interpolative decomposition
  394. Returns
  395. -------
  396. A : 2-D array
  397. An array representing the left-hand side of a system of equations
  398. rhs : 1-D array
  399. An array representing the right-hand side of a system of equations
  400. status: int
  401. An integer indicating the status of the system
  402. 0: No infeasibility identified
  403. 2: Trivially infeasible
  404. message : str
  405. A string descriptor of the exit status of the optimization.
  406. """
  407. status = 0
  408. message = ""
  409. inconsistent = ("There is a linear combination of rows of A_eq that "
  410. "results in zero, suggesting a redundant constraint. "
  411. "However the same linear combination of b_eq is "
  412. "nonzero, suggesting that the constraints conflict "
  413. "and the problem is infeasible.")
  414. A, rhs, status, message = _remove_zero_rows(A, rhs)
  415. if status != 0:
  416. return A, rhs, status, message
  417. m, n = A.shape
  418. k = rank
  419. if rank is None:
  420. k = np.linalg.matrix_rank(A)
  421. idx, proj = interp_decomp(A.T, k, rand=randomized)
  422. # first k entries in idx are indices of the independent rows
  423. # remaining entries are the indices of the m-k dependent rows
  424. # proj provides a linear combinations of rows of A2 that form the
  425. # remaining m-k (dependent) rows. The same linear combination of entries
  426. # in rhs2 must give the remaining m-k entries. If not, the system is
  427. # inconsistent, and the problem is infeasible.
  428. if not np.allclose(rhs[idx[:k]] @ proj, rhs[idx[k:]]):
  429. status = 2
  430. message = inconsistent
  431. # sort indices because the other redundancy removal routines leave rows
  432. # in original order and tests were written with that in mind
  433. idx = sorted(idx[:k])
  434. A2 = A[idx, :]
  435. rhs2 = rhs[idx]
  436. return A2, rhs2, status, message