inverse.py 11 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475
  1. from sympy.core.numbers import mod_inverse
  2. from .common import MatrixError, NonSquareMatrixError, NonInvertibleMatrixError
  3. from .utilities import _iszero
  4. def _pinv_full_rank(M):
  5. """Subroutine for full row or column rank matrices.
  6. For full row rank matrices, inverse of ``A * A.H`` Exists.
  7. For full column rank matrices, inverse of ``A.H * A`` Exists.
  8. This routine can apply for both cases by checking the shape
  9. and have small decision.
  10. """
  11. if M.is_zero_matrix:
  12. return M.H
  13. if M.rows >= M.cols:
  14. return M.H.multiply(M).inv().multiply(M.H)
  15. else:
  16. return M.H.multiply(M.multiply(M.H).inv())
  17. def _pinv_rank_decomposition(M):
  18. """Subroutine for rank decomposition
  19. With rank decompositions, `A` can be decomposed into two full-
  20. rank matrices, and each matrix can take pseudoinverse
  21. individually.
  22. """
  23. if M.is_zero_matrix:
  24. return M.H
  25. B, C = M.rank_decomposition()
  26. Bp = _pinv_full_rank(B)
  27. Cp = _pinv_full_rank(C)
  28. return Cp.multiply(Bp)
  29. def _pinv_diagonalization(M):
  30. """Subroutine using diagonalization
  31. This routine can sometimes fail if SymPy's eigenvalue
  32. computation is not reliable.
  33. """
  34. if M.is_zero_matrix:
  35. return M.H
  36. A = M
  37. AH = M.H
  38. try:
  39. if M.rows >= M.cols:
  40. P, D = AH.multiply(A).diagonalize(normalize=True)
  41. D_pinv = D.applyfunc(lambda x: 0 if _iszero(x) else 1 / x)
  42. return P.multiply(D_pinv).multiply(P.H).multiply(AH)
  43. else:
  44. P, D = A.multiply(AH).diagonalize(
  45. normalize=True)
  46. D_pinv = D.applyfunc(lambda x: 0 if _iszero(x) else 1 / x)
  47. return AH.multiply(P).multiply(D_pinv).multiply(P.H)
  48. except MatrixError:
  49. raise NotImplementedError(
  50. 'pinv for rank-deficient matrices where '
  51. 'diagonalization of A.H*A fails is not supported yet.')
  52. def _pinv(M, method='RD'):
  53. """Calculate the Moore-Penrose pseudoinverse of the matrix.
  54. The Moore-Penrose pseudoinverse exists and is unique for any matrix.
  55. If the matrix is invertible, the pseudoinverse is the same as the
  56. inverse.
  57. Parameters
  58. ==========
  59. method : String, optional
  60. Specifies the method for computing the pseudoinverse.
  61. If ``'RD'``, Rank-Decomposition will be used.
  62. If ``'ED'``, Diagonalization will be used.
  63. Examples
  64. ========
  65. Computing pseudoinverse by rank decomposition :
  66. >>> from sympy import Matrix
  67. >>> A = Matrix([[1, 2, 3], [4, 5, 6]])
  68. >>> A.pinv()
  69. Matrix([
  70. [-17/18, 4/9],
  71. [ -1/9, 1/9],
  72. [ 13/18, -2/9]])
  73. Computing pseudoinverse by diagonalization :
  74. >>> B = A.pinv(method='ED')
  75. >>> B.simplify()
  76. >>> B
  77. Matrix([
  78. [-17/18, 4/9],
  79. [ -1/9, 1/9],
  80. [ 13/18, -2/9]])
  81. See Also
  82. ========
  83. inv
  84. pinv_solve
  85. References
  86. ==========
  87. .. [1] https://en.wikipedia.org/wiki/Moore-Penrose_pseudoinverse
  88. """
  89. # Trivial case: pseudoinverse of all-zero matrix is its transpose.
  90. if M.is_zero_matrix:
  91. return M.H
  92. if method == 'RD':
  93. return _pinv_rank_decomposition(M)
  94. elif method == 'ED':
  95. return _pinv_diagonalization(M)
  96. else:
  97. raise ValueError('invalid pinv method %s' % repr(method))
  98. def _inv_mod(M, m):
  99. r"""
  100. Returns the inverse of the matrix `K` (mod `m`), if it exists.
  101. Method to find the matrix inverse of `K` (mod `m`) implemented in this function:
  102. * Compute `\mathrm{adj}(K) = \mathrm{cof}(K)^t`, the adjoint matrix of `K`.
  103. * Compute `r = 1/\mathrm{det}(K) \pmod m`.
  104. * `K^{-1} = r\cdot \mathrm{adj}(K) \pmod m`.
  105. Examples
  106. ========
  107. >>> from sympy import Matrix
  108. >>> A = Matrix(2, 2, [1, 2, 3, 4])
  109. >>> A.inv_mod(5)
  110. Matrix([
  111. [3, 1],
  112. [4, 2]])
  113. >>> A.inv_mod(3)
  114. Matrix([
  115. [1, 1],
  116. [0, 1]])
  117. """
  118. if not M.is_square:
  119. raise NonSquareMatrixError()
  120. N = M.cols
  121. det_K = M.det()
  122. det_inv = None
  123. try:
  124. det_inv = mod_inverse(det_K, m)
  125. except ValueError:
  126. raise NonInvertibleMatrixError('Matrix is not invertible (mod %d)' % m)
  127. K_adj = M.adjugate()
  128. K_inv = M.__class__(N, N,
  129. [det_inv * K_adj[i, j] % m for i in range(N) for j in range(N)])
  130. return K_inv
  131. def _verify_invertible(M, iszerofunc=_iszero):
  132. """Initial check to see if a matrix is invertible. Raises or returns
  133. determinant for use in _inv_ADJ."""
  134. if not M.is_square:
  135. raise NonSquareMatrixError("A Matrix must be square to invert.")
  136. d = M.det(method='berkowitz')
  137. zero = d.equals(0)
  138. if zero is None: # if equals() can't decide, will rref be able to?
  139. ok = M.rref(simplify=True)[0]
  140. zero = any(iszerofunc(ok[j, j]) for j in range(ok.rows))
  141. if zero:
  142. raise NonInvertibleMatrixError("Matrix det == 0; not invertible.")
  143. return d
  144. def _inv_ADJ(M, iszerofunc=_iszero):
  145. """Calculates the inverse using the adjugate matrix and a determinant.
  146. See Also
  147. ========
  148. inv
  149. inverse_GE
  150. inverse_LU
  151. inverse_CH
  152. inverse_LDL
  153. """
  154. d = _verify_invertible(M, iszerofunc=iszerofunc)
  155. return M.adjugate() / d
  156. def _inv_GE(M, iszerofunc=_iszero):
  157. """Calculates the inverse using Gaussian elimination.
  158. See Also
  159. ========
  160. inv
  161. inverse_ADJ
  162. inverse_LU
  163. inverse_CH
  164. inverse_LDL
  165. """
  166. from .dense import Matrix
  167. if not M.is_square:
  168. raise NonSquareMatrixError("A Matrix must be square to invert.")
  169. big = Matrix.hstack(M.as_mutable(), Matrix.eye(M.rows))
  170. red = big.rref(iszerofunc=iszerofunc, simplify=True)[0]
  171. if any(iszerofunc(red[j, j]) for j in range(red.rows)):
  172. raise NonInvertibleMatrixError("Matrix det == 0; not invertible.")
  173. return M._new(red[:, big.rows:])
  174. def _inv_LU(M, iszerofunc=_iszero):
  175. """Calculates the inverse using LU decomposition.
  176. See Also
  177. ========
  178. inv
  179. inverse_ADJ
  180. inverse_GE
  181. inverse_CH
  182. inverse_LDL
  183. """
  184. if not M.is_square:
  185. raise NonSquareMatrixError("A Matrix must be square to invert.")
  186. if M.free_symbols:
  187. _verify_invertible(M, iszerofunc=iszerofunc)
  188. return M.LUsolve(M.eye(M.rows), iszerofunc=_iszero)
  189. def _inv_CH(M, iszerofunc=_iszero):
  190. """Calculates the inverse using cholesky decomposition.
  191. See Also
  192. ========
  193. inv
  194. inverse_ADJ
  195. inverse_GE
  196. inverse_LU
  197. inverse_LDL
  198. """
  199. _verify_invertible(M, iszerofunc=iszerofunc)
  200. return M.cholesky_solve(M.eye(M.rows))
  201. def _inv_LDL(M, iszerofunc=_iszero):
  202. """Calculates the inverse using LDL decomposition.
  203. See Also
  204. ========
  205. inv
  206. inverse_ADJ
  207. inverse_GE
  208. inverse_LU
  209. inverse_CH
  210. """
  211. _verify_invertible(M, iszerofunc=iszerofunc)
  212. return M.LDLsolve(M.eye(M.rows))
  213. def _inv_QR(M, iszerofunc=_iszero):
  214. """Calculates the inverse using QR decomposition.
  215. See Also
  216. ========
  217. inv
  218. inverse_ADJ
  219. inverse_GE
  220. inverse_CH
  221. inverse_LDL
  222. """
  223. _verify_invertible(M, iszerofunc=iszerofunc)
  224. return M.QRsolve(M.eye(M.rows))
  225. def _inv_block(M, iszerofunc=_iszero):
  226. """Calculates the inverse using BLOCKWISE inversion.
  227. See Also
  228. ========
  229. inv
  230. inverse_ADJ
  231. inverse_GE
  232. inverse_CH
  233. inverse_LDL
  234. """
  235. from sympy.matrices.expressions.blockmatrix import BlockMatrix
  236. i = M.shape[0]
  237. if i <= 20 :
  238. return M.inv(method="LU", iszerofunc=_iszero)
  239. A = M[:i // 2, :i //2]
  240. B = M[:i // 2, i // 2:]
  241. C = M[i // 2:, :i // 2]
  242. D = M[i // 2:, i // 2:]
  243. try:
  244. D_inv = _inv_block(D)
  245. except NonInvertibleMatrixError:
  246. return M.inv(method="LU", iszerofunc=_iszero)
  247. B_D_i = B*D_inv
  248. BDC = B_D_i*C
  249. A_n = A - BDC
  250. try:
  251. A_n = _inv_block(A_n)
  252. except NonInvertibleMatrixError:
  253. return M.inv(method="LU", iszerofunc=_iszero)
  254. B_n = -A_n*B_D_i
  255. dc = D_inv*C
  256. C_n = -dc*A_n
  257. D_n = D_inv + dc*-B_n
  258. nn = BlockMatrix([[A_n, B_n], [C_n, D_n]]).as_explicit()
  259. return nn
  260. def _inv(M, method=None, iszerofunc=_iszero, try_block_diag=False):
  261. """
  262. Return the inverse of a matrix using the method indicated. Default for
  263. dense matrices is is Gauss elimination, default for sparse matrices is LDL.
  264. Parameters
  265. ==========
  266. method : ('GE', 'LU', 'ADJ', 'CH', 'LDL')
  267. iszerofunc : function, optional
  268. Zero-testing function to use.
  269. try_block_diag : bool, optional
  270. If True then will try to form block diagonal matrices using the
  271. method get_diag_blocks(), invert these individually, and then
  272. reconstruct the full inverse matrix.
  273. Examples
  274. ========
  275. >>> from sympy import SparseMatrix, Matrix
  276. >>> A = SparseMatrix([
  277. ... [ 2, -1, 0],
  278. ... [-1, 2, -1],
  279. ... [ 0, 0, 2]])
  280. >>> A.inv('CH')
  281. Matrix([
  282. [2/3, 1/3, 1/6],
  283. [1/3, 2/3, 1/3],
  284. [ 0, 0, 1/2]])
  285. >>> A.inv(method='LDL') # use of 'method=' is optional
  286. Matrix([
  287. [2/3, 1/3, 1/6],
  288. [1/3, 2/3, 1/3],
  289. [ 0, 0, 1/2]])
  290. >>> A * _
  291. Matrix([
  292. [1, 0, 0],
  293. [0, 1, 0],
  294. [0, 0, 1]])
  295. >>> A = Matrix(A)
  296. >>> A.inv('CH')
  297. Matrix([
  298. [2/3, 1/3, 1/6],
  299. [1/3, 2/3, 1/3],
  300. [ 0, 0, 1/2]])
  301. >>> A.inv('ADJ') == A.inv('GE') == A.inv('LU') == A.inv('CH') == A.inv('LDL') == A.inv('QR')
  302. True
  303. Notes
  304. =====
  305. According to the ``method`` keyword, it calls the appropriate method:
  306. GE .... inverse_GE(); default for dense matrices
  307. LU .... inverse_LU()
  308. ADJ ... inverse_ADJ()
  309. CH ... inverse_CH()
  310. LDL ... inverse_LDL(); default for sparse matrices
  311. QR ... inverse_QR()
  312. Note, the GE and LU methods may require the matrix to be simplified
  313. before it is inverted in order to properly detect zeros during
  314. pivoting. In difficult cases a custom zero detection function can
  315. be provided by setting the ``iszerofunc`` argument to a function that
  316. should return True if its argument is zero. The ADJ routine computes
  317. the determinant and uses that to detect singular matrices in addition
  318. to testing for zeros on the diagonal.
  319. See Also
  320. ========
  321. inverse_ADJ
  322. inverse_GE
  323. inverse_LU
  324. inverse_CH
  325. inverse_LDL
  326. Raises
  327. ======
  328. ValueError
  329. If the determinant of the matrix is zero.
  330. """
  331. from sympy.matrices import diag, SparseMatrix
  332. if method is None:
  333. method = 'LDL' if isinstance(M, SparseMatrix) else 'GE'
  334. if try_block_diag:
  335. blocks = M.get_diag_blocks()
  336. r = []
  337. for block in blocks:
  338. r.append(block.inv(method=method, iszerofunc=iszerofunc))
  339. return diag(*r)
  340. if method == "GE":
  341. rv = M.inverse_GE(iszerofunc=iszerofunc)
  342. elif method == "LU":
  343. rv = M.inverse_LU(iszerofunc=iszerofunc)
  344. elif method == "ADJ":
  345. rv = M.inverse_ADJ(iszerofunc=iszerofunc)
  346. elif method == "CH":
  347. rv = M.inverse_CH(iszerofunc=iszerofunc)
  348. elif method == "LDL":
  349. rv = M.inverse_LDL(iszerofunc=iszerofunc)
  350. elif method == "QR":
  351. rv = M.inverse_QR(iszerofunc=iszerofunc)
  352. elif method == "BLOCK":
  353. rv = M.inverse_BLOCK(iszerofunc=iszerofunc)
  354. else:
  355. raise ValueError("Inversion method unrecognized")
  356. return M._new(rv)