_extract.py 4.5 KB

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