/
csc.py
227 lines (174 loc) · 6.93 KB
/
csc.py
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
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
try:
import scipy.sparse
_scipy_available = True
except ImportError:
_scipy_available = False
import cupy
from cupy import cusparse
from cupy.sparse import compressed
class csc_matrix(compressed._compressed_sparse_matrix):
"""Compressed Sparse Column matrix.
Now it has only part of initializer formats:
``csc_matrix(D)``
``D`` is a rank-2 :class:`cupy.ndarray`.
``csc_matrix(S)``
``S`` is another sparse matrix. It is equivalent to ``S.tocsc()``.
``csc_matrix((M, N), [dtype])``
It constructs an empty matrix whose shape is ``(M, N)``. Default dtype
is flat64.
``csc_matrix((data, indices, indptr))``
All ``data``, ``indices`` and ``indptr`` are one-dimenaional
:class:`cupy.ndarray`.
Args:
arg1: Arguments for the initializer.
shape (tuple): Shape of a matrix. Its length must be two.
dtype: Data type. It must be an argument of :class:`numpy.dtype`.
copy (bool): If ``True``, copies of given arrays are always used.
.. seealso::
:class:`scipy.sparse.csc_matrix`
"""
format = 'csc'
def get(self, stream=None):
"""Returns a copy of the array on host memory.
.. warning::
You need to install SciPy to use this method.
Args:
stream (cupy.cuda.Stream): CUDA stream object. If it is given, the
copy runs asynchronously. Otherwise, the copy is synchronous.
Returns:
scipy.sparse.csc_matrix: Copy of the array on host memory.
"""
if not _scipy_available:
raise RuntimeError('scipy is not available')
data = self.data.get(stream)
indices = self.indices.get(stream)
indptr = self.indptr.get(stream)
return scipy.sparse.csc_matrix(
(data, indices, indptr), shape=self._shape)
def _convert_dense(self, x):
m = cusparse.dense2csc(x)
return m.data, m.indices, m.indptr
def _swap(self, x, y):
return (y, x)
# TODO(unno): Implement __getitem__
def __mul__(self, other):
if cupy.isscalar(other):
return self._with_data(self.data * other)
elif cupy.sparse.isspmatrix_csr(other):
return cusparse.csrgemm(self.T, other, transa=True)
elif isspmatrix_csc(other):
return cusparse.csrgemm(self.T, other.T, transa=True, transb=True)
elif cupy.sparse.isspmatrix(other):
return cusparse.csrgemm(self.T, other.tocsr(), transa=True)
elif cupy.sparse.base.isdense(other):
if other.ndim == 0:
return self._with_data(self.data * other)
elif other.ndim == 1:
return cusparse.csrmv(
self.T, cupy.asfortranarray(other), transa=True)
elif other.ndim == 2:
return cusparse.csrmm2(
self.T, cupy.asfortranarray(other), transa=True)
else:
raise ValueError('could not interpret dimensions')
else:
return NotImplemented
# TODO(unno): Implement argmax
# TODO(unno): Implement argmin
# TODO(unno): Implement check_format
# TODO(unno): Implement diagonal
# TODO(unno): Implement eliminate_zeros
# TODO(unno): Implement max
# TODO(unno): Implement maximum
# TODO(unno): Implement min
# TODO(unno): Implement minimum
# TODO(unno): Implement multiply
# TODO(unno): Implement prune
# TODO(unno): Implement reshape
def sort_indices(self):
"""Sorts the indices of the matrix in place."""
cusparse.cscsort(self)
def toarray(self, order=None, out=None):
"""Returns a dense matrix representing the same value.
Args:
order ({'C', 'F', None}): Whether to store data in C (row-major)
order or F (column-major) order. Default is C-order.
out: Not supported.
Returns:
cupy.ndarray: Dense array representing the same matrix.
.. seealso:: :func:`cupy.sparse.csc_array.toarray`
"""
if order is None:
order = 'C'
if self.nnz == 0:
return cupy.zeros(shape=self.shape, dtype=self.dtype, order=order)
self.sum_duplicates()
# csc2dense and csr2dense returns F-contiguous array.
if order == 'C':
# To return C-contiguous array, it uses transpose.
return cusparse.csr2dense(self.T).T
elif order == 'F':
return cusparse.csc2dense(self)
else:
raise TypeError('order not understood')
def _add_sparse(self, other, alpha, beta):
self.sum_duplicates()
other.sum_duplicates()
return cusparse.csrgeam(self.T, other.tocsc().T, alpha, beta).T
# TODO(unno): Implement tobsr
def tocoo(self, copy=False):
"""Converts the matrix to COOdinate format.
Args:
copy (bool): If ``False``, it shares data arrays as much as
possible.
Returns:
cupy.sparse.coo_matrix: Converted matrix.
"""
return self.T.tocoo(copy).T
def tocsc(self, copy=None):
"""Converts the matrix to Compressed Sparse Column format.
Args:
copy (bool): If ``False``, the method returns itself.
Otherwise it makes a copy of the matrix.
Returns:
cupy.sparse.csc_matrix: Converted matrix.
"""
if copy:
return self.copy()
else:
return self
def tocsr(self, copy=False):
"""Converts the matrix to Compressed Sparse Row format.
Args:
copy (bool): If ``False``, it shares data arrays as much as
possible. Actually this option is ignored because all
arrays in a matrix cannot be shared in csr to csc conversion.
Returns:
cupy.sparse.csr_matrix: Converted matrix.
"""
return self.T.tocsc(copy=False).T
# TODO(unno): Implement todia
# TODO(unno): Implement todok
# TODO(unno): Implement tolil
def transpose(self, axes=None, copy=False):
"""Returns a transpose matrix.
Args:
axes: This option is not supported.
copy (bool): If ``True``, a returned matrix shares no data.
Otherwise, it shared data arrays as much as possible.
Returns:
cupy.sparse.spmatrix: Transpose matrix.
"""
if axes is not None:
raise ValueError(
'Sparse matrices do not support an \'axes\' parameter because '
'swapping dimensions is the only logical permutation.')
shape = self.shape[1], self.shape[0]
return cupy.sparse.csr_matrix(
(self.data, self.indices, self.indptr), shape=shape, copy=copy)
def isspmatrix_csc(x):
"""Checks if a given matrix is of CSC format.
Returns:
bool: Returns if ``x`` is :class:`cupy.sparse.csc_matrix`.
"""
return isinstance(x, csc_matrix)