_extract.py 4.9 KB

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