extract.py 4.6 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171
  1. """Functions to extract parts of sparse matrices
  2. """
  3. from __future__ import division, print_function, absolute_import
  4. __docformat__ = "restructuredtext en"
  5. __all__ = ['find', 'tril', 'triu']
  6. from .coo import coo_matrix
  7. def find(A):
  8. """Return the indices and values of the nonzero elements of a matrix
  9. Parameters
  10. ----------
  11. A : dense or sparse matrix
  12. Matrix whose nonzero elements are desired.
  13. Returns
  14. -------
  15. (I,J,V) : tuple of arrays
  16. I,J, and V contain the row indices, column indices, and values
  17. of the nonzero matrix entries.
  18. Examples
  19. --------
  20. >>> from scipy.sparse import csr_matrix, find
  21. >>> A = csr_matrix([[7.0, 8.0, 0],[0, 0, 9.0]])
  22. >>> find(A)
  23. (array([0, 0, 1], dtype=int32), array([0, 1, 2], dtype=int32), array([ 7., 8., 9.]))
  24. """
  25. A = coo_matrix(A, copy=True)
  26. A.sum_duplicates()
  27. # remove explicit zeros
  28. nz_mask = A.data != 0
  29. return A.row[nz_mask], A.col[nz_mask], A.data[nz_mask]
  30. def tril(A, k=0, format=None):
  31. """Return the lower triangular portion of a matrix in sparse format
  32. Returns the elements on or below the k-th diagonal of the matrix A.
  33. - k = 0 corresponds to the main diagonal
  34. - k > 0 is above the main diagonal
  35. - k < 0 is below the main diagonal
  36. Parameters
  37. ----------
  38. A : dense or sparse matrix
  39. Matrix whose lower trianglar portion is desired.
  40. k : integer : optional
  41. The top-most diagonal of the lower triangle.
  42. format : string
  43. Sparse format of the result, e.g. format="csr", etc.
  44. Returns
  45. -------
  46. L : sparse matrix
  47. Lower triangular portion of A in sparse format.
  48. See Also
  49. --------
  50. triu : upper triangle in sparse format
  51. Examples
  52. --------
  53. >>> from scipy.sparse import csr_matrix, tril
  54. >>> A = csr_matrix([[1, 2, 0, 0, 3], [4, 5, 0, 6, 7], [0, 0, 8, 9, 0]],
  55. ... dtype='int32')
  56. >>> A.toarray()
  57. array([[1, 2, 0, 0, 3],
  58. [4, 5, 0, 6, 7],
  59. [0, 0, 8, 9, 0]])
  60. >>> tril(A).toarray()
  61. array([[1, 0, 0, 0, 0],
  62. [4, 5, 0, 0, 0],
  63. [0, 0, 8, 0, 0]])
  64. >>> tril(A).nnz
  65. 4
  66. >>> tril(A, k=1).toarray()
  67. array([[1, 2, 0, 0, 0],
  68. [4, 5, 0, 0, 0],
  69. [0, 0, 8, 9, 0]])
  70. >>> tril(A, k=-1).toarray()
  71. array([[0, 0, 0, 0, 0],
  72. [4, 0, 0, 0, 0],
  73. [0, 0, 0, 0, 0]])
  74. >>> tril(A, format='csc')
  75. <3x5 sparse matrix of type '<class 'numpy.int32'>'
  76. with 4 stored elements in Compressed Sparse Column format>
  77. """
  78. # convert to COOrdinate format where things are easy
  79. A = coo_matrix(A, copy=False)
  80. mask = A.row + k >= A.col
  81. return _masked_coo(A, mask).asformat(format)
  82. def triu(A, k=0, format=None):
  83. """Return the upper triangular portion of a matrix in sparse format
  84. Returns the elements on or above the k-th diagonal of the matrix A.
  85. - k = 0 corresponds to the main diagonal
  86. - k > 0 is above the main diagonal
  87. - k < 0 is below the main diagonal
  88. Parameters
  89. ----------
  90. A : dense or sparse matrix
  91. Matrix whose upper trianglar portion is desired.
  92. k : integer : optional
  93. The bottom-most diagonal of the upper triangle.
  94. format : string
  95. Sparse format of the result, e.g. format="csr", etc.
  96. Returns
  97. -------
  98. L : sparse matrix
  99. Upper triangular portion of A in sparse format.
  100. See Also
  101. --------
  102. tril : lower triangle in sparse format
  103. Examples
  104. --------
  105. >>> from scipy.sparse import csr_matrix, triu
  106. >>> A = csr_matrix([[1, 2, 0, 0, 3], [4, 5, 0, 6, 7], [0, 0, 8, 9, 0]],
  107. ... dtype='int32')
  108. >>> A.toarray()
  109. array([[1, 2, 0, 0, 3],
  110. [4, 5, 0, 6, 7],
  111. [0, 0, 8, 9, 0]])
  112. >>> triu(A).toarray()
  113. array([[1, 2, 0, 0, 3],
  114. [0, 5, 0, 6, 7],
  115. [0, 0, 8, 9, 0]])
  116. >>> triu(A).nnz
  117. 8
  118. >>> triu(A, k=1).toarray()
  119. array([[0, 2, 0, 0, 3],
  120. [0, 0, 0, 6, 7],
  121. [0, 0, 0, 9, 0]])
  122. >>> triu(A, k=-1).toarray()
  123. array([[1, 2, 0, 0, 3],
  124. [4, 5, 0, 6, 7],
  125. [0, 0, 8, 9, 0]])
  126. >>> triu(A, format='csc')
  127. <3x5 sparse matrix of type '<class 'numpy.int32'>'
  128. with 8 stored elements in Compressed Sparse Column format>
  129. """
  130. # convert to COOrdinate format where things are easy
  131. A = coo_matrix(A, copy=False)
  132. mask = A.row + k <= A.col
  133. return _masked_coo(A, mask).asformat(format)
  134. def _masked_coo(A, mask):
  135. row = A.row[mask]
  136. col = A.col[mask]
  137. data = A.data[mask]
  138. return coo_matrix((data, (row, col)), shape=A.shape, dtype=A.dtype)