sparse.construct

Functions to construct sparse matrices

Module Contents

Functions

spdiags(data,diags,m,n,format=None) Return a sparse matrix from diagonals.
diags(diagonals,offsets=0,shape=None,format=None,dtype=None) Construct a sparse matrix from diagonals.
identity(n,dtype=”d”,format=None) Identity matrix in sparse format
eye(m,n=None,k=0,dtype=float,format=None) Sparse matrix with ones on diagonal
kron(A,B,format=None) kronecker product of sparse matrices A and B
kronsum(A,B,format=None) kronecker sum of sparse matrices A and B
_compressed_sparse_stack(blocks,axis) Stacking fast path for CSR/CSC matrices
hstack(blocks,format=None,dtype=None) Stack sparse matrices horizontally (column wise)
vstack(blocks,format=None,dtype=None) Stack sparse matrices vertically (row wise)
bmat(blocks,format=None,dtype=None) Build a sparse matrix from sparse sub-blocks
block_diag(mats,format=None,dtype=None) Build a block diagonal sparse matrix from provided matrices.
random(m,n,density=0.01,format=”coo”,dtype=None,random_state=None,data_rvs=None) Generate a sparse matrix of the given shape and density with randomly
rand(m,n,density=0.01,format=”coo”,dtype=None,random_state=None) Generate a sparse matrix of the given shape and density with uniformly
spdiags(data, diags, m, n, format=None)

Return a sparse matrix from diagonals.

data : array_like
matrix diagonals stored row-wise
diags : diagonals to set
  • k = 0 the main diagonal
  • k > 0 the k-th upper diagonal
  • k < 0 the k-th lower diagonal
m, n : int
shape of the result
format : str, optional
Format of the result. By default (format=None) an appropriate sparse matrix format is returned. This choice is subject to change.

diags : more convenient form of this function dia_matrix : the sparse DIAgonal format.

>>> from scipy.sparse import spdiags
>>> data = np.array([[1, 2, 3, 4], [1, 2, 3, 4], [1, 2, 3, 4]])
>>> diags = np.array([0, -1, 2])
>>> spdiags(data, diags, 4, 4).toarray()
array([[1, 0, 3, 0],
       [1, 2, 0, 4],
       [0, 2, 3, 0],
       [0, 0, 3, 4]])
diags(diagonals, offsets=0, shape=None, format=None, dtype=None)

Construct a sparse matrix from diagonals.

diagonals : sequence of array_like
Sequence of arrays containing the matrix diagonals, corresponding to offsets.
offsets : sequence of int or an int, optional
Diagonals to set:
  • k = 0 the main diagonal (default)
  • k > 0 the k-th upper diagonal
  • k < 0 the k-th lower diagonal
shape : tuple of int, optional
Shape of the result. If omitted, a square matrix large enough to contain the diagonals is returned.
format : {“dia”, “csr”, “csc”, “lil”, …}, optional
Matrix format of the result. By default (format=None) an appropriate sparse matrix format is returned. This choice is subject to change.
dtype : dtype, optional
Data type of the matrix.

spdiags : construct matrix from diagonals

This function differs from spdiags in the way it handles off-diagonals.

The result from diags is the sparse equivalent of:

np.diag(diagonals[0], offsets[0])
+ ...
+ np.diag(diagonals[k], offsets[k])

Repeated diagonal offsets are disallowed.

New in version 0.11.

>>> from scipy.sparse import diags
>>> diagonals = [[1, 2, 3, 4], [1, 2, 3], [1, 2]]
>>> diags(diagonals, [0, -1, 2]).toarray()
array([[1, 0, 1, 0],
       [1, 2, 0, 2],
       [0, 2, 3, 0],
       [0, 0, 3, 4]])

Broadcasting of scalars is supported (but shape needs to be specified):

>>> diags([1, -2, 1], [-1, 0, 1], shape=(4, 4)).toarray()
array([[-2.,  1.,  0.,  0.],
       [ 1., -2.,  1.,  0.],
       [ 0.,  1., -2.,  1.],
       [ 0.,  0.,  1., -2.]])

If only one diagonal is wanted (as in numpy.diag), the following works as well:

>>> diags([1, 2, 3], 1).toarray()
array([[ 0.,  1.,  0.,  0.],
       [ 0.,  0.,  2.,  0.],
       [ 0.,  0.,  0.,  3.],
       [ 0.,  0.,  0.,  0.]])
identity(n, dtype="d", format=None)

Identity matrix in sparse format

Returns an identity matrix with shape (n,n) using a given sparse format and dtype.

