sisl.SparseOrbital
- class sisl.SparseOrbital(geometry: Geometry, dim: int = 1, dtype=None, nnzpr: int | None = None, **kwargs)
Bases:
_SparseGeometry
Sparse object with number of rows equal to the total number of orbitals in the
Geometry
Plotting
Plotting functions for the
SparseOrbital
class.plot.atomicmatrix
([dim, isc, ...])Builds a
AtomicMatrixPlot
by setting the value of "matrix" to the current object.Methods
Rij
([what, dtype])Create a sparse matrix with the vectors between atoms/orbitals
add
(other[, axis, offset])Add two sparse matrices by adding the parameters to one set.
append
(other, axis[, atol, scale])Append other along axis to construct a new connected sparse matrix
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
construct
function.edges
([atoms, exclude, orbitals])Retrieve edges (connections) for all atoms
eliminate_zeros
(*args, **kwargs)Removes all zero elements from the sparse matrix
empty
([keep_nnz])See
empty
for detailsfinalize
()Finalizes the model
fromsp
(geometry, P, **kwargs)Create a sparse model from a preset
Geometry
and a list of sparse matricesiter_nnz
([atoms, orbitals])Iterations of the non-zero elements
nonzero
([atoms, only_cols])Indices row and column indices where non-zero elements exists
prepend
(other, axis[, atol, scale])See
append
for detailsremove
(atoms)Create a subset of this sparse matrix by removing the atoms corresponding to atoms
remove_orbital
(atoms, orbitals)Remove a subset of orbitals on atoms according to orbitals
repeat
(reps, axis)Create a repeated sparse orbital object, equivalent to
Geometry.repeat
replace
(atoms, other[, other_atoms, atol, scale])Replace atoms in self with other_atoms in other and retain couplings between them
reset
([dim, dtype, nnzpr])The sparsity pattern has all elements removed and everything is reset.
rij
([what, dtype])Create a sparse matrix with the distance between atoms/orbitals
set_nsc
(*args, **kwargs)Reset the number of allowed supercells in the sparse orbital
spalign
(other)See
align
for 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 atoms corresponding to atoms
sub_orbital
(atoms, orbitals)Retain only a subset of the orbitals on atoms according to orbitals
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 orbital object, equivalent to
Geometry.tile
toSparseAtom
([dim, dtype])Convert the sparse object (without data) to a new sparse object with equivalent but reduced sparse pattern
tocsr
([dim, isc])Return a
csr_matrix
for 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(what: str = 'orbital', dtype=np.float64)[source]
Create a sparse matrix with the vectors between atoms/orbitals
- Parameters:
what (
{'orbital', 'atom'}
) – which kind of sparse vector matrix to return, either an atomic vector matrix or an orbital vector matrix. The orbital matrix is equivalent to the atomic one with the same vectors repeated for the same atomic orbitals. The default is the same type as the parent class.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.
- add(other, axis: int | None = None, offset: sisl.typing.Coord = (0, 0, 0))[source]
Add two sparse matrices by adding the parameters to one set. The final matrix will have no couplings between self and other
The final sparse matrix will not have any couplings between self and other. Not even if they have commensurate overlapping regions. If you want to create couplings you have to use
append
but that requires the structures are commensurate in the coupling region.- Parameters:
other (
SparseGeometry
) – the other sparse matrix to be added, all atoms will be appendedaxis – whether a specific axis of the cell will be added to the final geometry. For
None
the final cell will be that of self, otherwise the lattice vector corresponding to axis will be appended.offset – offset in geometry of other when adding the atoms.
- append(other, axis: int, atol: float = 0.005, scale: sisl.typing.SeqOrScalarFloat = 1)[source]
Append other along axis to construct a new connected sparse matrix
This method tries to append two sparse geometry objects together by the following these steps:
Create the new extended geometry
Use neighbor cell couplings from self as the couplings to other This may cause problems if the coupling atoms are not exactly equi-positioned. If the coupling coordinates and the coordinates in other differ by more than 0.01 Ang, a warning will be issued. If this difference is above atol the couplings will be removed.
When appending sparse matrices made up of atoms, this method assumes that the orbitals on the overlapping atoms have the same orbitals, as well as the same orbital ordering.
Examples
>>> sporb = SparseOrbital(....) >>> sporb2 = sporb.append(sporb, 0) >>> sporbt = sporb.tile(2, 0) >>> sporb2.spsame(sporbt) True
To retain couplings only from the left sparse matrix, do:
>>> sporb = left.append(right, 0, scale=(2, 0)) >>> sporb = (sporb + sporb.transpose()) / 2
To retain couplings only from the right sparse matrix, do:
>>> sporb = left.append(right, 0, scale=(0, 2.)) >>> sporb = (sporb + sporb.transpose()) / 2
Notes
The current implementation does not preserve the hermiticity of the matrix. If you want to preserve hermiticity of the matrix you have to do the following:
>>> sm = (sm + sm.transpose()) / 2
- Parameters:
other (
object
) – must be an object of the same type as selfaxis – axis to append the two sparse geometries along
atol – tolerance that all coordinates must be within to allow an append. It is important that this value is smaller than half the distance between the two closests atoms such that there is no ambiguity in selecting equivalent atoms. An internal stricter tolerance is used as a baseline, see above.
scale (
float
orarray_like
, optional) – the scale used for the overlapping region. For scalar values it corresponds to passing:(scale, scale)
. For array-like inputscale[0]
refers to the scale of the matrix elements coupling from self, whilescale[1]
is the scale of the matrix elements in other.
See also
prepend
equivalent scheme as this method
add
merge two matrices without considering overlap or commensurability
transpose
ensure hermiticity by using this routine
replace
replace a sub-set of atoms with another sparse matrix
Geometry.append
,Geometry.prepend
SparseCSR.scale_columns
method used to scale the two matrix elements values
- Raises:
ValueError – if the two geometries are not compatible for either coordinate, orbital or supercell errors
- Returns:
object
– a new instance with two sparse matrices joined and appended together
- construct(func, na_iR: int = 1000, method: str = '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_construct
which does the setting up.Pass a tuple/list in func which consists of two elements, one is
R
the 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 (
callable
orarray_like
) – 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_construct
a generic function used to create a generic function which this routine requires
tile
tiling after construct is much faster for very large systems
repeat
repeating after construct is much faster for very large systems
- copy(dtype=None) _SparseGeometry
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.dtype
- create_construct(R, params)
Create a simple function for passing to the
construct
function.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:
R (
array_like
) – radii parameters for different shells. Must have same length as params or one less. If one less it will be extended withR[0]/100
params (
array_like
) – coupling constants corresponding to the R ranges.params[0, :]
are the elements for the all atoms withinR[0]
of each atom.
See also
construct
routine to create the sparse matrix from a generic function (as returned from
create_construct
)
- edges(atoms: sisl.typing.AtomsIndex = None, exclude: sisl.typing.AtomsIndex = None, orbitals=None)[source]
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.no_s
).- Parameters:
atoms – the edges are returned only for the given atom (but by using all orbitals of the requested atom). The returned edges are also atoms.
exclude – remove edges which are in the exclude list, this list refers to orbitals.
orbitals (
int
orlist
ofint
) – the edges are returned only for the given orbital. The returned edges are orbitals.
See also
SparseCSR.edges
the 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_zeros
method called, see there for parameters
- finalize()
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.
Note that adding more elements to the sparse matrix is more time-consuming than for a non-finalized sparse matrix due to the internal data-representation.
- classmethod fromsp(geometry: Geometry, P, **kwargs)
Create a sparse model from a preset
Geometry
and a list of sparse matricesThe passed sparse matrices are in one of
scipy.sparse
formats.- Parameters:
- Returns:
SparseGeometry
– a new sparse matrix that holds the passed geometry and the elements of P
- iter_nnz(atoms: sisl.typing.AtomsIndex = None, orbitals=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 – only loop on the non-zero elements coinciding with the orbitals on these atoms (not compatible with the orbitals keyword)
orbitals (
int
orarray_like
) – only loop on the non-zero elements coinciding with the orbital (not compatible with the atoms keyword)
- nonzero(atoms: sisl.typing.AtomsIndex = None, only_cols: bool = False)[source]
Indices row and column indices where non-zero elements exists
- Parameters:
atoms – only return the tuples for the requested atoms, default is all atoms But for all orbitals.
only_cols – only return then non-zero columns
See also
SparseCSR.nonzero
the equivalent function call
- plot.atomicmatrix(dim: int = 0, isc: int | None = None, fill_value: float | None = None, geometry: Geometry | None = None, atom_lines: bool | dict = False, orbital_lines: bool | dict = False, sc_lines: bool | dict = False, color_pixels: bool = True, colorscale: Colorscale | None = 'RdBu', crange: tuple[float, float] | None = None, cmid: float | None = None, text: str | None = None, textfont: dict | None = {}, set_labels: bool = False, constrain_axes: bool = True, arrows: list[dict] = [], backend: str = 'plotly') AtomicMatrixPlot
Builds a
AtomicMatrixPlot
by 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 – If the matrix has a third dimension (e.g. spin), which index to plot in that third dimension.
isc – 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 – If the matrix is sparse, the value to use for the missing entries.
geometry – 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 – 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 – 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 – 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 – Whether to color the pixels of the matrix according to the colorscale.
colorscale – The colorscale to use to color the pixels.
crange – The minimum and maximum values of the colorscale.
cmid – The midpoint of the colorscale. If
crange
is 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 – If provided, show text of pixel value with the specified format. E.g. text=”.3f” shows the value with three decimal places.
textfont – The font to use for the text. This is a dictionary that may contain the keys “family”, “size”, “color”.
set_labels – 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
Atom
is the index of the atom and l and m are the quantum numbers of the orbital.constrain_axes – Whether to set the ranges of the axes to exactly fit the matrix.
backend – The backend to use for plotting.
- prepend(other, axis: int, atol: float = 0.005, scale: sisl.typing.SeqOrScalarFloat = 1)[source]
See
append
for detailsThis is currently equivalent to:
>>> other.append(self, axis, atol, scale)
- remove(atoms: sisl.typing.AtomsIndex) _SparseGeometry
Create a subset of this sparse matrix by removing the atoms corresponding to atoms
Negative indices are wrapped and thus works.
- Parameters:
atoms – indices of removed atoms
See also
Geometry.remove
equivalent to the resulting
Geometry
from this routineGeometry.sub
the negative of
Geometry.remove
sub
the opposite of
remove
, i.e. retain a subset of atoms
- remove_orbital(atoms: sisl.typing.AtomsIndex, orbitals)[source]
Remove a subset of orbitals on atoms according to orbitals
For more detailed examples, please see the equivalent (but opposite) method
sub_orbital
.- Parameters:
See also
sub_orbital
retaining a set of orbitals (see here for examples)
- repeat(reps: int, axis: int) SparseOrbital
Create a repeated sparse orbital object, equivalent to
Geometry.repeat
The already existing sparse elements are extrapolated to the new supercell by repeating them in blocks like the coordinates.
- Parameters:
reps – number of repetitions along cell-vector axis
axis – 0, 1, 2 according to the cell-direction
See also
Geometry.repeat
the same ordering as the final geometry
Geometry.tile
a different ordering of the final geometry
SparseOrbital.tile
a different ordering of the final geometry
- replace(atoms: sisl.typing.AtomsIndex, other, other_atoms: sisl.typing.AtomsIndex = None, atol: float = 0.005, scale: sisl.typing.SeqOrScalarFloat = 1.0)[source]
Replace atoms in self with other_atoms in other and retain couplings between them
This method replaces a subset of atoms in self with another sparse geometry retaining any couplings between them. The algorithm checks whether the coupling atoms have the same number of orbitals. Meaning that atoms in the overlapping region should have the same connections and number of orbitals per atom. It will _not_ check whether the orbitals or atoms _are_ the same, nor the order of the orbitals.
The replacement algorithm takes the couplings from
self -> other
on atoms belonging toself
andother -> self
fromother
. This will in some cases mean that the matrix becomes non-symmetric. See in Notes for details on symmetrizing the matrices.Examples
>>> minimal = SparseOrbital(....) >>> big = minimal.tile(2, 0) >>> big2 = big.replace(np.arange(big.na), minimal) >>> big.spsame(big2) True
To ensure hermiticity and using the average of the couplings from
big
andminimal
one can do:>>> big2 = big.replace(np.arange(big.na), minimal) >>> big2 = (big2 + big2.transpose()) / 2
To retain couplings only from the
big
sparse matrix, one should do the following (note the subsequent transposing which ensures hermiticy and is effectively copying couplings frombig
to the replaced region.>>> big2 = big.replace(np.arange(big.na), minimal, scale=(2, 0)) >>> big2 = (big2 + big2.transpose()) / 2
To only retain couplings from the
minimal
sparse matrix:>>> big2 = big.replace(np.arange(big.na), minimal, scale=(0, 2)) >>> big2 = (big2 + big2.transpose()) / 2
Notes
The current implementation does not preserve the hermiticity of the matrix. If you want to preserve hermiticity of the matrix you have to do the following:
>>> sm = (sm + sm.transpose()) / 2
Also note that the ordering of the atoms will be
range(atoms.min()), range(len(other_atoms)), <rest>
.Algorithms that utilizes atomic indices should be careful.
When the tolerance atol is high, the elements may be more prone to differences in the symmetry elements. A good idea would be to check the difference between the couplings. The below variable
diff
will contain the difference(self -> other) - (other -> self)
>>> diff = sm - sm.transpose()
- Parameters:
atoms – which atoms in self that are removed and replaced with
other.sub(other_atoms)
other (
object
) – must be an object of the same type as self, a subset is taken from this sparse matrix and combined with self to create a new sparse matrixother_atoms – to select a subset of atoms in other that are taken out. Defaults to all atoms in other.
atol – coordinate tolerance for allowing replacement. It is important that this value is at least smaller than half the distance between the two closests atoms such that there is no ambiguity in selecting equivalent atoms.
scale – the scale used for the overlapping region. For scalar values it corresponds to passing:
(scale, scale)
. For array-like inputscale[0]
refers to the scale of the matrix elements coupling from self, whilescale[1]
is the scale of the matrix elements in other.
See also
prepend
prepend two sparse matrices, see
append
for detailsadd
merge two matrices without considering overlap or commensurability
transpose
may be used to ensure hermiticity (symmetrization of the matrix elements)
append
append two sparse matrices
Geometry.append
,Geometry.prepend
SparseCSR.scale_columns
method used to scale the two matrix elements values
- Raises:
ValueError – if the two geometries are not compatible for either coordinate, orbital or supercell errors
AssertionError – if the two geometries are not compatible for either coordinate, orbital or supercell errors
- Warns:
SislWarning – in case the overlapping atoms are not comprising the same atomic specie. In some cases this may not be a problem. However, care must be taken by the user if this warning is issued.
- Returns:
object
– a new instance with two sparse matrices merged together by replacing some atoms
- reset(dim: int | None = None, dtype=np.float64, nnzpr: int | None = 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.
- Parameters:
dim – number of dimensions per element, default to the current number of elements per matrix element.
dtype (
numpy.dtype
, optional) – the datatype of the sparse elementsnnzpr – number of non-zero elements per row
- rij(what: str = 'orbital', dtype=np.float64)[source]
Create a sparse matrix with the distance between atoms/orbitals
- Parameters:
what (
{'orbital', 'atom'}
) – which kind of sparse distance matrix to return, either an atomic distance matrix or an orbital distance matrix. The orbital matrix is equivalent to the atomic one with the same distance repeated for the same atomic orbitals. The default is the same type as the parent class.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 orbital
If one reduces the number of supercells any sparse element that references the supercell will be deleted.
See
Lattice.set_nsc
for allowed parameters.See also
Lattice.set_nsc
the underlying called method
- spalign(other)
See
align
for details
- 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: sisl.typing.AtomsIndex) SparseOrbital
Create a subset of this sparse matrix by only retaining the atoms corresponding to atoms
Negative indices are wrapped and thus works, supercell atoms are also wrapped to the unit-cell.
- Parameters:
atoms – indices of retained atoms or
Atom
for retaining only that atom
Examples
>>> obj = SparseOrbital(...) >>> obj.sub(1) # only retain the second atom in the SparseGeometry >>> obj.sub(obj.atoms.atom[0]) # retain all atoms which is equivalent to >>> # the first atomic specie
See also
Geometry.remove
the negative of
Geometry.sub
Geometry.sub
equivalent to the resulting
Geometry
from this routineSparseOrbital.remove
the negative of
sub
, i.e. remove a subset of atoms
- sub_orbital(atoms: sisl.typing.AtomsIndex, orbitals)[source]
Retain only a subset of the orbitals on atoms according to orbitals
This allows one to retain only a given subset of the sparse matrix elements.
- Parameters:
Notes
Future implementations may allow one to re-arange orbitals using this method.
When using this method the internal species list will be populated by another species that is named after the orbitals removed. This is to distinguish different atoms.
Examples
>>> # a Carbon atom with 2 orbitals >>> C = sisl.Atom('C', [1., 2.]) >>> # an oxygen atom with 3 orbitals >>> O = sisl.Atom('O', [1., 2., 2.4]) >>> geometry = sisl.Geometry([[0] * 3, [1] * 3]], 2, [C, O]) >>> obj = SparseOrbital(geometry).tile(3, 0) >>> # fill in obj data...
Now
obj
is a sparse geometry with 2 different species and 6 atoms (3 of each). They are ordered[C, O, C, O, C, O]
. In the following we will note species that are different from the original by a'
in the list.Retain 2nd orbital on the 2nd atom:
[C, O', C, O, C, O]
>>> new_obj = obj.sub_orbital(1, 1)
Retain 2nd orbital on 1st and 2nd atom:
[C', O', C, O, C, O]
>>> new_obj = obj.sub_orbital([0, 1], 1)
Retain 2nd orbital on the 1st atom and 3rd orbital on 4th atom:
[C', O, C, O', C, O]
>>> new_obj = obj.sub_orbital(0, 1).sub_orbital(3, 2)
Retain 2nd orbital on all atoms equivalent to the first atom:
[C', O, C', O, C', O]
>>> new_obj = obj.sub_orbital(obj.geometry.atoms[0], 1)
Retain 1st orbital on 1st atom, and 2nd orbital on 3rd and 5th atom:
[C', O, C'', O, C'', O]
>>> new_obj = obj.sub_orbital(0, 0).sub_orbital([2, 4], 1)
See also
remove_orbital
removing a set of orbitals (opposite of this)
- swap(atoms_a: sisl.typing.AtomsIndex, atoms_b: sisl.typing.AtomsIndex) _SparseGeometry
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 – the first list of atomic coordinates
atoms_b – the second list of atomic coordinates
- tile(reps: int, axis: int) SparseOrbital
Create a tiled sparse orbital object, equivalent to
Geometry.tile
The already existing sparse elements are extrapolated to the new supercell by repeating them in blocks like the coordinates.
- Parameters:
reps – number of repetitions along cell-vector axis
axis – 0, 1, 2 according to the cell-direction
See also
SparseOrbital.repeat
a different ordering of the final geometry
SparseOrbital.untile
opposite of this method
Geometry.tile
the same ordering as the final geometry
Geometry.repeat
a different ordering of the final geometry
- toSparseAtom(dim: int = None, dtype=None)[source]
Convert the sparse object (without data) to a new sparse object with equivalent but reduced sparse pattern
This converts the orbital sparse pattern to an atomic sparse pattern.
- Parameters:
dim – number of dimensions allocated in the SparseAtom object, default to the same
dtype (
numpy.dtype
, optional) – used data-type for the sparse object. Defaults to the same.
- tocsr(dim: int = 0, isc=None, **kwargs)
Return a
csr_matrix
for the specified dimension- Parameters:
dim – the dimension in the sparse matrix (for non-orthogonal cases the last dimension is the overlap matrix)
isc (
int
, optional) – the supercell index, or all (ifisc=None
)
- translate2uc(atoms: sisl.typing.AtomsIndex = None, axes: sisl.typing.CellAxes | None = 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 – only translate the specified atoms. If not specified, all atoms will be translated.
axes – only translate certain lattice directions, None specifies only the periodic directions
- Returns:
SparseOrbital
orSparseAtom
– A new sparse matrix with the updated connections and a new associated geometry.
- transpose(sort: bool = 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 – 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.
Examples
Force a sparse geometry to be Hermitian:
>>> sp = SparseOrbital(...) >>> sp = (sp + sp.transpose()) / 2
- Returns:
object
– an equivalent sparse geometry with transposed matrix elements
- unrepeat(reps: int, axis: int, segment: int = 0, *args, sym: bool = True, **kwargs)
Unrepeats the sparse model into different parts (retaining couplings)
Please see
untile
for details, the algorithm and arguments are the same however, this is the opposite ofrepeat
.
- untile(reps: int, axis: int, segment: int = 0, *args, sym: bool = 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 – number of repetitions the tiling function created (opposite meaning as in
untile
)axis – which axis to untile along
segment – 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.untile
sym – 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
Notes
Untiling structures with
nsc == 1
along axis are assumed to have periodic boundary conditions.When untiling structures with
nsc == 1
along 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.
See also
tile
opposite of this method
Geometry.untile
same as this method, see details about parameters here
- property dim
Number of components per element
- property dkind
Data type of sparse elements (in str)
- property dtype
Data type of sparse elements
- property finalized
Whether the contained data is finalized and non-used elements have been removed
- property geometry
Associated geometry
- property nnz
Number of non-zero elements
- plot
Plotting functions for the
SparseOrbital
class.
- property shape
Shape of sparse matrix