cupy.sparse.csr_matrix

class cupy.sparse.csr_matrix(arg1, shape=None, dtype=None, copy=False)[source]

Compressed Sparse Row matrix.

Now it has only part of initializer formats:

csr_matrix(D)
D is a rank-2 cupy.ndarray.
csr_matrix(S)
S is another sparse matrix. It is equivalent to S.tocsr().
csr_matrix((M, N), [dtype])
It constructs an empty matrix whose shape is (M, N). Default dtype is float64.
csr_matrix((data, indices, indptr))
All data, indices and indptr are one-dimenaional cupy.ndarray.
Parameters:
  • 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 numpy.dtype.
  • copy (bool) – If True, copies of given arrays are always used.

Methods

__getitem__(slices)[source]
__len__()[source]
__iter__()[source]
arcsin()[source]

Elementwise arcsin.

arcsinh()[source]

Elementwise arcsinh.

arctan()[source]

Elementwise arctan.

arctanh()[source]

Elementwise arctanh.

asformat(format)[source]

Return this matrix in a given sparse format.

Parameters:format (str or None) – Format you need.
asfptype()[source]

Upcasts matrix to a floating point format.

When the matrix has floating point type, the method returns itself. Otherwise it makes a copy with floating point type and the same format.

Returns:A matrix with float type.
Return type:cupy.sparse.spmatrix
astype(t)[source]

Casts the array to given data type.

Parameters:dtype – Type specifier.
Returns:A copy of the array with a given type.
ceil()[source]

Elementwise ceil.

conj()[source]
conjugate()[source]
copy()[source]

Returns a copy of this matrix.

count_nonzero()[source]

Returns number of non-zero entries.

Note

This method counts the actual number of non-zero entories, which does not include explicit zero entries. Instead nnz returns the number of entries including explicit zeros.

Returns:Number of non-zero entries.
deg2rad()[source]

Elementwise deg2rad.

diagonal()[source]

Returns the main diagonal of the matrix

dot(other)[source]

Ordinary dot product

eliminate_zeros()[source]

Removes zero entories in place.

expm1()[source]

Elementwise expm1.

floor()[source]

Elementwise floor.

get(stream=None)[source]

Returns a copy of the array on host memory.

Parameters:stream (cupy.cuda.Stream) – CUDA stream object. If it is given, the copy runs asynchronously. Otherwise, the copy is synchronous.
Returns:Copy of the array on host memory.
Return type:scipy.sparse.csr_matrix
getH()[source]
get_shape()[source]

Returns the shape of the matrix.

Returns:Shape of the matrix.
Return type:tuple
getformat()[source]
getmaxprint()[source]
getnnz(axis=None)[source]

Returns the number of stored values, including explicit zeros.

Parameters:axis – Not supported yet.
Returns:The number of stored values.
Return type:int
log1p()[source]

Elementwise log1p.

maximum(other)[source]
minimum(other)[source]
multiply(other)[source]

Point-wise multiplication by another matrix

power(n, dtype=None)[source]

Elementwise power function.

Parameters:
  • n – Exponent.
  • dtype – Type specifier.
rad2deg()[source]

Elementwise rad2deg.

reshape(shape, order='C')[source]

Gives a new shape to a sparse matrix without changing its data.

rint()[source]

Elementwise rint.

set_shape(shape)[source]
sign()[source]

Elementwise sign.

sin()[source]

Elementwise sin.

sinh()[source]

Elementwise sinh.

sort_indices()[source]

Sorts the indices of the matrix in place.

sqrt()[source]

Elementwise sqrt.

sum(axis=None, dtype=None, out=None)[source]

Sums the matrix elements over a given axis.

Parameters:
  • axis (int or None) – Axis along which the sum is comuted. If it is None, it computes the sum of all the elements. Select from {None, 0, 1, -2, -1}.
  • dtype – The type of returned matrix. If it is not specified, type of the array is used.
  • out (cupy.ndarray) – Output matrix.
Returns:

Summed array.

Return type:

cupy.ndarray

sum_duplicates()[source]
tan()[source]

Elementwise tan.

tanh()[source]

Elementwise tanh.

toarray(order=None, out=None)[source]

Returns a dense matrix representing the same value.

Parameters:
  • 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:

Dense array representing the same matrix.

Return type:

cupy.ndarray

tobsr(blocksize=None, copy=False)[source]

Convert this matrix to Block Sparse Row format.

tocoo(copy=False)[source]

Converts the matrix to COOdinate format.

Parameters:copy (bool) – If False, it shares data arrays as much as possible.
Returns:Converted matrix.
Return type:cupy.sparse.coo_matrix
tocsc(copy=False)[source]

Converts the matrix to Compressed Sparse Column format.

Parameters: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:Converted matrix.
Return type:cupy.sparse.csc_matrix
tocsr(copy=False)[source]

Converts the matrix to Compressed Sparse Row format.

Parameters:copy (bool) – If False, the method returns itself. Otherwise it makes a copy of the matrix.
Returns:Converted matrix.
Return type:cupy.sparse.csr_matrix
todense(order=None, out=None)[source]

Return a dense matrix representation of this matrix.

todia(copy=False)[source]

Convert this matrix to sparse DIAgonal format.

todok(copy=False)[source]

Convert this matrix to Dictionary Of Keys format.

tolil(copy=False)[source]

Convert this matrix to LInked List format.

transpose(axes=None, copy=False)[source]

Returns a transpose matrix.

Parameters:
  • 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:

Transpose matrix.

Return type:

cupy.sparse.spmatrix

trunc()[source]

Elementwise trunc.

__eq__(other)[source]

Return self==value.

__ne__(other)[source]

Return self!=value.

__lt__(other)[source]

Return self<value.

__le__(other)[source]

Return self<=value.

__gt__(other)[source]

Return self>value.

__ge__(other)[source]

Return self>=value.

__nonzero__()[source]
__bool__()[source]

Attributes

A

Dense ndarray representation of this matrix.

This property is equivalent to toarray() method.

H
T
device

CUDA device on which this array resides.

dtype

Data type of the matrix.

format = 'csr'
has_canonical_format
ndim
nnz
shape
size