Sparse matrices (scipy.sparse)

Sparse Matrices

SciPy 2-D sparse matrix package.

Original code by Travis Oliphant. Modified and extended by Ed Schofield, Robert Cimrman, and Nathan Bell.

There are seven available sparse matrix types:
  1. csc_matrix: Compressed Sparse Column format
  2. csr_matrix: Compressed Sparse Row format
  3. bsr_matrix: Block Sparse Row format
  4. lil_matrix: List of Lists format
  5. dok_matrix: Dictionary of Keys format
  6. coo_matrix: COOrdinate format (aka IJV, triplet format)
  7. dia_matrix: DIAgonal format

To construct a matrix efficiently, use either lil_matrix (recommended) or dok_matrix. The lil_matrix class supports basic slicing and fancy indexing with a similar syntax to NumPy arrays. As illustrated below, the COO format may also be used to efficiently construct matrices.

To perform manipulations such as multiplication or inversion, first convert the matrix to either CSC or CSR format. The lil_matrix format is row-based, so conversion to CSR is efficient, whereas conversion to CSC is less so.

All conversions among the CSR, CSC, and COO formats are efficient, linear-time operations.

Example 1

Construct a 1000x1000 lil_matrix and add some values to it:

>>> from scipy.sparse import lil_matrix
>>> from scipy.sparse.linalg import spsolve
>>> from numpy.linalg import solve, norm
>>> from numpy.random import rand
>>> A = lil_matrix((1000, 1000))
>>> A[0, :100] = rand(100)
>>> A[1, 100:200] = A[0, :100]
>>> A.setdiag(rand(1000))

Now convert it to CSR format and solve A x = b for x:

>>> A = A.tocsr()
>>> b = rand(1000)
>>> x = spsolve(A, b)

Convert it to a dense matrix and solve, and check that the result is the same:

>>> x_ = solve(A.todense(), b)

Now we can compute norm of the error with:

>>> err = norm(x-x_)
>>> err < 1e-10
True

It should be small :)

Example 2

Construct a matrix in COO format:

>>> from scipy import sparse
>>> from numpy import array
>>> I = array([0,3,1,0])
>>> J = array([0,3,1,2])
>>> V = array([4,5,7,9])
>>> A = sparse.coo_matrix((V,(I,J)),shape=(4,4))

Notice that the indices do not need to be sorted.

Duplicate (i,j) entries are summed when converting to CSR or CSC.

>>> I = array([0,0,1,3,1,0,0])
>>> J = array([0,2,1,3,1,0,0])
>>> V = array([1,1,1,1,1,1,1])
>>> B = sparse.coo_matrix((V,(I,J)),shape=(4,4)).tocsr()

This is useful for constructing finite-element stiffness and mass matrices.

Further Details

CSR column indices are not necessarily sorted. Likewise for CSC row indices. Use the .sorted_indices() and .sort_indices() methods when sorted indices are required (e.g. when passing data to other libraries).

Package Contents

Modules

base Base class for sparse matrices
bsr Compressed Block Sparse Row matrix format
compressed Base class for sparse matrix formats using compressed storage
construct Functions to construct sparse matrices
coo A sparse matrix in COOrdinate or ‘triplet’ format
csc Compressed Sparse Column matrix format
csgraph Compressed Sparse graph algorithms
csr Compressed Sparse Row matrix format
data Base class for sparse matrice with a .data attribute
dia Sparse DIAgonal format
dok Dictionary Of Keys based matrix
extract Functions to extract parts of sparse matrices
lil LInked List sparse matrix class
linalg Sparse Linear Algebra
sparsetools sparsetools - a collection of routines for sparse matrix operations
spfuncs Functions that operate on sparse matrices
sputils Utility functions for sparse matrix module

Classes

