sisl.physics.DynamicalMatrix
- class sisl.physics.DynamicalMatrix
Bases:
SparseOrbitalBZDynamical matrix of a geometry
The dynamical matrix is defined as the mass-reduced quantity. Hence the quantities stored in this matrix are expected to contain a factor:
\[\frac1{\sqrt{M_IM_J}}\]for the elements that contains couplings between atoms \(I\) and \(J\).
Plotting
Plotting functions for the SparseOrbital class.
plot.atomicmatrix([dim, ...])Builds a
AtomicMatrixPlotby setting the value of "matrix" to the current object.Methods
Dk([k, dtype, gauge, format])Setup the dynamical matrix for a given k-point
Rij([what, dtype])Create a sparse matrix with the vectors between atoms/orbitals
Sk([k, dtype, gauge, format])Setup the overlap matrix for a given k-point
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
Sometimes the dynamical matrix does not obey Newtons 3rd law.
astype(dtype[, copy])Convert the stored data-type to something else
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.dDk([k, dtype, gauge, format])Setup the dynamical matrix derivative for a given k-point
dSk([k, dtype, gauge, format])Setup the \(\mathbf k\)-derivatie of the overlap matrix for a given k-point
ddDk([k, dtype, gauge, format])Setup the dynamical matrix double derivative for a given k-point
ddSk([k, dtype, gauge, format])Setup the double \(\mathbf k\)-derivatie of the overlap matrix for a given k-point
edges([atoms, exclude, orbitals])Retrieve edges (connections) for all atoms
eig([k, gauge, eigvals_only])Returns the eigenvalues of the physical quantity (using the non-Hermitian solver)
eigenmode([k, gauge])Calculate the eigenmodes at k and return an
EigenmodePhononobject containing all eigenmodeseigenvalue([k, gauge])Calculate the eigenvalues at k and return an
EigenvaluePhononobject containing all eigenvalues for a given keigh([k, gauge, eigvals_only])Returns the eigenvalues of the physical quantity
eigsh([k, n, gauge, eigvals_only])Calculates a subset of eigenvalues of the physical quantity using sparse matrices
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[, S])Create a sparse model from a preset Geometry and a list of sparse matrices
iter_nnz([atoms, orbitals])Iterations of the non-zero elements
iter_orbitals([atoms, local])Iterations of the orbital space in the geometry, two indices from loop
nonzero([atoms, only_cols])Indices row and column indices where non-zero elements exists
prepend(other, axis[, atol, scale])See
appendfor detailsprune_range(*, R[, atoms, atoms_to])Prune elements coupling further than R.
read(sile, *args, **kwargs)Reads dynamical matrix from Sile using read_dynamical_matrix.
remove(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
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 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_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)
write(sile, *args, **kwargs)Writes a dynamical matrix to the Sile as implemented in the
Sile.write_dynamical_matrixmethodAttributes
Access the dynamical matrix elements
Access the overlap elements associated with the sparse matrix
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
True if the object is using a non-orthogonal basis
True if the object is using an orthogonal basis
Shape of sparse matrix
- Dk(k=(0, 0, 0), dtype=None, gauge='lattice', format='csr', *args, **kwargs)[source]
Setup the dynamical matrix for a given k-point
Creation and return of the dynamical matrix for a given k-point (default to Gamma).
Notes
Currently the implemented gauge for the k-point is the cell vector gauge:
\[\mathbf D(\mathbf k) = \mathbf D_{I_\alpha J_\beta} e^{i\mathbf k\cdot\mathbf R}\]where \(\mathbf R\) is an integer times the cell vector and \(\alpha\), \(\beta\) are Cartesian directions.
Another possible gauge is the atomic distance which can be written as
\[\mathbf D(\mathbf k) = \mathbf D_{I_\alpha J_\beta} e^{i\mathbf k\cdot\mathbf r}\]where \(\mathbf r\) is the distance between the atoms.
- Parameters:
k (Sequence[float]) – the k-point to setup the dynamical matrix at
dtype (numpy.dtype , optional) – the data type of the returned matrix. Do NOT request non-complex data-type for non-Gamma k. The default data-type is
numpy.complex128gauge (Literal['lattice', 'atomic']) – the chosen gauge, lattice for lattice vector gauge, and atomic for atomic distance gauge.
format (
{'csr', 'array', 'dense', 'coo', ...}) – the returned format of the matrix, defaulting to thescipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return innumpy.ndarray(‘array’/’dense’/’matrix’). Prefixing with ‘sc:’, or simply ‘sc’ returns the matrix in supercell format with phases.
See also
- Returns:
matrix – the dynamical matrix at \(\mathbf k\). The returned object depends on format.
- Return type:
numpy.ndarrayorscipy.sparse.*_matrix- Parameters:
- Rij(what='orbital', dtype=np.float64)
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.
- Sk(k=(0, 0, 0), dtype=None, gauge='lattice', format='csr', *args, **kwargs)
Setup the overlap matrix for a given k-point
Creation and return of the overlap matrix for a given k-point (default to Gamma).
Notes
Currently the implemented gauge for the k-point is the lattice vector gauge:
\[\mathbf S(\mathbf k) = \mathbf S_{ij} e^{i\mathbf k\cdot\mathbf R}\]where \(\mathbf R\) is an integer times the cell vector and \(i\), \(j\) are orbital indices.
Another possible gauge is the atomic distance which can be written as
\[\mathbf S(\mathbf k) = \mathbf S_{ij} e^{i\mathbf k\cdot\mathbf r}\]where \(\mathbf r\) is the distance between the orbitals.
- Parameters:
k (Sequence[float]) – the k-point to setup the overlap at (default Gamma point)
dtype (
numpy.dtype, optional) – the data type of the returned matrix. Do NOT request non-complex data-type for non-Gamma k. The default data-type isnumpy.complex128gauge (Literal['lattice', 'atomic']) – the chosen gauge
format (
{"csr", "array", "matrix", "coo", ...}) – the returned format of the matrix, defaulting to thescipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return innumpy.ndarray(“array”/”dense”/”matrix”). Prefixing with “sc:”, or simply “sc” returns the matrix in supercell format with phases. This is useful for e.g. bond-current calculations where individual hopping + phases are required.
See also
- Returns:
matrix – the overlap matrix at \(\mathbf k\). The returned object depends on format.
- Return type:
numpy.ndarrayorscipy.sparse.*_matrix- Parameters:
- add(other, axis=None, offset=(0, 0, 0))
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
appendbut 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 (int | None) – whether a specific axis of the cell will be added to the final geometry. For
Nonethe final cell will be that of self, otherwise the lattice vector corresponding to axis will be appended.offset (Sequence[float]) – offset in geometry of other when adding the atoms.
- append(other, axis, atol=0.005, scale=1)
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 (int) – axis to append the two sparse geometries along
atol (float) – 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 (
floatorndarray, 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.
- Return type:
Self
See also
prependequivalent scheme as this method
addmerge two matrices without considering overlap or commensurability
transposeensure hermiticity by using this routine
replacereplace a sub-set of atoms with another sparse matrix
Geometry.append,Geometry.prependSparseCSR.scale_columnsmethod 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:
a new instance with two sparse matrices joined and appended together
- Return type:
- Parameters:
axis (int)
atol (float)
scale (SeqOrScalarFloat)
- apply_newton()[source]
Sometimes the dynamical matrix does not obey Newtons 3rd law.
We correct the dynamical matrix by imposing zero force.
Correcting for Newton forces the matrix to be finalized.
Notes
This is an in-place operation.
- Return type:
- astype(dtype, copy=True)
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
- 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 in Geometry.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)
- dDk(k=(0, 0, 0), dtype=None, gauge='lattice', format='csr', *args, **kwargs)[source]
Setup the dynamical matrix derivative for a given k-point
Creation and return of the dynamical matrix derivative for a given k-point (default to Gamma).
Notes
Currently the implemented gauge for the k-point is the cell vector gauge:
\[\nabla_{\mathbf k} \mathbf D_\gamma(\mathbf k) = i \mathbf R_\gamma \mathbf D_{I_\alpha J_\beta} e^{i\mathbf k\cdot\mathbf R}\]where \(\mathbf R\) is an integer times the cell vector and \(\alpha\), \(\beta\) are atomic indices. And \(\gamma\) is one of the Cartesian directions.
Another possible gauge is the atomic distance which can be written as
\[\nabla_{\mathbf k} \mathbf D_\gamma(\mathbf k) = i \mathbf r_\gamma \mathbf D_{I_\alpha J_\beta} e^{i\mathbf k\cdot\mathbf r}\]where \(\mathbf r\) is the distance between the atoms.
- Parameters:
k (Sequence[float]) – the k-point to setup the dynamical matrix at
dtype (numpy.dtype , optional) – the data type of the returned matrix. Do NOT request non-complex data-type for non-Gamma k. The default data-type is
numpy.complex128gauge (Literal['lattice', 'atomic']) – the chosen gauge, lattice for lattice vector gauge, and atomic for atomic distance gauge.
format (
{'csr', 'array', 'dense', 'coo', ...}) – the returned format of the matrix, defaulting to thescipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return innumpy.ndarray(‘array’/’dense’/’matrix’).
See also
- dSk(k=(0, 0, 0), dtype=None, gauge='lattice', format='csr', *args, **kwargs)
Setup the \(\mathbf k\)-derivatie of the overlap matrix for a given k-point
Creation and return of the derivative of the overlap matrix for a given k-point (default to Gamma).
Notes
Currently the implemented gauge for the k-point is the lattice vector gauge:
\[\nabla_{\mathbf k} \mathbf S_\alpha(\mathbf k) = i \mathbf R_\alpha \mathbf S_{ij} e^{i\mathbf k\cdot\mathbf R}\]where \(\mathbf R\) is an integer times the cell vector and \(i\), \(j\) are orbital indices. And \(\alpha\) is one of the Cartesian directions.
Another possible gauge is the atomic distance which can be written as
\[\nabla_{\mathbf k} \mathbf S_\alpha(\mathbf k) = i \mathbf r_\alpha \mathbf S_{ij} e^{i\mathbf k\cdot\mathbf r}\]where \(\mathbf r\) is the distance between the orbitals.
- Parameters:
k (Sequence[float]) – the k-point to setup the overlap at (default Gamma point)
dtype (
numpy.dtype, optional) – the data type of the returned matrix. Do NOT request non-complex data-type for non-Gamma k. The default data-type isnumpy.complex128gauge (Literal['lattice', 'atomic']) – the chosen gauge.
format (
{"csr", "array", "matrix", "coo", ...}) – the returned format of the matrix, defaulting to thescipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return innumpy.ndarray(“array”/”dense”/”matrix”).
- ddDk(k=(0, 0, 0), dtype=None, gauge='lattice', format='csr', *args, **kwargs)[source]
Setup the dynamical matrix double derivative for a given k-point
Creation and return of the dynamical matrix double derivative for a given k-point (default to Gamma).
Notes
Currently the implemented gauge for the k-point is the cell vector gauge:
\[\nabla_{\mathbf k^2} \mathbf D_{\gamma\sigma}(\mathbf k) = - \mathbf R_\gamma \mathbf R_\sigma \mathbf D_{I_\alpha J_\beta} e^{i\mathbf k\cdot\mathbf R}\]where \(\mathbf R\) is an integer times the cell vector and \(\alpha\), \(\beta\) are Cartesian directions. And \(\gamma\), \(\sigma\) are one of the Cartesian directions.
Another possible gauge is the atomic distance which can be written as
\[\nabla_{\mathbf k^2} \mathbf D_{\gamma\sigma}(\mathbf k) = - \mathbf r_\gamma \mathbf r_\sigma \mathbf D_{I_\alpha J_\beta} e^{i\mathbf k\cdot\mathbf r}\]where \(\mathbf r\) is atomic distance.
- Parameters:
k (Sequence[float]) – the k-point to setup the dynamical matrix at
dtype (numpy.dtype , optional) – the data type of the returned matrix. Do NOT request non-complex data-type for non-Gamma k. The default data-type is
numpy.complex128gauge (Literal['lattice', 'atomic']) – the chosen gauge,
latticefor cell vector gauge, andatomicfor atomic distance gauge.format (
{'csr', 'array', 'dense', 'coo', ...}) – the returned format of the matrix, defaulting to thescipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return innumpy.ndarray(‘array’/’dense’/’matrix’).
- ddSk(k=(0, 0, 0), dtype=None, gauge='lattice', format='csr', *args, **kwargs)
Setup the double \(\mathbf k\)-derivatie of the overlap matrix for a given k-point
Creation and return of the double derivative of the overlap matrix for a given k-point (default to Gamma).
Notes
Currently the implemented gauge for the k-point is the cell vector gauge:
\[\nabla_{\mathbf k^2} \mathbf S_{\alpha\beta}(\mathbf k) = - \mathbf R_\alpha \mathbf R_\beta \mathbf S_{ij} e^{i\mathbf k\cdot\mathbf R}\]where \(\mathbf R\) is an integer times the cell vector and \(i\), \(j\) are orbital indices. And \(\alpha\) and \(\beta\) are one of the Cartesian directions.
Another possible gauge is the atomic distance which can be written as
\[\nabla_{\mathbf k^2} \mathbf S_{\alpha\beta}(\mathbf k) = - \mathbf r_\alpha \mathbf r_\beta \mathbf S_{ij} e^{i\mathbf k\cdot\mathbf r}\]where \(\mathbf r\) is the distance between the orbitals.
- Parameters:
k (Sequence[float]) – the k-point to setup the overlap at (default Gamma point)
dtype (
numpy.dtype, optional) – the data type of the returned matrix. Do NOT request non-complex data-type for non-Gamma k. The default data-type isnumpy.complex128gauge (Literal['lattice', 'atomic']) – the chosen gauge,
cellfor cell vector gauge, andatomfor atomic distance gauge.format (
{"csr", "array", "matrix", "coo", ...}) – the returned format of the matrix, defaulting to thescipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return innumpy.ndarray(“array”/”dense”/”matrix”).
- edges(atoms=None, exclude=None, orbitals=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.no_s).- Parameters:
atoms (AtomsIndex) – 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 (AtomsIndex) – remove edges which are in the exclude list, this list refers to orbitals.
orbitals (
intorlistofint) – the edges are returned only for the given orbital. The returned edges are orbitals.
- Returns:
If orbitals is None, the returned indices are atomic indices. Otherwise it will be orbital indices.
- Return type:
indices
See also
SparseCSR.edgesthe underlying routine used for extracting the edges
- eig(k=(0, 0, 0), gauge='lattice', eigvals_only=True, **kwargs)
Returns the eigenvalues of the physical quantity (using the non-Hermitian solver)
Setup the system and overlap matrix with respect to the given k-point and calculate the eigenvalues.
All subsequent arguments gets passed directly to
scipy.linalg.eig
- eigenmode(k=(0, 0, 0), gauge='lattice', **kwargs)[source]
Calculate the eigenmodes at k and return an
EigenmodePhononobject containing all eigenmodesNotes
Note that the phonon modes are not mass-scaled.
- Parameters:
k (Sequence[float]) – the k-point at which to evaluate the eigenmodes at
gauge (Literal['lattice', 'atomic']) – the gauge used for calculating the eigenmodes
sparse (
bool, optional) – ifTrue,eigshwill be called, elseeighwill be called (default).**kwargs (
dict, optional) – passed arguments to the eigenvalue calculator routine
- Return type:
See also
- eigenvalue(k=(0, 0, 0), gauge='lattice', **kwargs)[source]
Calculate the eigenvalues at k and return an
EigenvaluePhononobject containing all eigenvalues for a given k- Parameters:
k (Sequence[float]) – the k-point at which to evaluate the eigenvalues at
gauge (Literal['lattice', 'atomic']) – the gauge used for calculating the eigenvalues
sparse (
bool, optional) – ifTrue,eigshwill be called, elseeighwill be called (default).**kwargs (
dict, optional) – passed arguments to the eigenvalue calculator routine
- Return type:
- eigh(k=(0, 0, 0), gauge='lattice', eigvals_only=True, **kwargs)
Returns the eigenvalues of the physical quantity
Setup the system and overlap matrix with respect to the given k-point and calculate the eigenvalues.
All subsequent arguments gets passed directly to
scipy.linalg.eigh
- eigsh(k=(0, 0, 0), n=1, gauge='lattice', eigvals_only=True, **kwargs)
Calculates a subset of eigenvalues of the physical quantity using sparse matrices
Setup the quantity and overlap matrix with respect to the given k-point and calculate a subset of the eigenvalues using the sparse algorithms.
All subsequent arguments gets passed directly to
scipy.sparse.linalg.eigsh.- Parameters:
- Return type:
Notes
The performance and accuracy of this method depends heavily on kwargs. Playing around with a small test example before doing large scale calculations is adviced!
- 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, S=None, **kwargs)
Create a sparse model from a preset Geometry and a list of sparse matrices
The passed sparse matrices are in one of
scipy.sparseformats.- Parameters:
geometry (Geometry) – geometry to describe the new sparse geometry
P (Union[OrSequence[SparseMatrix], SparseMatrixPhysical]) – the new sparse matrices that are to be populated in the sparse matrix. If P contains a
sislsparse matrix with an overlap matrix, that part of the matrix will be omitted. UseSfor included the overlap.S (Optional[Union[SparseMatrix, SparseMatrixPhysical]]) – if provided this refers to the overlap matrix and will force the returned sparse matrix to be non-orthogonal. If the passed matrix is a non-orthogonal
sislmatrix object (e.g. aHamiltonian), then it will take the overlap part of the object and pass that along. See examples for details.**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 and optionally being non-orthogonal if
Sis not none- Return type:
SparseGeometry
Examples
Merging two Hamiltonians, for instance a spin-up/down Hamiltonian
>>> H1 = si.Hamiltonian(...) >>> H2 = si.Hamiltonian(...) >>> H = H1.fromsp([H1, H2])
Adding an overlap from another matrix.
H, will now only contain theH1data and the overlap matrix fromH2(the Hamiltonian values inH2will be neglected)>>> H1 = si.Hamiltonian(..., orthogonal=True) >>> H2 = si.Hamiltonian(..., orthogonal=False) >>> H = H1.fromsp(H1, S=H2)
If one wishes to construct a merged Hamiltonian with the overlap parts in the final matrix, then it should be added explicitly.
>>> H1 = si.Hamiltonian(..., orthogonal=False) >>> s = H1.shape >>> assert H1.fromsp([H1, H1]).shape == (s[0], s[1], s[2] * 2 - 2) >>> assert H1.fromsp([H1, H1], S=H1).shape == (s[0], s[1], s[2] * 2 - 1)
- iter_nnz(atoms=None, orbitals=None)
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 orbitals on these atoms (not compatible with the orbitals keyword)
orbitals (
intorndarray) – only loop on the non-zero elements coinciding with the orbital (not compatible with the atoms keyword)
- iter_orbitals(atoms=None, local=False)
Iterations of the orbital space in the geometry, two indices from loop
An iterator returning the current atomic index and the corresponding orbital index.
>>> for ia, io in self.iter_orbitals():
In the above case
ioalways belongs to atom ia and ia may be repeated according to the number of orbitals associated with the atom ia.- Parameters:
- Yields:
ia– atomic indexio– orbital index
See also
Geometry.iter_orbitalsmethod used to iterate orbitals
- nonzero(atoms=None, only_cols=False)
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 But for all orbitals.
only_cols (bool) – only return then 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 Atom is 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:
- prepend(other, axis, atol=0.005, scale=1)
See
appendfor detailsThis is currently equivalent to:
>>> other.append(self, axis, atol, scale)
- Parameters:
axis (int)
atol (float)
scale (SeqOrScalarFloat)
- Return type:
Self
- prune_range(*, R, atoms=None, atoms_to=None)
Prune elements coupling further than R.
Search for connections from atoms to atoms_to where the distances are further than R. If such couplings exists, remove them.
- Parameters:
- Return type:
Self
Notes
The transposed couplings are also deleted to ensure a symmetric matrix.
Currently, one cannot select subset of atoms.
- static read(sile, *args, **kwargs)[source]
Reads dynamical matrix from Sile using read_dynamical_matrix.
- Parameters:
sile (
Sile,strorpathlib.Path) – a Sile object which will be used to read the dynamical matrix. If it is a string it will create a new sile using get_sile.* (
args passed directlytoread_dynamical_matrix(,**))
- 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
- remove_orbital(atoms, orbitals)
Remove a subset of orbitals on atoms according to orbitals
For more detailed examples, please see the equivalent (but opposite) method
sub_orbital.- Parameters:
atoms (AtomsIndex) – indices of atoms or Atom that will be reduced in size according to orbitals
orbitals (
ndarrayofintorOrbital) – indices of the orbitals on atoms that are removed from the sparse matrix.
See also
sub_orbitalretaining a set of orbitals (see here for examples)
- repeat(reps, axis)
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 (int) – number of repetitions along cell-vector axis
axis (int) – 0, 1, 2 according to the cell-direction
SO (SparseOrbital)
- Return type:
See also
Geometry.repeatthe same ordering as the final geometry
Geometry.tilea different ordering of the final geometry
SparseOrbital.tilea different ordering of the final geometry
- replace(atoms, other, other_atoms=None, atol=0.005, scale=1.0)
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 -> otheron atoms belonging toselfandother -> selffromother. 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
bigandminimalone can do:>>> big2 = big.replace(np.arange(big.na), minimal) >>> big2 = (big2 + big2.transpose()) / 2
To retain couplings only from the
bigsparse matrix, one should do the following (note the subsequent transposing which ensures hermiticy and is effectively copying couplings frombigto 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
minimalsparse 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
diffwill contain the difference(self -> other) - (other -> self)>>> diff = sm - sm.transpose()
- Parameters:
atoms (AtomsIndex) – 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 (AtomsIndex) – to select a subset of atoms in other that are taken out. Defaults to all atoms in other.
atol (float) – 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 (SeqOrScalarFloat) – 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.
- Return type:
Self
See also
prependprepend two sparse matrices, see
appendfor detailsaddmerge two matrices without considering overlap or commensurability
transposemay be used to ensure hermiticity (symmetrization of the matrix elements)
appendappend two sparse matrices
Geometry.append,Geometry.prependSparseCSR.scale_columnsmethod 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:
a new instance with two sparse matrices merged together by replacing some atoms
- Return type:
- Parameters:
atoms (AtomsIndex)
other_atoms (AtomsIndex)
atol (float)
scale (SeqOrScalarFloat)
- 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(what='orbital', dtype=np.float64)
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)
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_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 atoms corresponding to atoms
Negative indices are wrapped and thus works, supercell atoms are also wrapped to the unit-cell.
- Parameters:
atoms (AtomsIndex) – indices of retained atoms or Atom for retaining only that atom
SO (SparseOrbital)
- Return type:
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.removethe negative of Geometry.sub
Geometry.subequivalent to the resulting Geometry from this routine
SparseOrbital.removethe negative of
sub, i.e. remove a subset of atoms
- sub_orbital(atoms, orbitals)
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:
atoms (AtomsIndex) – indices of atoms or Atom that will be reduced in size according to orbitals
orbitals (Union[SeqOrScalarInt, Orbital]) – indices of the orbitals on atoms that are retained in the sparse matrix, the list of orbitals will be sorted. One cannot re-arrange matrix elements currently.
Notes
Future implementations may allow one to re-arrange 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
objis 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_orbitalremoving a set of orbitals (opposite of this)
- 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 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 (int) – number of repetitions along cell-vector axis
axis (int) – 0, 1, 2 according to the cell-direction
SO (SparseOrbital)
- Return type:
See also
SparseOrbital.repeata different ordering of the final geometry
SparseOrbital.untileopposite of this method
Geometry.tilethe same ordering as the final geometry
Geometry.repeata different ordering of the final geometry
- toSparseAtom(dim=None, dtype=None)
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 (int | None) – 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=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:
SparseOrbitalorSparseAtom
- 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)
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.untile
sym (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
- write(sile, *args, **kwargs)
Writes a dynamical matrix to the Sile as implemented in the
Sile.write_dynamical_matrixmethod- Parameters:
dyn (DynamicalMatrix)
sile (SileLike)
- Return type:
- property D
Access the dynamical matrix elements
- property S: Self
Access the overlap elements associated with the sparse matrix
- 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 SparseOrbital class.
- um = <module 'numpy._core.umath' from '/opt/hostedtoolcache/Python/3.12.12/x64/lib/python3.12/site-packages/numpy/_core/umath.py'>