scipy.sparse.dok_matrix#

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

Dictionary Of Keys based sparse matrix.

This is an efficient structure for constructing sparse matrices incrementally.

This can be instantiated in several ways:
dok_matrix(D)

where D is a 2-D ndarray

dok_matrix(S)

with another sparse array or matrix S (equivalent to S.todok())

dok_matrix((M,N), [dtype])

create the matrix with initial shape (M,N) dtype is optional, defaulting to dtype=’d’

Notes

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

  • Allows for efficient O(1) access of individual elements.

  • Duplicates are not allowed.

  • Can be efficiently converted to a coo_matrix once constructed.

Examples

>>> import numpy as np
>>> from scipy.sparse import dok_matrix
>>> S = dok_matrix((5, 5), dtype=np.float32)
>>> for i in range(5):
...     for j in range(5):
...         S[i, j] = i + j    # Update element
Attributes:
dtypedtype

Data type of the matrix

shape2-tuple

Get shape of a sparse matrix.

ndimint

Number of dimensions (this is always 2)

nnz

Number of stored values, including explicit zeros.

size

Number of stored values.

T

Transpose.

Methods

__getitem__(key)

x.__getitem__(y) <==> x[y]

__len__()

Return len(self).

__mul__(other)

asformat(format[, copy])

Return this array/matrix in the passed format.

asfptype()

Upcast matrix to a floating point format (if necessary)

astype(dtype[, casting, copy])

Cast the array/matrix elements to a specified type.

clear()

conj([copy])

Element-wise complex conjugation.

conjtransp()

Return the conjugate transpose.

conjugate([copy])

Element-wise complex conjugation.

copy()

Returns a copy of this array/matrix.

count_nonzero()

Number of non-zero entries, equivalent to

diagonal([k])

Returns the kth diagonal of the array/matrix.

dot(other)

Ordinary dot product

fromkeys(iterable[, value])

Create a new dictionary with keys from iterable and values set to value.

get(key[, default])

This provides dict.get method functionality with type checking

getH()

Return the Hermitian transpose of this matrix.

get_shape()

Get shape of a sparse matrix.

getcol(j)

Returns a copy of column j of the matrix, as an (m x 1) sparse matrix (column vector).

getformat()

Matrix storage format

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).

items()

keys()

maximum(other)

Element-wise maximum between this and another array/matrix.

mean([axis, dtype, out])

Compute the arithmetic mean along the specified axis.

minimum(other)

Element-wise minimum between this and another array/matrix.

multiply(other)

Point-wise multiplication by another array/matrix.

nonzero()

Nonzero indices of the array/matrix.

pop(k[,d])

If the key is not found, return the default if given; otherwise, raise a KeyError.

popitem()

Remove and return a (key, value) pair as a 2-tuple.

power(n[, dtype])

Element-wise power.

reshape(self, shape[, order, copy])

Gives a new shape to a sparse array/matrix without changing its data.

resize(*shape)

Resize the array/matrix in-place to dimensions given by shape

set_shape(shape)

Set the shape of the matrix in-place

setdefault(key[, default])

Insert key with a value of default if key is not in the dictionary.

setdiag(values[, k])

Set diagonal or off-diagonal elements of the array/matrix.

sum([axis, dtype, out])

Sum the array/matrix elements over a given axis.

toarray([order, out])

Return a dense ndarray representation of this sparse array/matrix.

tobsr([blocksize, copy])

Convert this array/matrix to Block Sparse Row format.

tocoo([copy])

Convert this array/matrix to COOrdinate format.

tocsc([copy])

Convert this array/matrix to Compressed Sparse Column format.

tocsr([copy])

Convert this array/matrix to Compressed Sparse Row format.

todense([order, out])

Return a dense representation of this sparse array/matrix.

todia([copy])

Convert this array/matrix to sparse DIAgonal format.

todok([copy])

Convert this array/matrix to Dictionary Of Keys format.

tolil([copy])

Convert this array/matrix to List of Lists format.

trace([offset])

Returns the sum along diagonals of the sparse array/matrix.

transpose([axes, copy])

Reverses the dimensions of the sparse array/matrix.

update([E, ]**F)

If E is present and has a .keys() method, then does: for k in E: D[k] = E[k] If E is present and lacks a .keys() method, then does: for k, v in E: D[k] = v In either case, this is followed by: for k in F: D[k] = F[k]

values()