SparseEfficiencyWarning
SparseWarning
bsr_matrix(arg1[, shape, dtype, copy, blocksize]) Block Sparse Row matrix
coo_matrix(arg1[, shape, dtype, copy]) A sparse matrix in COOrdinate format.
csc_matrix(arg1[, shape, dtype, copy]) Compressed Sparse Column matrix
csr_matrix(arg1[, shape, dtype, copy]) Compressed Sparse Row matrix
dia_matrix(arg1[, shape, dtype, copy]) Sparse matrix with DIAgonal storage
dok_matrix(arg1[, shape, dtype, copy]) Dictionary Of Keys based sparse matrix.
lil_matrix(arg1[, shape, dtype, copy]) Row-based linked list sparse matrix

Functions

bmat(blocks[, format, dtype]) Build a sparse matrix from sparse sub-blocks
cs_graph_components(x) Determine connected compoments of a graph stored as a compressed sparse row or column matrix.
eye(m, n[, k, dtype, format]) eye(m, n) returns a sparse (m x n) matrix where the k-th diagonal
find(A) Return the indices and values of the nonzero elements of a matrix
hstack(blocks[, format, dtype]) Stack sparse matrices horizontally (column wise)
identity(n[, dtype, format]) Identity matrix in sparse format
issparse(x)
isspmatrix(x)
isspmatrix_bsr(x)
isspmatrix_coo(x)
isspmatrix_csc(x)
isspmatrix_csr(x)
isspmatrix_dia(x)
isspmatrix_dok(x)
isspmatrix_lil(x)
kron(A, B[, format]) kronecker product of sparse matrices A and B
kronsum(A, B[, format]) kronecker sum of sparse matrices A and B
lil_diags
lil_eye
rand(m, n[, density, format, dtype]) Generate a sparse matrix of the given shape and density with uniformely distributed values.
spdiags(data, diags, m, n[, format]) Return a sparse matrix from diagonals.
tril(A[, k, format]) Return the lower triangular portion of a matrix in sparse format
triu(A[, k, format]) Return the upper triangular portion of a matrix in sparse format
vstack(blocks[, format, dtype]) Stack sparse matrices vertically (row wise)

Sparse matrix classes

csc_matrix(arg1[, shape, dtype, copy]) Compressed Sparse Column matrix
csr_matrix(arg1[, shape, dtype, copy]) Compressed Sparse Row matrix
bsr_matrix(arg1[, shape, dtype, copy, blocksize]) Block Sparse Row matrix
lil_matrix(arg1[, shape, dtype, copy]) Row-based linked list sparse matrix
dok_matrix(arg1[, shape, dtype, copy]) Dictionary Of Keys based sparse matrix.
coo_matrix(arg1[, shape, dtype, copy]) A sparse matrix in COOrdinate format.
dia_matrix(arg1[, shape, dtype, copy]) Sparse matrix with DIAgonal storage

Functions

Building sparse matrices:

eye(m, n[, k, dtype, format]) eye(m, n) returns a sparse (m x n) matrix where the k-th diagonal
identity(n[, dtype, format]) Identity matrix in sparse format
kron(A, B[, format]) kronecker product of sparse matrices A and B
kronsum(A, B[, format]) kronecker sum of sparse matrices A and B
spdiags(data, diags, m, n[, format]) Return a sparse matrix from diagonals.
tril(A[, k, format]) Return the lower triangular portion of a matrix in sparse format
triu(A[, k, format]) Return the upper triangular portion of a matrix in sparse format
bmat(blocks[, format, dtype]) Build a sparse matrix from sparse sub-blocks
hstack(blocks[, format, dtype]) Stack sparse matrices horizontally (column wise)
vstack(blocks[, format, dtype]) Stack sparse matrices vertically (row wise)

Identifying sparse matrices:

issparse(x)
isspmatrix(x)
isspmatrix_csc(x)
isspmatrix_csr(x)
isspmatrix_bsr(x)
isspmatrix_lil(x)
isspmatrix_dok(x)
isspmatrix_coo(x)
isspmatrix_dia(x)

Exceptions

exception scipy.sparse.SparseEfficiencyWarning
exception scipy.sparse.SparseWarning