n : int
Shape of the identity matrix.
dtype : dtype, optional
Data type of the matrix
format : str, optional
Sparse format of the result, e.g. format=”csr”, etc.
>>> from scipy.sparse import identity
>>> identity(3).toarray()
array([[ 1.,  0.,  0.],
       [ 0.,  1.,  0.],
       [ 0.,  0.,  1.]])
>>> identity(3, dtype='int8', format='dia')
<3x3 sparse matrix of type '<class 'numpy.int8'>'
        with 3 stored elements (1 diagonals) in DIAgonal format>
eye(m, n=None, k=0, dtype=float, format=None)

Sparse matrix with ones on diagonal

Returns a sparse (m x n) matrix where the k-th diagonal is all ones and everything else is zeros.

m : int
Number of rows in the matrix.
n : int, optional
Number of columns. Default: m.
k : int, optional
Diagonal to place ones on. Default: 0 (main diagonal).
dtype : dtype, optional
Data type of the matrix.
format : str, optional
Sparse format of the result, e.g. format=”csr”, etc.
>>> from scipy import sparse
>>> sparse.eye(3).toarray()
array([[ 1.,  0.,  0.],
       [ 0.,  1.,  0.],
       [ 0.,  0.,  1.]])
>>> sparse.eye(3, dtype=np.int8)
<3x3 sparse matrix of type '<class 'numpy.int8'>'
    with 3 stored elements (1 diagonals) in DIAgonal format>
kron(A, B, format=None)

kronecker product of sparse matrices A and B

A : sparse or dense matrix
first matrix of the product
B : sparse or dense matrix
second matrix of the product
format : str, optional
format of the result (e.g. “csr”)

kronecker product in a sparse matrix format

>>> from scipy import sparse
>>> A = sparse.csr_matrix(np.array([[0, 2], [5, 0]]))
>>> B = sparse.csr_matrix(np.array([[1, 2], [3, 4]]))
>>> sparse.kron(A, B).toarray()
array([[ 0,  0,  2,  4],
       [ 0,  0,  6,  8],
       [ 5, 10,  0,  0],
       [15, 20,  0,  0]])
>>> sparse.kron(A, [[1, 2], [3, 4]]).toarray()
array([[ 0,  0,  2,  4],
       [ 0,  0,  6,  8],
       [ 5, 10,  0,  0],
       [15, 20,  0,  0]])
kronsum(A, B, format=None)

kronecker sum of sparse matrices A and B

Kronecker sum of two sparse matrices is a sum of two Kronecker products kron(I_n,A) + kron(B,I_m) where A has shape (m,m) and B has shape (n,n) and I_m and I_n are identity matrices of shape (m,m) and (n,n) respectively.

A
square matrix
B
square matrix
format : str
format of the result (e.g. “csr”)

kronecker sum in a sparse matrix format

_compressed_sparse_stack(blocks, axis)

Stacking fast path for CSR/CSC matrices (i) vstack for CSR, (ii) hstack for CSC.

hstack(blocks, format=None, dtype=None)

Stack sparse matrices horizontally (column wise)

blocks
sequence of sparse matrices with compatible shapes
format : str
sparse format of the result (e.g. “csr”) by default an appropriate sparse matrix format is returned. This choice is subject to change.
dtype : dtype, optional
The data-type of the output matrix. If not given, the dtype is determined from that of blocks.

vstack : stack sparse matrices vertically (row wise)

>>> from scipy.sparse import coo_matrix, hstack
>>> A = coo_matrix([[1, 2], [3, 4]])
>>> B = coo_matrix([[5], [6]])
>>> hstack([A,B]).toarray()
array([[1, 2, 5],
       [3, 4, 6]])
vstack(blocks, format=None, dtype=None)

Stack sparse matrices vertically (row wise)

blocks
sequence of sparse matrices with compatible shapes
format : str, optional
sparse format of the result (e.g. “csr”) by default an appropriate sparse matrix format is returned. This choice is subject to change.
dtype : dtype, optional
The data-type of the output matrix. If not given, the dtype is determined from that of blocks.

hstack : stack sparse matrices horizontally (column wise)

>>> from scipy.sparse import coo_matrix, vstack
>>> A = coo_matrix([[1, 2], [3, 4]])
>>> B = coo_matrix([[5, 6]])
>>> vstack([A, B]).toarray()
array([[1, 2],
       [3, 4],
       [5, 6]])
bmat(blocks, format=None, dtype=None)

Build a sparse matrix from sparse sub-blocks

