SciPy

scipy.sparse.coo_matrix

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

A sparse matrix in COOrdinate format.

Also known as the ‘ijv’ or ‘triplet’ format.

This can be instantiated in several ways:
coo_matrix(D)
with a dense matrix D
coo_matrix(S)
with another sparse matrix S (equivalent to S.tocoo())
coo_matrix((M, N), [dtype])
to construct an empty matrix with shape (M, N) dtype is optional, defaulting to dtype=’d’.
coo_matrix((data, (i, j)), [shape=(M, N)])
to construct from three arrays:
  1. data[:] the entries of the matrix, in any order
  2. i[:] the row indices of the matrix entries
  3. j[:] the column indices of the matrix entries

Where A[i[k], j[k]] = data[k]. When shape is not specified, it is inferred from the index arrays

Notes

Sparse matrices can be used in arithmetic operations: they support addition, subtraction, multiplication, division, and matrix power.

Advantages of the COO format
  • facilitates fast conversion among sparse formats
  • permits duplicate entries (see example)
  • very fast conversion to and from CSR/CSC formats
Disadvantages of the COO format
  • does not directly support:
    • arithmetic operations
    • slicing
Intended Usage
  • COO is a fast format for constructing sparse matrices
  • Once a matrix has been constructed, convert to CSR or CSC format for fast arithmetic and matrix vector operations
  • By default when converting to CSR or CSC format, duplicate (i,j) entries will be summed together. This facilitates efficient construction of finite element matrices and the like. (see example)

Examples

>>> # Constructing an empty matrix
>>> from scipy.sparse import coo_matrix
>>> coo_matrix((3, 4), dtype=np.int8).toarray()
array([[0, 0, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 0]], dtype=int8)
>>> # Constructing a matrix using ijv format
>>> row  = np.array([0, 3, 1, 0])
>>> col  = np.array([0, 3, 1, 2])
>>> data = np.array([4, 5, 7, 9])
>>> coo_matrix((data, (row, col)), shape=(4, 4)).toarray()
array([[4, 0, 9, 0],
       [0, 7, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 5]])
>>> # Constructing a matrix with duplicate indices
>>> row  = np.array([0, 0, 1, 3, 1, 0, 0])
>>> col  = np.array([0, 2, 1, 3, 1, 0, 0])
>>> data = np.array([1, 1, 1, 1, 1, 1, 1])
>>> coo = coo_matrix((data, (row, col)), shape=(4, 4))
>>> # Duplicate indices are maintained until implicitly or explicitly summed
>>> np.max(coo.data)
1
>>> coo.toarray()
array([[3, 0, 1, 0],
       [0, 2, 0, 0],
       [0, 0, 0, 0],
       [0, 0, 0, 1]])
Attributes:
dtype : dtype

Data type of the matrix

shape : 2-tuple

Get shape of a matrix.

ndim : int

Number of dimensions (this is always 2)

nnz

Number of stored values, including explicit zeros.

data

COO format data array of the matrix

row

COO format row index array of the matrix

col

COO format column index array of the matrix

Methods

arcsin() Element-wise arcsin.
arcsinh() Element-wise arcsinh.
arctan() Element-wise arctan.
arctanh() Element-wise arctanh.
argmax([axis, out]) Return indices of maximum elements along an axis.
argmin([axis, out]) Return indices of minimum elements along an axis.
asformat(format[, copy]) Return this matrix in the passed format.
asfptype() Upcast matrix to a floating point format (if necessary)
astype(dtype[, casting, copy]) Cast the matrix elements to a specified type.
ceil() Element-wise ceil.
conj([copy]) Element-wise complex conjugation.
conjugate([copy]) Element-wise complex conjugation.
copy() Returns a copy of this matrix.
count_nonzero() Number of non-zero entries, equivalent to
deg2rad() Element-wise deg2rad.
diagonal([k]) Returns the k-th diagonal of the matrix.
dot(other) Ordinary dot product
eliminate_zeros() Remove zero entries from the matrix
expm1() Element-wise expm1.
floor() Element-wise floor.
getH() Return the Hermitian transpose of this matrix.
get_shape() Get shape of a matrix.
getcol(j) Returns a copy of column j of the matrix, as an (m x 1) sparse matrix (column vector).
getformat() Format of a matrix representation as a string.
getmaxprint() Maximum number of elements to display when printed.
getnnz([axis]) Number of stored values, including explicit zeros.
getrow(i) Returns a copy of row i of the matrix, as a (1 x n) sparse matrix (row vector).
log1p() Element-wise log1p.
max([axis, out]) Return the maximum of the matrix or maximum along an axis.
maximum(other) Element-wise maximum between this and another matrix.
mean([axis, dtype, out]) Compute the arithmetic mean along the specified axis.
min([axis, out]) Return the minimum of the matrix or maximum along an axis.
minimum(other) Element-wise minimum between this and another matrix.
multiply(other) Point-wise multiplication by another matrix
nonzero() nonzero indices
power(n[, dtype]) This function performs element-wise power.
rad2deg() Element-wise rad2deg.
reshape(self, shape[, order, copy]) Gives a new shape to a sparse matrix without changing its data.
resize(*shape) Resize the matrix in-place to dimensions given by shape
rint() Element-wise rint.
set_shape(shape) See reshape.
setdiag(values[, k]) Set diagonal or off-diagonal elements of the array.
sign() Element-wise sign.
sin() Element-wise sin.
sinh() Element-wise sinh.
sqrt() Element-wise sqrt.
sum([axis, dtype, out]) Sum the matrix elements over a given axis.
sum_duplicates() Eliminate duplicate matrix entries by adding them together
tan() Element-wise tan.
tanh() Element-wise tanh.
toarray([order, out]) See the docstring for spmatrix.toarray.
tobsr([blocksize, copy]) Convert this matrix to Block Sparse Row format.
tocoo([copy]) Convert this matrix to COOrdinate format.
tocsc([copy]) Convert this matrix to Compressed Sparse Column format
tocsr([copy]) Convert this matrix to Compressed Sparse Row format
todense([order, out]) Return a dense matrix representation of this matrix.
todia([copy]) Convert this matrix to sparse DIAgonal format.
todok([copy]) Convert this matrix to Dictionary Of Keys format.
tolil([copy]) Convert this matrix to LInked List format.
transpose([axes, copy]) Reverses the dimensions of the sparse matrix.
trunc() Element-wise trunc.

Previous topic

scipy.sparse.bsr_matrix.trunc

Next topic

scipy.sparse.coo_matrix.arcsin