SparseAtom

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

Sparse object with number of rows equal to the total number of atoms in the Geometry

Attributes

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

shape

Shape of sparse matrix

Methods

Rij(self[, dtype])

Create a sparse matrix with the vectors between atoms

__init__(self, geometry[, dim, dtype, nnzpr])

Create sparse object with element between orbitals

construct(self, func[, na_iR, method, eta])

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

copy(self[, dtype])

A copy of this object

create_construct(self, R, param)

Create a simple function for passing to the construct function.

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

Cuts the sparse atom model into different parts.

edges(self, atom[, exclude])

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

eliminate_zeros(self[, atol])

Removes all zero elements from the sparse matrix

empty(self[, keep_nnz])

See empty for details

finalize(self)

Finalizes the model

fromsp(geom, \*sp)

Create a sparse model from a preset Geometry and a list of sparse matrices

iter_nnz(self[, atom])

Iterations of the non-zero elements

nonzero(self[, atom, only_col])

Indices row and column indices where non-zero elements exists

remove(self, atom)

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

repeat(self, reps, axis)

Create a repeated sparse atom object, equivalent to Geometry.repeat

reset(self[, dim, dtype, nnzpr])

The sparsity pattern has all elements removed and everything is reset.

rij(self[, dtype])

Create a sparse matrix with the distance between atoms

set_nsc(self, \*args, \*\*kwargs)

Reset the number of allowed supercells in the sparse atom

spalign(self, other)

See align for details

spsame(self, other)

Compare two sparse objects and check whether they have the same entries.

sub(self, atom)

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

swap(self, a, b)

Swaps atoms in the sparse geometry to obtain a new order of atoms

tile(self, reps, axis)

Create a tiled sparse atom object, equivalent to Geometry.tile

tocsr(self[, dim, isc])

Return a csr_matrix for the specified dimension

transpose(self)

Create the transposed sparse geometry by interchanging supercell indices

Rij(self, dtype=<class 'numpy.float64'>)[source]

Create a sparse matrix with the vectors between atoms

Parameters
dtypenumpy.dtype, optional

the data-type of the sparse matrix.

Notes

The returned sparse matrix with vectors are taken from the current sparse pattern. I.e. a subsequent addition of sparse elements will make them inequivalent. It is thus important to only create the sparse vector matrix when the sparse structure is completed.

construct(self, 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
funccallable 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_iRint, 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

etabool, 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(self, dtype=None)

A copy of this object

Parameters
dtypenumpy.dtype, optional

it is possible to convert the data to a different data-type If not specified, it will use self.dtype

create_construct(self, 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
Rarray_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

paramarray_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(self, seps, axis, *args, **kwargs)[source]

Cuts the sparse atom model into different parts.

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

Cutting is the opposite of tiling.

Parameters
sepsint

number of times the structure will be cut

axisint

the axis that will be cut

property dim

Number of components per element

property dkind

Data type of sparse elements (in str)

property dtype

Data type of sparse elements

edges(self, atom, exclude=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.na_s).

Parameters
atomint or list of int

the edges are returned only for the given atom

excludeint or list of int, optional

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

See also

SparseCSR.edges

the underlying routine used for extracting the edges

eliminate_zeros(self, atol=0.0)

Removes all zero elements from the sparse matrix

This is an in-place operation.

Parameters
atolfloat, optional

absolute tolerance equal or below this value will be considered 0.

empty(self, keep_nnz=False)

See empty for details

finalize(self)

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.

property finalized

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

classmethod fromsp(geom, *sp)

Create a sparse model from a preset Geometry and a list of sparse matrices

property geom

Associated geometry

property geometry

Associated geometry

iter_nnz(self, atom=None)[source]

Iterations of the non-zero elements

An iterator on the sparse matrix with, row and column

Parameters
atomint or array_like

only loop on the non-zero elements coinciding with the atoms

Examples

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

Number of non-zero elements

nonzero(self, atom=None, only_col=False)[source]

Indices row and column indices where non-zero elements exists

Parameters
atomint or array_like of int, optional

only return the tuples for the requested atoms, default is all atoms

only_colbool, optional

only return then non-zero columns

See also

SparseCSR.nonzero

the equivalent function call

remove(self, atom)

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

Negative indices are wrapped and thus works.

Parameters
atomarray_like of int

indices of removed atoms

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(self, reps, axis)[source]

Create a repeated sparse atom 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
repsint

number of repetitions along cell-vector axis

axisint

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(self, 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
dimint, optional

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

dtypenumpy.dtype, optional

the datatype of the sparse elements

nnzprint, optional

number of non-zero elements per row

rij(self, dtype=<class 'numpy.float64'>)[source]

Create a sparse matrix with the distance between atoms

Parameters
dtypenumpy.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(self, *args, **kwargs)[source]

Reset the number of allowed supercells in the sparse atom

If one reduces the number of supercells any sparse element that references the supercell will be deleted.

See SuperCell.set_nsc for allowed parameters.

See also

SuperCell.set_nsc

the underlying called method

property shape

Shape of sparse matrix

spalign(self, other)

See align for details

spsame(self, 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(self, atom)[source]

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

Indices passed MUST be unique.

Negative indices are wrapped and thus works.

Parameters
atomarray_like of int

indices of retained atoms

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

swap(self, 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
aarray_like

the first list of atomic coordinates

barray_like

the second list of atomic coordinates

tile(self, reps, axis)[source]

Create a tiled sparse atom 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
repsint

number of repetitions along cell-vector axis

axisint

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

Notes

Calling this routine will automatically finalize the SparseAtom. This is required to greatly increase performance.

tocsr(self, dim=0, isc=None, **kwargs)

Return a csr_matrix for the specified dimension

Parameters
dimint, optional

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

iscint, optional

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

transpose(self)

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

Returns
object

an equivalent sparse geometry with transposed matrix elements

Notes

For Hamiltonians with non-collinear or spin-orbit there is no transposing of the sub-spin matrix box. This needs to be done manually.

Examples

Force a sparse geometry to be Hermitian:

>>> sp = SparseOrbital(...)
>>> sp = (sp + sp.transpose()) * 0.5