blocks : array_like
Grid of sparse matrices with compatible shapes. An entry of None implies an all-zero matrix.
format : {‘bsr’, ‘coo’, ‘csc’, ‘csr’, ‘dia’, ‘dok’, ‘lil’}, optional
The sparse format of the result (e.g. “csr”). By default an appropriate sparse matrix format is returned. This choice is subject to change.
dtype : dtype, optional
The data-type of the output matrix. If not given, the dtype is determined from that of blocks.

bmat : sparse matrix

block_diag, diags

>>> from scipy.sparse import coo_matrix, bmat
>>> A = coo_matrix([[1, 2], [3, 4]])
>>> B = coo_matrix([[5], [6]])
>>> C = coo_matrix([[7]])
>>> bmat([[A, B], [None, C]]).toarray()
array([[1, 2, 5],
       [3, 4, 6],
       [0, 0, 7]])
>>> bmat([[A, None], [None, C]]).toarray()
array([[1, 2, 0],
       [3, 4, 0],
       [0, 0, 7]])
block_diag(mats, format=None, dtype=None)

Build a block diagonal sparse matrix from provided matrices.

mats : sequence of matrices
Input matrices.
format : str, optional
The sparse format of the result (e.g. “csr”). If not given, the matrix is returned in “coo” format.
dtype : dtype specifier, optional
The data-type of the output matrix. If not given, the dtype is determined from that of blocks.

res : sparse matrix

New in version 0.11.0.

bmat, diags

>>> from scipy.sparse import coo_matrix, block_diag
>>> A = coo_matrix([[1, 2], [3, 4]])
>>> B = coo_matrix([[5], [6]])
>>> C = coo_matrix([[7]])
>>> block_diag((A, B, C)).toarray()
array([[1, 2, 0, 0],
       [3, 4, 0, 0],
       [0, 0, 5, 0],
       [0, 0, 6, 0],
       [0, 0, 0, 7]])
random(m, n, density=0.01, format="coo", dtype=None, random_state=None, data_rvs=None)

Generate a sparse matrix of the given shape and density with randomly distributed values.

m, n : int
shape of the matrix
density : real, optional
density of the generated matrix: density equal to one means a full matrix, density of 0 means a matrix with no non-zero items.
format : str, optional
sparse matrix format.
dtype : dtype, optional
type of the returned matrix values.
random_state : {numpy.random.RandomState, int}, optional
Random number generator or random seed. If not given, the singleton numpy.random will be used. This random state will be used for sampling the sparsity structure, but not necessarily for sampling the values of the structurally nonzero entries of the matrix.
data_rvs : callable, optional
Samples a requested number of random values. This function should take a single argument specifying the length of the ndarray that it will return. The structurally nonzero entries of the sparse random matrix will be taken from the array sampled by this function. By default, uniform [0, 1) random values will be sampled using the same random state as is used for sampling the sparsity structure.

res : sparse matrix

>>> from scipy.sparse import random
>>> from scipy import stats
>>> class CustomRandomState(object):
...     def randint(self, k):
...         i = np.random.randint(k)
...         return i - i % 2
>>> rs = CustomRandomState()
>>> rvs = stats.poisson(25, loc=10).rvs
>>> S = random(3, 4, density=0.25, random_state=rs, data_rvs=rvs)
>>> S.A
array([[ 36.,   0.,  33.,   0.],   # random
       [  0.,   0.,   0.,   0.],
       [  0.,   0.,  36.,   0.]])

Only float types are supported for now.

rand(m, n, density=0.01, format="coo", dtype=None, random_state=None)

Generate a sparse matrix of the given shape and density with uniformly distributed values.

m, n : int
shape of the matrix
density : real, optional
density of the generated matrix: density equal to one means a full matrix, density of 0 means a matrix with no non-zero items.
format : str, optional
sparse matrix format.
dtype : dtype, optional
type of the returned matrix values.
random_state : {numpy.random.RandomState, int}, optional
Random number generator or random seed. If not given, the singleton numpy.random will be used.

res : sparse matrix

Only float types are supported for now.

scipy.sparse.random : Similar function that allows a user-specified random
data source.
>>> from scipy.sparse import rand
>>> matrix = rand(3, 4, density=0.25, format="csr", random_state=42)
>>> matrix
<3x4 sparse matrix of type '<class 'numpy.float64'>'
   with 3 stored elements in Compressed Sparse Row format>
>>> matrix.todense()
matrix([[ 0.        ,  0.59685016,  0.779691  ,  0.        ],
        [ 0.        ,  0.        ,  0.        ,  0.44583275],
        [ 0.        ,  0.        ,  0.        ,  0.        ]])