DynamicalMatrix

class sisl.physics.DynamicalMatrix(geometry, dim=1, dtype=None, nnzpr=None, **kwargs)[source]

Dynamical matrix of a geometry

Attributes

D
S
dim Number of components per element
dkind Data type of sparse elements (in str)
dtype Data type of sparse elements
finalized Whether the contained data is finalized and non-used elements have been removed
geom Associated geometry
geometry Associated geometry
nnz Number of non-zero elements
non_orthogonal True if the object is using a non-orthogonal basis
orthogonal True if the object is using an orthogonal basis
shape Shape of sparse matrix

Methods

DOS(E[, k, distribution]) Calculate the DOS at the given energies for a specific k point
Dk([k, dtype, gauge, format]) Setup the dynamical matrix for a given k-point
PDOS(E[, k, distribution]) Calculate the projected DOS at the given energies for a specific 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
__init__(geometry[, dim, dtype, nnzpr]) Create sparse object with element between orbitals
apply_newton() Sometimes the dynamical matrix does not obey Newtons 3rd law.
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, param) Create a simple function for passing to the construct function.
cut(seps, axis, *args, **kwargs) Cuts the sparse orbital model into different parts.
dDk([k, dtype, gauge, format]) Setup the dynamical matrix derivative for a given k-point
dSk([k, dtype, gauge, format]) Setup the \(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 \(k\)-derivatie of the overlap matrix for a given k-point
displacement([k]) Calculate the displacement for the eigenmodes for a given k point
edges([atom, exclude, orbital]) Retrieve edges (connections) of a given atom or list of atom’s
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 EigenmodePhonon object containing all eigenmodes
eigenvalue([k, gauge]) Calculate the eigenvalues at k and return an EigenvaluePhonon object containing all eigenvalues for a given k
eigh([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 (default 10)
eliminate_zeros([atol]) Removes all zero elements from the sparse matrix
empty([keep_nnz]) See empty for details
finalize() Finalizes the model
fromsp(geometry, P[, S]) Read and return the object with possible overlap
iter([local]) Iterations of the orbital space in the geometry, two indices from loop
iter_nnz([atom, orbital]) Iterations of the non-zero elements
make_hermitian() Ensures the matrix is Hermitian by doing an in-place symmetrization
nonzero([atom, only_col]) Indices row and column indices where non-zero elements exists
read(sile, *args, **kwargs) Reads dynamical matrix from Sile using read_dynamical_matrix.
remove(atom[, orb_index]) Remove a subset of this sparse matrix by only retaining the atoms corresponding to atom and optionally a subset of the atom orbitals
repeat(reps, axis) Create a repeated sparse orbital object, equivalent to Geometry.repeat
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 details
spsame(other) Compare two sparse objects and check whether they have the same entries.
sub(atom[, orb_index]) Create a subset of this sparse matrix by only retaining the atoms corresponding to atom
swap(a, 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([dtype]) Convert the sparse object (without data) to a new sparse object with equivalent but reduced sparse pattern
tocsr(index[, isc]) Return a scipy.sparse.csr_matrix of the specified index
velocity([k]) Calculate the velocity for the eigenmodes for a given k point
write(sile, *args, **kwargs) Writes a dynamical matrix to the Sile as implemented in the Sile.write_dynamical_matrix method
D
DOS(E, k=(0, 0, 0), distribution='gaussian', **kwargs)[source]

Calculate the DOS at the given energies for a specific k point

Parameters:
E : array_like

energies to calculate the DOS at

k : array_like, optional

k-point at which the DOS is calculated

distribution : func or str, optional

a function that accepts \(E-\hbar\omega\) as argument and calculates the distribution function.

**kwargs: optional

additional parameters passed to the eigenvalue routine

See also

sisl.physics.distribution
setup a distribution function, see details regarding the distribution argument
eigenvalue
method used to calculate the eigenvalues
PDOS
Calculate projected DOS
EigenvaluePhonon.DOS
Underlying method used to calculate the DOS
Dk(k=(0, 0, 0), dtype=None, gauge='R', 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).

Parameters:
k : array_like

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

gauge : {‘R’, ‘r’}

the chosen gauge, R for cell vector gauge, and r for atomic distance gauge.

format : {‘csr’, ‘array’, ‘dense’, ‘coo’, …}

the returned format of the matrix, defaulting to the scipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return in numpy.ndarray (‘array’) or numpy.matrix (‘dense’).

Returns:
object : the dynamical matrix at \(k\). The returned object depends on format.

See also

dDk
dynamical matrix derivative with respect to k
ddDk
dynamical matrix double derivative with respect to k

Notes

Currently the implemented gauge for the k-point is the cell vector gauge:

\[\mathbf D(k) = \mathbf D_{i_\alpha j_\beta} e^{i q R}\]

where \(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(k) = \mathbf D_{\nu\mu} e^{i k r}\]

where \(r\) is the distance between the atoms.

PDOS(E, k=(0, 0, 0), distribution='gaussian', **kwargs)[source]

Calculate the projected DOS at the given energies for a specific k point

Parameters:
E : array_like

energies to calculate the projected DOS at

k : array_like, optional

k-point at which the projected DOS is calculated

distribution : func or str, optional

a function that accepts \(E-\epsilon\) as argument and calculates the distribution function.

**kwargs: optional

additional parameters passed to the eigenmode routine

See also

sisl.physics.distribution
setup a distribution function, see details regarding the distribution argument
eigenmode
method used to calculate the eigenmodes
DOS
Calculate total DOS
EigenmodePhonon.PDOS
Underlying method used to calculate the projected DOS
Rij(what='orbital', dtype=<class 'numpy.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.

S
Sk(k=(0, 0, 0), dtype=None, gauge='R', 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).

Parameters:
k : array_like, optional

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 is numpy.complex128

gauge : {‘R’, ‘r’}

the chosen gauge, R for cell vector gauge, and r for orbital distance gauge.

format : {‘csr’, ‘array’, ‘matrix’, ‘coo’, …}

the returned format of the matrix, defaulting to the scipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return in numpy.ndarray (‘array’) or numpy.matrix (‘matrix’).

Returns:
object : the overlap matrix for the \(k\)-point, format determines the object type.

See also

dSk
Overlap matrix derivative with respect to k
ddSk
Overlap matrix double derivative with respect to k

Notes

Currently the implemented gauge for the k-point is the cell vector gauge:

\[\mathbf S(k) = \mathbf S_{\nu\mu} e^{i k R}\]

where \(R\) is an integer times the cell vector and \(\nu\), \(\mu\) are orbital indices.

Another possible gauge is the orbital distance which can be written as

\[\mathbf S(k) = \mathbf S_{\nu\mu} e^{i k r}\]

where \(r\) is the distance between the orbitals.

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.

construct(func, na_iR=1000, method='rand', eta=False)

Automatically construct the sparse model based on a function that does the setting up of the elements

This may be called in two variants.

  1. Pass a function (func), see e.g. create_construct which does the setting up.
  2. 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 the R[i] elements. In this second case all atoms must only have one orbital.
Parameters:
func: callable or array_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, idxs, idxs_xyz=None): 
...     idx = self.geometry.close(ia, R=[0.1, 1.44], idx=idxs, idx_xyz=idxs_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 details

eta: 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)

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 use self.dtype

create_construct(R, param)

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, idxs, idxs_xyz=None): 
...     idx = self.geometry.close(ia, R=R, idx=idxs) 
...     for ix, p in zip(idx, param): 
...         self[ia, ix] = p 
Parameters:
R : array_like

radii parameters for different shells. Must have same length as param or one less. If one less it will be extended with R[0]/100

param : array_like

coupling constants corresponding to the R ranges. param[0,:] are the elements for the all atoms within R[0] of each atom.

See also

construct
routine to create the sparse matrix from a generic function (as returned from create_construct)

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.

cut(seps, axis, *args, **kwargs)

Cuts the sparse orbital model into different parts.

Recreates a new sparse orbital object with only the cutted atoms in the structure.

Cutting is the opposite of tiling.

Parameters:
seps : int

number of times the structure will be cut

axis : int

the axis that will be cut

dDk(k=(0, 0, 0), dtype=None, gauge='R', 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).

Parameters:
k : array_like

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

gauge : {‘R’, ‘r’}

the chosen gauge, R for cell vector gauge, and r for atomic distance gauge.

format : {‘csr’, ‘array’, ‘dense’, ‘coo’, …}

the returned format of the matrix, defaulting to the scipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return in numpy.ndarray (‘array’) or numpy.matrix (‘dense’).

Returns:
tuple : for each of the Cartesian directions a \(\partial \mathbf D(k)/\partial k_\gamma\) is returned.

See also

Dk
dynamical matrix with respect to k
ddDk
dynamical matrix double derivative with respect to k

Notes

Currently the implemented gauge for the k-point is the cell vector gauge:

\[\nabla_k \mathbf D_\gamma(k) = i R_\gamma \mathbf D_{i_\alpha j_\beta} e^{i q R}\]

where \(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_k \mathbf D_\gamma(k) = i r_\gamma \mathbf D_{i_\alpha j_\beta} e^{i k r}\]

where \(r\) is the distance between the atoms.

dSk(k=(0, 0, 0), dtype=None, gauge='R', format='csr', *args, **kwargs)

Setup the \(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).

Parameters:
k : array_like, optional

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 is numpy.complex128

gauge : {‘R’, ‘r’}

the chosen gauge, R for cell vector gauge, and r for orbital distance gauge.

format : {‘csr’, ‘array’, ‘matrix’, ‘coo’, …}

the returned format of the matrix, defaulting to the scipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return in numpy.ndarray (‘array’) or numpy.matrix (‘matrix’).

Returns:
tuple : for each of the Cartesian directions a \(\partial \mathbf S(k)/\partial k\) is returned.

See also

Sk
Overlap matrix at k
ddSk
Overlap matrix double derivative at k

Notes

Currently the implemented gauge for the k-point is the cell vector gauge:

\[\nabla_k \mathbf S_\alpha(k) = i R_\alpha \mathbf S_{\nu\mu} e^{i k R}\]

where \(R\) is an integer times the cell vector and \(\nu\), \(\mu\) are orbital indices. And \(\alpha\) is one of the Cartesian directions.

Another possible gauge is the orbital distance which can be written as

\[\nabla_k \mathbf S_\alpha(k) = i r_\alpha \mathbf S_{ij} e^{i k r}\]

where \(r\) is the distance between the orbitals.

ddDk(k=(0, 0, 0), dtype=None, gauge='R', 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).

Parameters:
k : array_like

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

gauge : {‘R’, ‘r’}

the chosen gauge, R for cell vector gauge, and r for orbital distance gauge.

format : {‘csr’, ‘array’, ‘dense’, ‘coo’, …}

the returned format of the matrix, defaulting to the scipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return in numpy.ndarray (‘array’) or numpy.matrix (‘dense’).

Returns:
tuple of tuples : for each of the Cartesian directions

See also

Dk
dynamical matrix with respect to k
dDk
dynamical matrix derivative with respect to k

Notes

Currently the implemented gauge for the k-point is the cell vector gauge:

\[\nabla_k^2 \mathbf D_{\gamma\sigma}(k) = - R_\gamma R_\sigma \mathbf D_{i_\alpha j_\beta} e^{i q R}\]

where \(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_k^2 \mathbf D_{\gamma\sigma}(k) = - r_\gamma r_\sigma \mathbf D_{i_\alpha j_\beta} e^{i k r}\]

where \(r\) is atomic distance.

ddSk(k=(0, 0, 0), dtype=None, gauge='R', format='csr', *args, **kwargs)

Setup the double \(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).

Parameters:
k : array_like, optional

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 is numpy.complex128

gauge : {‘R’, ‘r’}

the chosen gauge, R for cell vector gauge, and r for orbital distance gauge.

format : {‘csr’, ‘array’, ‘matrix’, ‘coo’, …}

the returned format of the matrix, defaulting to the scipy.sparse.csr_matrix, however if one always requires operations on dense matrices, one can always return in numpy.ndarray (‘array’) or numpy.matrix (‘matrix’).

Returns:
tuple of tuples : for each of the Cartesian directions

See also

Sk
Overlap matrix at k
dSk
Overlap matrix derivative at k

Notes

Currently the implemented gauge for the k-point is the cell vector gauge:

\[\nabla_k^2 \mathbf S_{\alpha\beta}(k) = - R_\alpha R_\beta \mathbf S_{\nu\mu} e^{i k R}\]

where \(R\) is an integer times the cell vector and \(\nu\), \(\mu\) are orbital indices. And \(\alpha\) and \(\beta\) are one of the Cartesian directions.

Another possible gauge is the orbital distance which can be written as

\[\nabla_k^2 \mathbf S_{\alpha\beta}(k) = - r_\alpha r_\beta \mathbf S_{ij} e^{i k r}\]

where \(r\) is the distance between the orbitals.

dim

Number of components per element

displacement(k=(0, 0, 0), **kwargs)[source]

Calculate the displacement for the eigenmodes for a given k point

Parameters:
k : array_like, optional

k-point at which the displacement are calculated

**kwargs: optional

additional parameters passed to the eigenmode routine

See also

eigenmode
method used to calculate the eigenmodes
velocity
Calculate mode velocity
EigenmodePhonon.displacement
Underlying method used to calculate the velocity
dkind

Data type of sparse elements (in str)

dtype

Data type of sparse elements

edges(atom=None, exclude=None, orbital=None)

Retrieve edges (connections) of a given atom or list of atom’s

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:
atom : int or list of int

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 : int or list of int, optional

remove edges which are in the exclude list. Default to atom.

orbital : int or list of int

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
eig(k=(0, 0, 0), gauge='R', 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='R', **kwargs)[source]

Calculate the eigenmodes at k and return an EigenmodePhonon object containing all eigenmodes

Parameters:
k : array_like*3, optional

the k-point at which to evaluate the eigenmodes at

gauge : str, optional

the gauge used for calculating the eigenmodes

sparse : bool, optional

if True, eigsh will be called, else eigh will be called (default).

**kwargs : dict, optional

passed arguments to the eigh routine

Returns:
EigenmodePhonon

See also

eigh
eigenvalue routine
eigsh
eigenvalue routine
eigenvalue(k=(0, 0, 0), gauge='R', **kwargs)[source]

Calculate the eigenvalues at k and return an EigenvaluePhonon object containing all eigenvalues for a given k

Parameters:
k : array_like*3, optional

the k-point at which to evaluate the eigenvalues at

gauge : str, optional

the gauge used for calculating the eigenvalues

sparse : bool, optional

if True, eigsh will be called, else eigh will be called (default).

**kwargs : dict, optional

passed arguments to the eigh routine

Returns:
EigenvaluePhonon

See also

eigh
eigenvalue routine
eigsh
eigenvalue routine
eigh(k=(0, 0, 0), gauge='R', 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=10, gauge='R', eigvals_only=True, **kwargs)

Calculates a subset of eigenvalues of the physical quantity (default 10)

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.linalg.eigsh

eliminate_zeros(atol=0.0)

Removes all zero elements from the sparse matrix

This is an in-place operation.

Parameters:
atol : float, optional

absolute tolerance below this value will be considered 0.

empty(keep_nnz=False)

See empty for details

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.

finalized

Whether the contained data is finalized and non-used elements have been removed

classmethod fromsp(geometry, P, S=None)

Read and return the object with possible overlap

geom

Associated geometry

geometry

Associated geometry

iter(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: 

In the above case io always belongs to atom ia and ia may be repeated according to the number of orbitals associated with the atom ia.

Parameters:
local : bool, optional

whether the orbital index is the global index, or the local index relative to the atom it resides on.

iter_nnz(atom=None, orbital=None)

Iterations of the non-zero elements

An iterator on the sparse matrix with, row and column

Parameters:
atom : int or array_like

only loop on the non-zero elements coinciding with the orbitals on these atoms (not compatible with the orbital keyword)

orbital : int or array_like

only loop on the non-zero elements coinciding with the orbital (not compatible with the atom keyword)

Examples

>>> for i, j in self.iter_nnz(): 
...    self[i, j] # is then the non-zero value 
make_hermitian()

Ensures the matrix is Hermitian by doing an in-place symmetrization

nnz

Number of non-zero elements

non_orthogonal

True if the object is using a non-orthogonal basis

nonzero(atom=None, only_col=False)

Indices row and column indices where non-zero elements exists

Parameters:
atom : int or array_like of int, optional

only return the tuples for the requested atoms, default is all atoms But for all orbitals.

only_col : bool, optional

only return then non-zero columns

See also

SparseCSR.nonzero
the equivalent function call
orthogonal

True if the object is using an orthogonal basis

static read(sile, *args, **kwargs)[source]

Reads dynamical matrix from Sile using read_dynamical_matrix.

Parameters:
sile : Sile, str

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 directly to read_dynamical_matrix(,**)
remove(atom, orb_index=None)

Remove a subset of this sparse matrix by only retaining the atoms corresponding to atom and optionally a subset of the atom orbitals

See sub for details regarding atom and orb_index arguments.

Parameters:
atom : array_like of int or Atom

indices of removed atoms or Atom for direct removal of all atoms

orb_index : array_like of int, optional

if atom is an instance of Atom, this variable correspond to the orbital indices for the atom to remove.

See also

Geometry.remove
equivalent to the resulting Geometry from this routine
Geometry.sub
the negative of Geometry.remove
sub
the opposite of remove, i.e. retain a subset of atoms
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

See also

Geometry.repeat
the same ordering as the final geometry
Geometry.tile
a different ordering of the final geometry
tile
a different ordering of the final geometry
reset(dim=None, dtype=<class 'numpy.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.

Parameters:
dim: int, optional

number of dimensions per element, default to the current number of elements per matrix element.

dtype: numpy.dtype, optional

the datatype of the sparse elements

nnzpr: int, optional

number of non-zero elements per row

rij(what='orbital', dtype=<class 'numpy.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 SuperCell.set_nsc for allowed parameters.

See also

SuperCell.set_nsc
the underlying called method
shape

Shape of sparse matrix

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(atom, orb_index=None)

Create a subset of this sparse matrix by only retaining the atoms corresponding to atom

Indices passed MUST be unique.

Negative indices are wrapped and thus works.

Parameters:
atom : array_like of int or Atom

indices of retained atoms or Atom for retaining only that atom

orb_index : array_like of int, optional

if atom is an instance of Atom, this variable correspond to the orbital indices for the atom to retain.

See also

Geometry.remove
the negative of Geometry.sub
Geometry.sub
equivalent to the resulting Geometry from this routine
remove
the negative of sub, i.e. remove a subset of atoms

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
>>> obj.sub(obj.atoms.atom[0], [1, 2]) # remove all but the 2nd and 3rd
>>>                                    # from the first atomic specie
>>>                                    # All other atomic species retain their orbitals.
swap(a, 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:
a : array_like

the first list of atomic coordinates

b : array_like

the second list of atomic coordinates

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

See also

Geometry.tile
the same ordering as the final geometry
Geometry.repeat
a different ordering of the final geometry
repeat
a different ordering of the final geometry
toSparseAtom(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:
dtype: numpy.dtype, optional

the data-container for the sparse object. Defaults to the same.

tocsr(index, isc=None, **kwargs)

Return a scipy.sparse.csr_matrix of the specified index

Parameters:
index : int

the index in the sparse matrix (for non-orthogonal cases the last dimension is the overlap matrix)

isc : int, optional

the supercell index, or all (if isc=None)

velocity(k=(0, 0, 0), **kwargs)[source]

Calculate the velocity for the eigenmodes for a given k point

Parameters:
k : array_like, optional

k-point at which the velocities are calculated

**kwargs: optional

additional parameters passed to the eigenmode routine

See also

eigenmode
method used to calculate the eigenmodes
displacement
Calculate mode displacements
EigenmodePhonon.velocity
Underlying method used to calculate the velocity
write(sile, *args, **kwargs)[source]

Writes a dynamical matrix to the Sile as implemented in the Sile.write_dynamical_matrix method