sisl.SparseAtom
- class sisl.SparseAtom
Bases:
_SparseGeometrySparse object with number of rows equal to the total number of atoms in the
GeometryPlotting
Plotting functions for the
SparseAtomclass.plot.atomicmatrix([dim, isc, ...])Builds a
AtomicMatrixPlotby setting the value of "matrix" to the current object.Methods
Rij([dtype])Create a sparse matrix with vectors between atoms
construct(func[, na_iR, method, eta])Automatically construct the sparse model based on a function that does the setting up of the elements
copy([dtype])A copy of this object
create_construct(R, params)Create a simple function for passing to the
constructfunction.edges(atoms[, exclude])Retrieve edges (connections) for all atoms
eliminate_zeros(*args, **kwargs)Removes all zero elements from the sparse matrix
empty([keep_nnz])See
emptyfor detailsfinalize(*args, **kwargs)Finalizes the model
fromsp(geometry, P, **kwargs)Create a sparse model from a preset
Geometryand a list of sparse matricesiter_nnz([atoms])Iterations of the non-zero elements
nonzero([atoms, only_cols])Indices row and column indices where non-zero elements exists
remove(atoms)Create a subset of this sparse matrix by removing the atoms corresponding to atoms
repeat(reps, axis)Create a repeated sparse atom object, equivalent to
Geometry.repeatreset([dim, dtype, nnzpr])The sparsity pattern has all elements removed and everything is reset.
rij([dtype])Create a sparse matrix with the distance between atoms
set_nsc(*args, **kwargs)Reset the number of allowed supercells in the sparse atom
spalign(other)See
alignfor detailsspsame(other)Compare two sparse objects and check whether they have the same entries.
sub(atoms)Create a subset of this sparse matrix by only retaining the elements corresponding to the atoms
swap(atoms_a, atoms_b)Swaps atoms in the sparse geometry to obtain a new order of atoms
tile(reps, axis)Create a tiled sparse atom object, equivalent to
Geometry.tiletocsr([dim, isc])Return a
csr_matrixfor the specified dimensiontranslate2uc([atoms, axes])Translates all primary atoms to the unit cell.
transpose([sort])Create the transposed sparse geometry by interchanging supercell indices
unrepeat(reps, axis[, segment, sym])Unrepeats the sparse model into different parts (retaining couplings)
untile(reps, axis[, segment, sym])Untiles the sparse model into different parts (retaining couplings)
Attributes
Number of components per element
Data type of sparse elements (in str)
Data type of sparse elements
Whether the contained data is finalized and non-used elements have been removed
Associated geometry
Number of non-zero elements
Shape of sparse matrix
- Rij(dtype=np.float64)[source]
Create a sparse matrix with vectors between atoms
- Parameters:
dtype (
numpy.dtype, optional) – the data-type of the sparse matrix.
Notes
The returned sparse matrix with vectors are taken from the current sparse pattern. I.e. a subsequent addition of sparse elements will make them inequivalent. It is thus important to only create the sparse vector matrix when the sparse structure is completed.
- construct(func, na_iR=1000, method='rand', eta=None)
Automatically construct the sparse model based on a function that does the setting up of the elements
This may be called in two variants.
Pass a function (func), see e.g.
create_constructwhich does the setting up.Pass a tuple/list in func which consists of two elements, one is
Rthe radii parameters for the corresponding parameters. The second is the parameters corresponding to theR[i]elements. In this second case all atoms must only have one orbital.
- Parameters:
func (
callableorndarray) –this function must take 4 arguments. 1. Is this object (
self) 2. Is the currently examined atom (ia) 3. Is the currently bounded indices (idxs) 4. Is the currently bounded indices atomic coordinates (idxs_xyz) An example func could be:>>> def func(self, ia, atoms, atoms_xyz=None): ... idx = self.geometry.close(ia, R=[0.1, 1.44], atoms=atoms, atoms_xyz=atoms_xyz) ... self[ia, idx[0]] = 0 ... self[ia, idx[1]] = -2.7
na_iR (
int, optional) – number of atoms within the sphere for speeding up the iter_block loop.method (
{'rand', str}) – method used inGeometry.iter_block, see there for detailseta (
bool, optional) – whether an ETA will be printed
See also
create_constructa generic function used to create a generic function which this routine requires
tiletiling after construct is much faster for very large systems
repeatrepeating after construct is much faster for very large systems
- copy(dtype=None)
A copy of this object
- Parameters:
dtype (
numpy.dtype, optional) – it is possible to convert the data to a different data-type If not specified, it will useself.dtypeS (_SparseGeometry)
- Return type:
_SparseGeometry
- create_construct(R, params)
Create a simple function for passing to the
constructfunction.This is simply to leviate the creation of simplistic functions needed for setting up the sparse elements.
Basically this returns a function:
>>> def func(self, ia, atoms, atoms_xyz=None): ... idx = self.geometry.close(ia, R=R, atoms=atoms, atoms_xyz=atoms_xyz) ... for ix, p in zip(idx, params): ... self[ia, ix] = p
Notes
This function only works for geometry sparse matrices (i.e. one element per atom). If you have more than one element per atom you have to implement the function your-self.
- Parameters:
See also
constructroutine to create the sparse matrix from a generic function (as returned from
create_construct)
- edges(atoms, exclude=None)
Retrieve edges (connections) for all atoms
The returned edges are unique and sorted (see
numpy.unique) and are returned in supercell indices (i.e.0 <= edge < self.geometry.na_s).- Parameters:
atoms (AtomsIndex) – the edges are returned only for the given atom
exclude (AtomsIndex) – remove edges which are in the exclude list.
See also
SparseCSR.edgesthe underlying routine used for extracting the edges
- eliminate_zeros(*args, **kwargs)
Removes all zero elements from the sparse matrix
This is an in-place operation.
See also
SparseCSR.eliminate_zerosmethod called, see there for parameters
- Return type:
- finalize(*args, **kwargs)
Finalizes the model
Finalizes the model so that all non-used elements are removed. I.e. this simply reduces the memory requirement for the sparse matrix.
Notes
Adding more elements to the sparse matrix is more time-consuming than for a non-finalized sparse matrix due to the internal data-representation.
- Return type:
- classmethod fromsp(geometry, P, **kwargs)
Create a sparse model from a preset
Geometryand a list of sparse matricesThe passed sparse matrices are in one of
scipy.sparseformats.- Parameters:
geometry (Geometry) – geometry to describe the new sparse geometry
P (OrSequence[SparseMatrix]) – the new sparse matrices that are to be populated in the sparse matrix
**kwargs – any arguments that are directly passed to the
__init__method of the class.
- Returns:
a new sparse matrix that holds the passed geometry and the elements of P
- Return type:
SparseGeometry
- iter_nnz(atoms=None)[source]
Iterations of the non-zero elements
An iterator on the sparse matrix with, row and column
Examples
>>> for i, j in self.iter_nnz(): ... self[i, j] # is then the non-zero value
- Parameters:
atoms (AtomsIndex) – only loop on the non-zero elements coinciding with the atoms
- nonzero(atoms=None, only_cols=False)[source]
Indices row and column indices where non-zero elements exists
- Parameters:
atoms (AtomsIndex) – only return the tuples for the requested atoms, default is all atoms
only_cols (bool) – only return the non-zero columns
See also
SparseCSR.nonzerothe equivalent function call
- 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(atoms)
Create a subset of this sparse matrix by removing the atoms corresponding to atoms
Negative indices are wrapped and thus works.
- Parameters:
atoms (AtomsIndex) – indices of removed atoms
S (_SparseGeometry)
- Return type:
_SparseGeometry
See also
Geometry.removeequivalent to the resulting
Geometryfrom this routineGeometry.subthe negative of
Geometry.removesubthe opposite of
remove, i.e. retain a subset of atoms
- repeat(reps, axis)
Create a repeated sparse atom object, equivalent to
Geometry.repeatThe already existing sparse elements are extrapolated to the new supercell by repeating them in blocks like the coordinates.
- Parameters:
reps (int) – number of repetitions along cell-vector axis
axis (int) – 0, 1, 2 according to the cell-direction
SA (SparseAtom)
- Return type:
See also
Geometry.repeatthe same ordering as the final geometry
Geometry.tilea different ordering of the final geometry
SparseAtom.tilea different ordering of the final geometry
- reset(dim=None, dtype=np.float64, nnzpr=None)
The sparsity pattern has all elements removed and everything is reset.
The object will be the same as if it had been initialized with the same geometry as it were created with.
- rij(dtype=np.float64)[source]
Create a sparse matrix with the distance between atoms
- Parameters:
dtype (
numpy.dtype, optional) – the data-type of the sparse matrix.
Notes
The returned sparse matrix with distances are taken from the current sparse pattern. I.e. a subsequent addition of sparse elements will make them inequivalent. It is thus important to only create the sparse distance when the sparse structure is completed.
- set_nsc(*args, **kwargs)[source]
Reset the number of allowed supercells in the sparse atom
If one reduces the number of supercells any sparse element that references the supercell will be deleted.
See
Lattice.set_nscfor allowed parameters.See also
Lattice.set_nscthe underlying called method
- spsame(other)
Compare two sparse objects and check whether they have the same entries.
This does not necessarily mean that the elements are the same
- sub(atoms)
Create a subset of this sparse matrix by only retaining the elements corresponding to the atoms
Indices passed MUST be unique.
Negative indices are wrapped and thus works.
- Parameters:
atoms (AtomsIndex) – indices of retained atoms
SA (SparseAtom)
- Return type:
See also
Geometry.removethe negative of
Geometry.subGeometry.subequivalent to the resulting
Geometryfrom this routineSparseAtom.removethe negative of
sub, i.e. remove a subset of atoms
- swap(atoms_a, atoms_b)
Swaps atoms in the sparse geometry to obtain a new order of atoms
This can be used to reorder elements of a geometry.
- Parameters:
atoms_a (AtomsIndex) – the first list of atomic coordinates
atoms_b (AtomsIndex) – the second list of atomic coordinates
S (_SparseGeometry)
- Return type:
_SparseGeometry
- tile(reps, axis)
Create a tiled sparse atom object, equivalent to
Geometry.tileThe already existing sparse elements are extrapolated to the new supercell by repeating them in blocks like the coordinates.
Notes
Calling this routine will automatically
finalizetheSparseAtom. This is required to greatly increase performance.- Parameters:
reps (int) – number of repetitions along cell-vector axis
axis (int) – 0, 1, 2 according to the cell-direction
SA (SparseAtom)
- Return type:
See also
SparseAtom.repeata different ordering of the final geometry
SparseAtom.untileopposite of this method
Geometry.tilethe same ordering as the final geometry
Geometry.repeata different ordering of the final geometry
- tocsr(dim=0, isc=None, **kwargs)
Return a
csr_matrixfor the specified dimension
- translate2uc(atoms=None, axes=None)
Translates all primary atoms to the unit cell.
With this, the coordinates of the geometry are translated to the unit cell and the supercell connections in the matrix are updated accordingly.
- Parameters:
atoms (AtomsIndex) – only translate the specified atoms. If not specified, all atoms will be translated.
axes (Optional[CellAxes]) – only translate certain lattice directions, None specifies only the periodic directions
- Returns:
A new sparse matrix with the updated connections and a new associated geometry.
- Return type:
- transpose(sort=True)
Create the transposed sparse geometry by interchanging supercell indices
Sparse geometries are (typically) relying on symmetry in the supercell picture. Thus when one transposes a sparse geometry one should ideally get the same matrix. This is true for the Hamiltonian, density matrix, etc.
This routine transposes all rows and columns such that any interaction between row, r, and column c in a given supercell (i,j,k) will be transposed into row c, column r in the supercell (-i,-j,-k).
- Parameters:
sort (bool) – the returned columns for the transposed structure will be sorted if this is true, default
- Return type:
Self
Notes
The components for each sparse element are not changed in this method.
Examples
Force a sparse geometry to be Hermitian:
>>> sp = SparseOrbital(...) >>> sp = (sp + sp.transpose()) / 2
- unrepeat(reps, axis, segment=0, *args, sym=True, **kwargs)
Unrepeats the sparse model into different parts (retaining couplings)
Please see
untilefor details, the algorithm and arguments are the same however, this is the opposite ofrepeat.
- untile(reps, axis, segment=0, *args, sym=True, **kwargs)[source]
Untiles the sparse model into different parts (retaining couplings)
Recreates a new sparse object with only the cutted atoms in the structure. This will preserve matrix elements in the supercell.
- Parameters:
reps (int) – number of repetitions the tiling function created (opposite meaning as in
untile)axis (int) – which axis to untile along
segment (int) – which segment to retain. Generally each segment should be equivalent, however requesting individiual segments can help uncover inconsistencies in the sparse matrix
*args – arguments passed directly to
Geometry.untilesym (bool) – if True, the algorithm will ensure the returned matrix is symmetrized (i.e. return
(M + M.transpose())/2, else return data as is. False should generally only be used for debugging precision of the matrix elements, or if one wishes to check the warnings.**kwargs – keyword arguments passed directly to
Geometry.untile
- Return type:
Self
Notes
Untiling structures with
nsc == 1along axis are assumed to have periodic boundary conditions.When untiling structures with
nsc == 1along axis it is important to untile as much as possible. This is because otherwise the algorithm cannot determine the correct couplings. Therefore to create a geometry of 3 times a unit-cell, one should untile to the unit-cell, and subsequently tile 3 times.Consider for example a system of 4 atoms, each atom connects to its 2 neighbors. Due to the PBC atom 0 will connect to 1 and 3. Untiling this structure in 2 will group couplings of atoms 0 and 1. As it will only see one coupling to the right it will halve the coupling and use the same coupling to the left, which is clearly wrong.
In the following the latter is the correct way to do it.
>>> SPM.untile(2, 0) != SPM.untile(4, 0).tile(2, 0)
- Raises:
ValueError : – in case the matrix elements are not conseuctive when determining the new supercell structure. This may often happen if untiling a matrix too few times, and then untiling it again.
- Parameters:
- Return type:
Self
See also
tileopposite of this method
Geometry.untilesame as this method, see details about parameters here
- property dkind
Data type of sparse elements (in str)
- property dtype
Data type of sparse elements
- property finalized: bool
Whether the contained data is finalized and non-used elements have been removed
- plot
Plotting functions for the
SparseAtomclass.
- um = <module 'numpy._core.umath' from '/opt/hostedtoolcache/Python/3.12.12/x64/lib/python3.12/site-packages/numpy/_core/umath.py'>