Scipy 2D sparse matrix module.
Original code by Travis Oliphant. Modified and extended by Ed Schofield, Robert Cimrman, and Nathan Bell.
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.
Construct a 1000x1000 lil_matrix and add some values to it:
>>> from scipy import sparse, linsolve
>>> from numpy import linalg
>>> from numpy.random import rand
>>> A = sparse.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 = linsolve.spsolve(A, b)
Convert it to a dense matrix and solve, and check that the result is the same:
>>> x_ = linalg.solve(A.todense(), b)
Now we can compute norm of the error with:
>>> err = linalg.norm(x-x_)
>>> err < 1e-10
True
It should be small :)
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.
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).
csc_matrix | Compressed Sparse Column matrix |
csr_matrix | Compressed Sparse Row matrix |
bsr_matrix | Block Sparse Row matrix |
lil_matrix | Row-based linked list sparse matrix |
dok_matrix | Dictionary Of Keys based sparse matrix. |
coo_matrix | A sparse matrix in COOrdinate format. |
dia_matrix | Sparse matrix with DIAgonal storage |
Building sparse matrices:
eye (m, n[, k, dtype, format]) | eye(m, n) returns a sparse (m x n) matrix where the k-th diagonal is all ones and everything else is zeros. |
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 |
lil_eye ((r, c)[, k, dtype]) | Generate a lil_matrix of dimensions (r,c) with the k-th diagonal set to 1. |
lil_diags (diags, offsets, (m, n)[, dtype]) | Generate a lil_matrix with the given diagonals. |
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) |