scipy.sparse.coo_matrix

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

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, ij), [shape=(M, N)])
The arguments ‘data’ and ‘ij’ represent three arrays:
  1. data[:] the entries of the matrix, in any order
  2. ij[0][:] the row indices of the matrix entries
  3. ij[1][:] the column indices of the matrix entries

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

Notes

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

>>> from scipy.sparse import *
>>> from scipy import *
>>> coo_matrix( (3,4), dtype=int8 ).todense()
matrix([[0, 0, 0, 0],
        [0, 0, 0, 0],
        [0, 0, 0, 0]], dtype=int8)
>>> row  = array([0,3,1,0])
>>> col  = array([0,3,1,2])
>>> data = array([4,5,7,9])
>>> coo_matrix( (data,(row,col)), shape=(4,4) ).todense()
matrix([[4, 0, 9, 0],
        [0, 7, 0, 0],
        [0, 0, 0, 0],
        [0, 0, 0, 5]])
>>> # example with duplicates
>>> row  = array([0,0,1,3,1,0,0])
>>> col  = array([0,2,1,3,1,0,0])
>>> data = array([1,1,1,1,1,1,1])
>>> coo_matrix( (data,(row,col)), shape=(4,4)).todense()
matrix([[3, 0, 1, 0],
        [0, 2, 0, 0],
        [0, 0, 0, 0],
        [0, 0, 0, 1]])

Methods

asformat
asfptype
astype
conj
conjugate
copy Generic (shallow and deep) copying operations.
diagonal
dot
getH
get_shape
getcol
getformat
getmaxprint
getnnz
getrow
mean
multiply
nonzero
reshape
set_shape
setdiag
sum
toarray
tobsr
tocoo
tocsc
tocsr
todense
todia
todok
tolil
transpose

This Page