sisl.SparseCSR
- class sisl.SparseCSR
Bases:
NDArrayOperatorsMixinA compressed sparse row matrix, slightly different than
csr_matrix.This class holds all required information regarding the CSR matrix format.
Note that this sparse matrix of data does not retain the number of columns in the matrix, i.e. it has no way of determining whether the input is correct.
This sparse matrix class tries to resemble the
csr_matrixas much as possible with the difference of this class being multi-dimensional.Creating a new sparse matrix is much similar to the
scipyequivalent.nnzis only used ifnnz > nr * nnzpr.This class may be instantiated by verious means.
SparseCSR(S)whereSis ascipy.sparsematrixSparseCSR((M,N)[, dtype])the shape of the sparse matrix (equivalent toSparseCSR((M,N,1)[, dtype]).SparseCSR((M,N), dim=K, [, dtype])the shape of the sparse matrix (equivalent toSparseCSR((M,N,K)[, dtype]).SparseCSR((M,N,K)[, dtype])creating a sparse matrix withMrows,Ncolumns andKelements per sparse element.SparseCSR((data, ptr, indices), [shape, dtype])creating a sparse matrix with specific data as would be used when creatingscipy.sparse.csr_matrix.
Additionally these parameters control the creation of the sparse matrix.
- Parameters:
arg1 (
tuple) – various initialization methods as described abovedim (
int, optional) – number of elements stored per sparse element, only used if (M,N) is passeddtype (
numpy.dtype, optional) – data type of the matrix, defaults tonumpy.float64nnzpr (
int, optional) – initial number of non-zero elements per row. Only used ifnnzis not suppliednnz (
int, optional) – initial total number of non-zero elements This quantity has precedence over nnzpr
Plotting
Plotting functions for the
SparseCSRclass.plot.atomicmatrix([dim, isc, ...])Builds a
AtomicMatrixPlotby setting the value of "matrix" to the current object.Methods
align(other)Aligns this sparse matrix with the sparse elements of the other sparse matrix
astype(dtype[, copy])Convert the stored data-type to something else
copy([dims, dtype])A deepcopy of the sparse matrix
delete_columns(cols[, keep_shape])Delete all columns in columns (in-place action)
diagonal()Return the diagonal elements from the matrix
diags(diagonals[, offsets, dim, dtype])Create a
SparseCSRwith diagonal elements with the same shape as the routineedges(rows[, exclude])Retrieve edges (connections) of given rows
eliminate_zeros([atol])Remove all zero elememts from the sparse matrix
empty([keep_nnz])Delete all sparse information from the sparsity pattern
finalize([sort])Finalizes the sparse matrix by removing all non-set elements
fromsp(sparse_matrices[, dtype])Combine multiple single-dimension sparse matrices into one SparseCSR matrix
iter_nnz([rows])Iterations of the non-zero elements, returns a tuple of row and column with non-zero elements
nonzero([rows, only_cols])Row and column indices where non-zero elements exists
remove(indices)Return a new sparse CSR matrix with all the indices removed
scale_columns(cols, scale[, rows])Scale all values with certain column values with a number
sparsity_union(sparse_matrices[, dtype, ...])Create a
SparseCSRwith constant fill value in all places that sparse_matrices have nonzerosspsame(other)Check whether two sparse matrices have the same non-zero elements
sub(indices)Create a new sparse CSR matrix with the data only for the given rows and columns
toarray()Return a dense
numpy.ndarraywhich has 3 dimensions (self.shape)tocsr([dim])Convert dimension
diminto acsr_matrixformattodense()Return a dense
numpy.ndarraywhich has 3 dimensions (self.shape)transform(matrix[, dtype])Apply a linear transformation \(R^n \rightarrow R^m\) to the \(n\)-dimensional elements of the sparse matrix
translate_columns(old, new[, rows, clean])Takes all old columns and translates them to new.
transpose([sort])Create the transposed sparse matrix
Attributes
Data contained in the sparse matrix (numpy array of elements)
The extra dimensionality of the sparse matrix (elements per matrix element)
The data-type in the sparse matrix (in str)
The data-type in the sparse matrix
Whether the contained data is finalized and non-used elements have been removed
Number of non-zero elements in the sparse matrix
The shape of the sparse matrix
- align(other)[source]
Aligns this sparse matrix with the sparse elements of the other sparse matrix
Routine for ensuring that all non-zero elements in other are also in this object.
I.e. this will, possibly, change the sparse elements in-place.
A
ValueErrorwill be raised if the shapes are not mergeable.- Parameters:
other (
SparseCSR) – the other sparse matrix to align.
- astype(dtype, copy=True)[source]
Convert the stored data-type to something else
- Parameters:
dtype – the new dtype for the sparse matrix
copy (bool) – copy when needed, or do not copy when not needed.
- Return type:
Self
- copy(dims=None, dtype=None)[source]
A deepcopy of the sparse matrix
- Parameters:
dims (Optional[SeqOrScalarInt]) – which dimensions to store in the copy, defaults to all.
dtype (
numpy.dtype) – this defaults to the dtype of the object, but one may change it if supplied.
- delete_columns(cols, keep_shape=False)[source]
Delete all columns in columns (in-place action)
- Parameters:
- diags(diagonals, offsets=0, dim=None, dtype=None)[source]
Create a
SparseCSRwith diagonal elements with the same shape as the routine- Parameters:
diagonals (
scalarorndarray) – the diagonal values, if scalar theshapemust be present.offsets (
scalarorndarray) – the offsets from the diagonal for each of the components (defaults to the diagonal)dim (Optional[int]) – the extra dimension of the new diagonal matrix (default to the current extra dimension)
dtype (
numpy.dtype, optional) – the data-type to create (default tonumpy.float64)
- edges(rows, exclude=None)[source]
Retrieve edges (connections) of given rows
The returned edges are unique and sorted (see
numpy.unique).- Parameters:
rows (SeqOrScalarInt) – the edges are returned only for the given row
exclude (Optional[SeqOrScalarInt]) – remove edges which are in the exclude list.
- eliminate_zeros(atol=0.0)[source]
Remove all zero elememts from the sparse matrix
This is an in-place operation
- empty(keep_nnz=False)[source]
Delete all sparse information from the sparsity pattern
Essentially this deletes all entries.
- finalize(sort=True)[source]
Finalizes the sparse matrix by removing all non-set elements
One may still interact with the sparse matrix as one would previously.
NOTE: This is mainly an internal used routine to ensure data structure when converting to
csr_matrix
- classmethod fromsp(sparse_matrices, dtype=None)[source]
Combine multiple single-dimension sparse matrices into one SparseCSR matrix
The different sparse matrices need not have the same sparsity pattern.
- Parameters:
sparse_matrices (OrSequence[SparseMatrix]) – any sparse matrix which can convert to a
scipy.sparse.csr_matrixmatrixdtype (
numpy.dtype, optional) – data-type to store in the matrix, default to largestdtypefor the passed sparse matrices
- iter_nnz(rows=None)[source]
Iterations of the non-zero elements, returns a tuple of row and column with non-zero elements
An iterator returning the current row index and the corresponding column index.
>>> for r, c in self:
In the above case
randcare rows and columns such that>>> self[r, c]
returns the non-zero element of the sparse matrix.
- Parameters:
row (
intorndarrayofint) – only loop on the given row(s) default to all rowsrows (Optional[SeqOrScalarInt])
- nonzero(rows=None, only_cols=False)[source]
Row and column indices where non-zero elements exists
- Parameters:
rows (Optional[SeqOrScalarInt]) – only return the tuples for the requested rows, default is all rows
only_cols (bool) – only return the non-zero columns
- plot.atomicmatrix(dim=0, isc=None, fill_value=None, geometry=None, atom_lines=False, orbital_lines=False, sc_lines=False, color_pixels=True, colorscale='RdBu', crange=None, cmid=None, text=None, textfont={}, set_labels=False, constrain_axes=True, arrows=[], backend='plotly')
Builds a
AtomicMatrixPlotby setting the value of “matrix” to the current object.Plots a (possibly sparse) matrix where rows and columns are either orbitals or atoms.
- Parameters:
dim (int) – If the matrix has a third dimension (e.g. spin), which index to plot in that third dimension.
isc (Optional[int]) – If the matrix contains data for an auxiliary supercell, the index of the cell to plot. If None, the whole matrix is plotted.
fill_value (Optional[float]) – If the matrix is sparse, the value to use for the missing entries.
geometry (Union[Geometry, None]) – Only needed if the matrix does not contain a geometry (e.g. it is a numpy array) and separator lines or labels are requested.
atom_lines (Union[bool, dict]) – If a boolean, whether to draw lines separating atom blocks, using default styles. If a dict, draws the lines with the specified plotly line styles.
orbital_lines (Union[bool, dict]) – If a boolean, whether to draw lines separating blocks of orbital sets, using default styles. If a dict, draws the lines with the specified plotly line styles.
sc_lines (Union[bool, dict]) – If a boolean, whether to draw lines separating the supercells, using default styles. If a dict, draws the lines with the specified plotly line styles.
color_pixels (bool) – Whether to color the pixels of the matrix according to the colorscale.
colorscale (Optional[Colorscale]) – The colorscale to use to color the pixels.
crange (Optional[tuple[float, float]]) – The minimum and maximum values of the colorscale.
cmid (Optional[float]) –
The midpoint of the colorscale. If
crangeis provided, this is ignored.If None and crange is also None, the midpoint is set to 0 if the data contains both positive and negative values.
text (Optional[str]) – If provided, show text of pixel value with the specified format. E.g. text=”.3f” shows the value with three decimal places.
textfont (Optional[dict]) – The font to use for the text. This is a dictionary that may contain the keys “family”, “size”, “color”.
set_labels (bool) – Whether to set the axes labels to the atom/orbital that each row and column corresponds to. For orbitals the labels will be of the form “Atom: (l, m)”, where
Atomis the index of the atom and l and m are the quantum numbers of the orbital.constrain_axes (bool) – Whether to set the ranges of the axes to exactly fit the matrix.
backend (str) – The backend to use for plotting.
- Return type:
- remove(indices)[source]
Return a new sparse CSR matrix with all the indices removed
- Parameters:
indices (SeqOrScalarInt) – the indices of the rows and columns that are removed in the sparse pattern
- Return type:
Self
- scale_columns(cols, scale, rows=None)[source]
Scale all values with certain column values with a number
This will multiply all values with certain column values with
scale\[\mathbf M\[\mathrm{rows}, \mathrm{cols}\] *= \mathrm{scale}\]This is an in-place operation.
- Parameters:
cols (SeqOrScalarInt) – column indices to scale
scale (SeqOrScalarFloat) – scale value for each value (if array-like it has to have the same dimension as the sparsity dimension)
rows (Optional[SeqOrScalarInt]) – only scale the column values that exists in these rows, default to all
- classmethod sparsity_union(sparse_matrices, dtype=None, dim=None, value=0)[source]
Create a
SparseCSRwith constant fill value in all places that sparse_matrices have nonzerosBy default the returned matrix will be sorted.
- Parameters:
sparse_matrices (OrSequence[SparseMatrix]) – SparseCSRs to find the sparsity pattern union of.
dtype (
dtype, optional) – Output dtype. If not given, use the result dtype of the spmats.dim (Optional[int]) – If given, the returned SparseCSR will have this as dim. By default the first sparse matrix in sparse_matrices determines the resulting 3rd dimension.
value (float) – The used fill value.
- Return type:
Self
- sub(indices)[source]
Create a new sparse CSR matrix with the data only for the given rows and columns
All rows and columns in indices are retained, everything else is removed.
- Parameters:
indices (SeqOrScalarInt) – the indices of the rows and columns that are retained in the sparse pattern
- Return type:
Self
- toarray()[source]
Return a dense
numpy.ndarraywhich has 3 dimensions (self.shape)
- tocsr(dim=0, **kwargs)[source]
Convert dimension
diminto acsr_matrixformat- Parameters:
dim (int) – dimension of the data returned in a scipy sparse matrix format
**kwargs – arguments passed to the
csr_matrixroutine
- Return type:
- todense()[source]
Return a dense
numpy.ndarraywhich has 3 dimensions (self.shape)
- transform(matrix, dtype=None)[source]
Apply a linear transformation \(R^n \rightarrow R^m\) to the \(n\)-dimensional elements of the sparse matrix
Notes
The transformation matrix does not act on the rows and columns, only on the final dimension of the matrix.
- Parameters:
matrix (
ndarray) – transformation matrix of shape \(m \times n\), \(n\) should correspond to the number of elements inself.shape[2]dtype (
numpy.dtype, optional) – defaults to the common dtype of the object and the transformation matrix
- translate_columns(old, new, rows=None, clean=True)[source]
Takes all old columns and translates them to new.
- Parameters:
old (SeqOrScalarInt) – old column indices
new (SeqOrScalarInt) – new column indices
rows (Optional[SeqOrScalarInt]) – only translate columns for the given rows
clean (bool) – whether the new translated columns, outside the shape, should be deleted or not (default delete)
- transpose(sort=True)[source]
Create the transposed sparse matrix
- Parameters:
sort (bool) – the returned columns for the transposed structure will be sorted if this is true, default
Notes
The components for each sparse element are not changed in this method.
- col
- property dkind
The data-type in the sparse matrix (in str)
- property dtype
The data-type in the sparse matrix
- property finalized: bool
Whether the contained data is finalized and non-used elements have been removed
- ncol
- ptr
- um = <module 'numpy._core.umath' from '/opt/hostedtoolcache/Python/3.12.12/x64/lib/python3.12/site-packages/numpy/_core/umath.py'>