inverse.py 13 KB

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