normalforms.py 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406
  1. '''Functions returning normal forms of matrices'''
  2. from collections import defaultdict
  3. from .domainmatrix import DomainMatrix
  4. from .exceptions import DMDomainError, DMShapeError
  5. from sympy.ntheory.modular import symmetric_residue
  6. from sympy.polys.domains import QQ, ZZ
  7. # TODO (future work):
  8. # There are faster algorithms for Smith and Hermite normal forms, which
  9. # we should implement. See e.g. the Kannan-Bachem algorithm:
  10. # <https://www.researchgate.net/publication/220617516_Polynomial_Algorithms_for_Computing_the_Smith_and_Hermite_Normal_Forms_of_an_Integer_Matrix>
  11. def smith_normal_form(m):
  12. '''
  13. Return the Smith Normal Form of a matrix `m` over the ring `domain`.
  14. This will only work if the ring is a principal ideal domain.
  15. Examples
  16. ========
  17. >>> from sympy import ZZ
  18. >>> from sympy.polys.matrices import DomainMatrix
  19. >>> from sympy.polys.matrices.normalforms import smith_normal_form
  20. >>> m = DomainMatrix([[ZZ(12), ZZ(6), ZZ(4)],
  21. ... [ZZ(3), ZZ(9), ZZ(6)],
  22. ... [ZZ(2), ZZ(16), ZZ(14)]], (3, 3), ZZ)
  23. >>> print(smith_normal_form(m).to_Matrix())
  24. Matrix([[1, 0, 0], [0, 10, 0], [0, 0, -30]])
  25. '''
  26. invs = invariant_factors(m)
  27. smf = DomainMatrix.diag(invs, m.domain, m.shape)
  28. return smf
  29. def add_columns(m, i, j, a, b, c, d):
  30. # replace m[:, i] by a*m[:, i] + b*m[:, j]
  31. # and m[:, j] by c*m[:, i] + d*m[:, j]
  32. for k in range(len(m)):
  33. e = m[k][i]
  34. m[k][i] = a*e + b*m[k][j]
  35. m[k][j] = c*e + d*m[k][j]
  36. def invariant_factors(m):
  37. '''
  38. Return the tuple of abelian invariants for a matrix `m`
  39. (as in the Smith-Normal form)
  40. References
  41. ==========
  42. [1] https://en.wikipedia.org/wiki/Smith_normal_form#Algorithm
  43. [2] https://web.archive.org/web/20200331143852/https://sierra.nmsu.edu/morandi/notes/SmithNormalForm.pdf
  44. '''
  45. domain = m.domain
  46. if not domain.is_PID:
  47. msg = "The matrix entries must be over a principal ideal domain"
  48. raise ValueError(msg)
  49. if 0 in m.shape:
  50. return ()
  51. rows, cols = shape = m.shape
  52. m = list(m.to_dense().rep)
  53. def add_rows(m, i, j, a, b, c, d):
  54. # replace m[i, :] by a*m[i, :] + b*m[j, :]
  55. # and m[j, :] by c*m[i, :] + d*m[j, :]
  56. for k in range(cols):
  57. e = m[i][k]
  58. m[i][k] = a*e + b*m[j][k]
  59. m[j][k] = c*e + d*m[j][k]
  60. def clear_column(m):
  61. # make m[1:, 0] zero by row and column operations
  62. if m[0][0] == 0:
  63. return m # pragma: nocover
  64. pivot = m[0][0]
  65. for j in range(1, rows):
  66. if m[j][0] == 0:
  67. continue
  68. d, r = domain.div(m[j][0], pivot)
  69. if r == 0:
  70. add_rows(m, 0, j, 1, 0, -d, 1)
  71. else:
  72. a, b, g = domain.gcdex(pivot, m[j][0])
  73. d_0 = domain.div(m[j][0], g)[0]
  74. d_j = domain.div(pivot, g)[0]
  75. add_rows(m, 0, j, a, b, d_0, -d_j)
  76. pivot = g
  77. return m
  78. def clear_row(m):
  79. # make m[0, 1:] zero by row and column operations
  80. if m[0][0] == 0:
  81. return m # pragma: nocover
  82. pivot = m[0][0]
  83. for j in range(1, cols):
  84. if m[0][j] == 0:
  85. continue
  86. d, r = domain.div(m[0][j], pivot)
  87. if r == 0:
  88. add_columns(m, 0, j, 1, 0, -d, 1)
  89. else:
  90. a, b, g = domain.gcdex(pivot, m[0][j])
  91. d_0 = domain.div(m[0][j], g)[0]
  92. d_j = domain.div(pivot, g)[0]
  93. add_columns(m, 0, j, a, b, d_0, -d_j)
  94. pivot = g
  95. return m
  96. # permute the rows and columns until m[0,0] is non-zero if possible
  97. ind = [i for i in range(rows) if m[i][0] != 0]
  98. if ind and ind[0] != 0:
  99. m[0], m[ind[0]] = m[ind[0]], m[0]
  100. else:
  101. ind = [j for j in range(cols) if m[0][j] != 0]
  102. if ind and ind[0] != 0:
  103. for row in m:
  104. row[0], row[ind[0]] = row[ind[0]], row[0]
  105. # make the first row and column except m[0,0] zero
  106. while (any(m[0][i] != 0 for i in range(1,cols)) or
  107. any(m[i][0] != 0 for i in range(1,rows))):
  108. m = clear_column(m)
  109. m = clear_row(m)
  110. if 1 in shape:
  111. invs = ()
  112. else:
  113. lower_right = DomainMatrix([r[1:] for r in m[1:]], (rows-1, cols-1), domain)
  114. invs = invariant_factors(lower_right)
  115. if m[0][0]:
  116. result = [m[0][0]]
  117. result.extend(invs)
  118. # in case m[0] doesn't divide the invariants of the rest of the matrix
  119. for i in range(len(result)-1):
  120. if result[i] and domain.div(result[i+1], result[i])[1] != 0:
  121. g = domain.gcd(result[i+1], result[i])
  122. result[i+1] = domain.div(result[i], g)[0]*result[i+1]
  123. result[i] = g
  124. else:
  125. break
  126. else:
  127. result = invs + (m[0][0],)
  128. return tuple(result)
  129. def _gcdex(a, b):
  130. r"""
  131. This supports the functions that compute Hermite Normal Form.
  132. Explanation
  133. ===========
  134. Let x, y be the coefficients returned by the extended Euclidean
  135. Algorithm, so that x*a + y*b = g. In the algorithms for computing HNF,
  136. it is critical that x, y not only satisfy the condition of being small
  137. in magnitude -- namely that |x| <= |b|/g, |y| <- |a|/g -- but also that
  138. y == 0 when a | b.
  139. """
  140. x, y, g = ZZ.gcdex(a, b)
  141. if a != 0 and b % a == 0:
  142. y = 0
  143. x = -1 if a < 0 else 1
  144. return x, y, g
  145. def _hermite_normal_form(A):
  146. r"""
  147. Compute the Hermite Normal Form of DomainMatrix *A* over :ref:`ZZ`.
  148. Parameters
  149. ==========
  150. A : :py:class:`~.DomainMatrix` over domain :ref:`ZZ`.
  151. Returns
  152. =======
  153. :py:class:`~.DomainMatrix`
  154. The HNF of matrix *A*.
  155. Raises
  156. ======
  157. DMDomainError
  158. If the domain of the matrix is not :ref:`ZZ`.
  159. References
  160. ==========
  161. .. [1] Cohen, H. *A Course in Computational Algebraic Number Theory.*
  162. (See Algorithm 2.4.5.)
  163. """
  164. if not A.domain.is_ZZ:
  165. raise DMDomainError('Matrix must be over domain ZZ.')
  166. # We work one row at a time, starting from the bottom row, and working our
  167. # way up.
  168. m, n = A.shape
  169. A = A.to_dense().rep.copy()
  170. # Our goal is to put pivot entries in the rightmost columns.
  171. # Invariant: Before processing each row, k should be the index of the
  172. # leftmost column in which we have so far put a pivot.
  173. k = n
  174. for i in range(m - 1, -1, -1):
  175. if k == 0:
  176. # This case can arise when n < m and we've already found n pivots.
  177. # We don't need to consider any more rows, because this is already
  178. # the maximum possible number of pivots.
  179. break
  180. k -= 1
  181. # k now points to the column in which we want to put a pivot.
  182. # We want zeros in all entries to the left of the pivot column.
  183. for j in range(k - 1, -1, -1):
  184. if A[i][j] != 0:
  185. # Replace cols j, k by lin combs of these cols such that, in row i,
  186. # col j has 0, while col k has the gcd of their row i entries. Note
  187. # that this ensures a nonzero entry in col k.
  188. u, v, d = _gcdex(A[i][k], A[i][j])
  189. r, s = A[i][k] // d, A[i][j] // d
  190. add_columns(A, k, j, u, v, -s, r)
  191. b = A[i][k]
  192. # Do not want the pivot entry to be negative.
  193. if b < 0:
  194. add_columns(A, k, k, -1, 0, -1, 0)
  195. b = -b
  196. # The pivot entry will be 0 iff the row was 0 from the pivot col all the
  197. # way to the left. In this case, we are still working on the same pivot
  198. # col for the next row. Therefore:
  199. if b == 0:
  200. k += 1
  201. # If the pivot entry is nonzero, then we want to reduce all entries to its
  202. # right in the sense of the division algorithm, i.e. make them all remainders
  203. # w.r.t. the pivot as divisor.
  204. else:
  205. for j in range(k + 1, n):
  206. q = A[i][j] // b
  207. add_columns(A, j, k, 1, -q, 0, 1)
  208. # Finally, the HNF consists of those columns of A in which we succeeded in making
  209. # a nonzero pivot.
  210. return DomainMatrix.from_rep(A)[:, k:]
  211. def _hermite_normal_form_modulo_D(A, D):
  212. r"""
  213. Perform the mod *D* Hermite Normal Form reduction algorithm on
  214. :py:class:`~.DomainMatrix` *A*.
  215. Explanation
  216. ===========
  217. If *A* is an $m \times n$ matrix of rank $m$, having Hermite Normal Form
  218. $W$, and if *D* is any positive integer known in advance to be a multiple
  219. of $\det(W)$, then the HNF of *A* can be computed by an algorithm that
  220. works mod *D* in order to prevent coefficient explosion.
  221. Parameters
  222. ==========
  223. A : :py:class:`~.DomainMatrix` over :ref:`ZZ`
  224. $m \times n$ matrix, having rank $m$.
  225. D : :ref:`ZZ`
  226. Positive integer, known to be a multiple of the determinant of the
  227. HNF of *A*.
  228. Returns
  229. =======
  230. :py:class:`~.DomainMatrix`
  231. The HNF of matrix *A*.
  232. Raises
  233. ======
  234. DMDomainError
  235. If the domain of the matrix is not :ref:`ZZ`, or
  236. if *D* is given but is not in :ref:`ZZ`.
  237. DMShapeError
  238. If the matrix has more rows than columns.
  239. References
  240. ==========
  241. .. [1] Cohen, H. *A Course in Computational Algebraic Number Theory.*
  242. (See Algorithm 2.4.8.)
  243. """
  244. if not A.domain.is_ZZ:
  245. raise DMDomainError('Matrix must be over domain ZZ.')
  246. if not ZZ.of_type(D) or D < 1:
  247. raise DMDomainError('Modulus D must be positive element of domain ZZ.')
  248. def add_columns_mod_R(m, R, i, j, a, b, c, d):
  249. # replace m[:, i] by (a*m[:, i] + b*m[:, j]) % R
  250. # and m[:, j] by (c*m[:, i] + d*m[:, j]) % R
  251. for k in range(len(m)):
  252. e = m[k][i]
  253. m[k][i] = symmetric_residue((a * e + b * m[k][j]) % R, R)
  254. m[k][j] = symmetric_residue((c * e + d * m[k][j]) % R, R)
  255. W = defaultdict(dict)
  256. m, n = A.shape
  257. if n < m:
  258. raise DMShapeError('Matrix must have at least as many columns as rows.')
  259. A = A.to_dense().rep.copy()
  260. k = n
  261. R = D
  262. for i in range(m - 1, -1, -1):
  263. k -= 1
  264. for j in range(k - 1, -1, -1):
  265. if A[i][j] != 0:
  266. u, v, d = _gcdex(A[i][k], A[i][j])
  267. r, s = A[i][k] // d, A[i][j] // d
  268. add_columns_mod_R(A, R, k, j, u, v, -s, r)
  269. b = A[i][k]
  270. if b == 0:
  271. A[i][k] = b = R
  272. u, v, d = _gcdex(b, R)
  273. for ii in range(m):
  274. W[ii][i] = u*A[ii][k] % R
  275. if W[i][i] == 0:
  276. W[i][i] = R
  277. for j in range(i + 1, m):
  278. q = W[i][j] // W[i][i]
  279. add_columns(W, j, i, 1, -q, 0, 1)
  280. R //= d
  281. return DomainMatrix(W, (m, m), ZZ).to_dense()
  282. def hermite_normal_form(A, *, D=None, check_rank=False):
  283. r"""
  284. Compute the Hermite Normal Form of :py:class:`~.DomainMatrix` *A* over
  285. :ref:`ZZ`.
  286. Examples
  287. ========
  288. >>> from sympy import ZZ
  289. >>> from sympy.polys.matrices import DomainMatrix
  290. >>> from sympy.polys.matrices.normalforms import hermite_normal_form
  291. >>> m = DomainMatrix([[ZZ(12), ZZ(6), ZZ(4)],
  292. ... [ZZ(3), ZZ(9), ZZ(6)],
  293. ... [ZZ(2), ZZ(16), ZZ(14)]], (3, 3), ZZ)
  294. >>> print(hermite_normal_form(m).to_Matrix())
  295. Matrix([[10, 0, 2], [0, 15, 3], [0, 0, 2]])
  296. Parameters
  297. ==========
  298. A : $m \times n$ ``DomainMatrix`` over :ref:`ZZ`.
  299. D : :ref:`ZZ`, optional
  300. Let $W$ be the HNF of *A*. If known in advance, a positive integer *D*
  301. being any multiple of $\det(W)$ may be provided. In this case, if *A*
  302. also has rank $m$, then we may use an alternative algorithm that works
  303. mod *D* in order to prevent coefficient explosion.
  304. check_rank : boolean, optional (default=False)
  305. The basic assumption is that, if you pass a value for *D*, then
  306. you already believe that *A* has rank $m$, so we do not waste time
  307. checking it for you. If you do want this to be checked (and the
  308. ordinary, non-modulo *D* algorithm to be used if the check fails), then
  309. set *check_rank* to ``True``.
  310. Returns
  311. =======
  312. :py:class:`~.DomainMatrix`
  313. The HNF of matrix *A*.
  314. Raises
  315. ======
  316. DMDomainError
  317. If the domain of the matrix is not :ref:`ZZ`, or
  318. if *D* is given but is not in :ref:`ZZ`.
  319. DMShapeError
  320. If the mod *D* algorithm is used but the matrix has more rows than
  321. columns.
  322. References
  323. ==========
  324. .. [1] Cohen, H. *A Course in Computational Algebraic Number Theory.*
  325. (See Algorithms 2.4.5 and 2.4.8.)
  326. """
  327. if not A.domain.is_ZZ:
  328. raise DMDomainError('Matrix must be over domain ZZ.')
  329. if D is not None and (not check_rank or A.convert_to(QQ).rank() == A.shape[0]):
  330. return _hermite_normal_form_modulo_D(A, D)
  331. else:
  332. return _hermite_normal_form(A)