1

2

3

4

5

6

7

8

9

10

11

12

13

14

15

16

17

18

19

20

21

22

23

24

25

26

27

28

29

30

31

32

33

34

35

36

37

38

39

40

41

42

43

44

45

46

47

48

49

50

51

52

53

54

55

56

57

58

59

60

61

62

63

64

65

66

67

68

69

70

71

72

73

74

75

76

77

78

79

80

81

82

83

84

85

86

87

88

89

90

91

92

93

94

95

96

97

98

99

100

101

102

103

104

105

106

107

108

109

110

111

112

113

114

115

116

117

118

119

120

121

122

123

124

125

126

127

128

129

130

131

132

133

134

135

136

137

138

139

140

141

142

143

144

145

146

147

148

149

150

151

152

153

154

155

156

157

158

159

160

161

162

163

164

165

166

167

168

169

170

171

"""Functions to extract parts of sparse matrices 

""" 

 

from __future__ import division, print_function, absolute_import 

 

__docformat__ = "restructuredtext en" 

 

__all__ = ['find', 'tril', 'triu'] 

 

 

from .coo import coo_matrix 

 

 

def find(A): 

"""Return the indices and values of the nonzero elements of a matrix 

 

Parameters 

---------- 

A : dense or sparse matrix 

Matrix whose nonzero elements are desired. 

 

Returns 

------- 

(I,J,V) : tuple of arrays 

I,J, and V contain the row indices, column indices, and values 

of the nonzero matrix entries. 

 

 

Examples 

-------- 

>>> from scipy.sparse import csr_matrix, find 

>>> A = csr_matrix([[7.0, 8.0, 0],[0, 0, 9.0]]) 

>>> find(A) 

(array([0, 0, 1], dtype=int32), array([0, 1, 2], dtype=int32), array([ 7., 8., 9.])) 

 

""" 

 

A = coo_matrix(A, copy=True) 

A.sum_duplicates() 

# remove explicit zeros 

nz_mask = A.data != 0 

return A.row[nz_mask], A.col[nz_mask], A.data[nz_mask] 

 

 

def tril(A, k=0, format=None): 

"""Return the lower triangular portion of a matrix in sparse format 

 

Returns the elements on or below the k-th diagonal of the matrix A. 

- k = 0 corresponds to the main diagonal 

- k > 0 is above the main diagonal 

- k < 0 is below the main diagonal 

 

Parameters 

---------- 

A : dense or sparse matrix 

Matrix whose lower trianglar portion is desired. 

k : integer : optional 

The top-most diagonal of the lower triangle. 

format : string 

Sparse format of the result, e.g. format="csr", etc. 

 

Returns 

------- 

L : sparse matrix 

Lower triangular portion of A in sparse format. 

 

See Also 

-------- 

triu : upper triangle in sparse format 

 

Examples 

-------- 

>>> from scipy.sparse import csr_matrix, tril 

>>> A = csr_matrix([[1, 2, 0, 0, 3], [4, 5, 0, 6, 7], [0, 0, 8, 9, 0]], 

... dtype='int32') 

>>> A.toarray() 

array([[1, 2, 0, 0, 3], 

[4, 5, 0, 6, 7], 

[0, 0, 8, 9, 0]]) 

>>> tril(A).toarray() 

array([[1, 0, 0, 0, 0], 

[4, 5, 0, 0, 0], 

[0, 0, 8, 0, 0]]) 

>>> tril(A).nnz 

4 

>>> tril(A, k=1).toarray() 

array([[1, 2, 0, 0, 0], 

[4, 5, 0, 0, 0], 

[0, 0, 8, 9, 0]]) 

>>> tril(A, k=-1).toarray() 

array([[0, 0, 0, 0, 0], 

[4, 0, 0, 0, 0], 

[0, 0, 0, 0, 0]]) 

>>> tril(A, format='csc') 

<3x5 sparse matrix of type '<class 'numpy.int32'>' 

with 4 stored elements in Compressed Sparse Column format> 

 

""" 

 

# convert to COOrdinate format where things are easy 

A = coo_matrix(A, copy=False) 

mask = A.row + k >= A.col 

return _masked_coo(A, mask).asformat(format) 

 

 

def triu(A, k=0, format=None): 

"""Return the upper triangular portion of a matrix in sparse format 

 

Returns the elements on or above the k-th diagonal of the matrix A. 

- k = 0 corresponds to the main diagonal 

- k > 0 is above the main diagonal 

- k < 0 is below the main diagonal 

 

Parameters 

---------- 

A : dense or sparse matrix 

Matrix whose upper trianglar portion is desired. 

k : integer : optional 

The bottom-most diagonal of the upper triangle. 

format : string 

Sparse format of the result, e.g. format="csr", etc. 

 

Returns 

------- 

L : sparse matrix 

Upper triangular portion of A in sparse format. 

 

See Also 

-------- 

tril : lower triangle in sparse format 

 

Examples 

-------- 

>>> from scipy.sparse import csr_matrix, triu 

>>> A = csr_matrix([[1, 2, 0, 0, 3], [4, 5, 0, 6, 7], [0, 0, 8, 9, 0]], 

... dtype='int32') 

>>> A.toarray() 

array([[1, 2, 0, 0, 3], 

[4, 5, 0, 6, 7], 

[0, 0, 8, 9, 0]]) 

>>> triu(A).toarray() 

array([[1, 2, 0, 0, 3], 

[0, 5, 0, 6, 7], 

[0, 0, 8, 9, 0]]) 

>>> triu(A).nnz 

8 

>>> triu(A, k=1).toarray() 

array([[0, 2, 0, 0, 3], 

[0, 0, 0, 6, 7], 

[0, 0, 0, 9, 0]]) 

>>> triu(A, k=-1).toarray() 

array([[1, 2, 0, 0, 3], 

[4, 5, 0, 6, 7], 

[0, 0, 8, 9, 0]]) 

>>> triu(A, format='csc') 

<3x5 sparse matrix of type '<class 'numpy.int32'>' 

with 8 stored elements in Compressed Sparse Column format> 

 

""" 

 

# convert to COOrdinate format where things are easy 

A = coo_matrix(A, copy=False) 

mask = A.row + k <= A.col 

return _masked_coo(A, mask).asformat(format) 

 

 

def _masked_coo(A, mask): 

row = A.row[mask] 

col = A.col[mask] 

data = A.data[mask] 

return coo_matrix((data, (row, col)), shape=A.shape, dtype=A.dtype)