SparseOrbitalBZSpin

class sisl.physics.SparseOrbitalBZSpin(geometry, dim=1, dtype=None, nnzpr=None, **kwargs)

Bases: sisl.physics.SparseOrbitalBZ

Sparse object containing the orbital connections in a Brillouin zone with possible spin-components

It contains an intrinsic sparse matrix of the physical elements.

Assigning or changing elements is as easy as with standard numpy assignments:

>>> S = SparseOrbitalBZSpin(...)
>>> S[1,2] = 0.1

which assigns 0.1 as the element between orbital 2 and 3. (remember that Python is 0-based elements).

Parameters
  • geometry (Geometry) – parent geometry to create a sparse matrix from. The matrix will have size equivalent to the number of orbitals in the geometry

  • dim (int or Spin, optional) – number of components per element, may be a Spin object

  • dtype (np.dtype, optional) – data type contained in the matrix. See details of Spin for default values.

  • nnzpr (int, optional) – number of initially allocated memory per orbital in the matrix. For increased performance this should be larger than the actual number of entries per orbital.

  • spin (Spin, optional) – equivalent to dim argument. This keyword-only argument has precedence over dim.

  • orthogonal (bool, optional) – whether the matrix corresponds to a non-orthogonal basis. In this case the dimensionality of the matrix is one more than dim. This is a keyword-only argument.

Attributes

S

Access the overlap elements associated with the sparse matrix

__array_priority__

__dict__

__doc__

__hash__

__module__

__weakref__

list of weak references to the object (if defined)

_size

The size of the sparse object

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

deprecated 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

spin

Associated spin class

Methods

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

_Pk([k, dtype, gauge, format, _dim])

Sparse matrix (scipy.sparse.csr_matrix) at k for a polarized system

_Pk_non_colinear([k, dtype, gauge, format])

Sparse matrix (scipy.sparse.csr_matrix) at k for a non-collinear system

_Pk_polarized([k, spin, dtype, gauge, format])

Sparse matrix (scipy.sparse.csr_matrix) at k for a polarized system

_Pk_spin_orbit([k, dtype, gauge, format])

Sparse matrix (scipy.sparse.csr_matrix) at k for a spin-orbit system

_Pk_unpolarized([k, dtype, gauge, format])

Sparse matrix (scipy.sparse.csr_matrix) at k

_Sk([k, dtype, gauge, format])

Overlap matrix in a scipy.sparse.csr_matrix at k.

_Sk_diagonal([k, dtype, gauge, format])

For an orthogonal case we always return the identity matrix

_Sk_non_colinear([k, dtype, gauge, format])

Overlap matrix (scipy.sparse.csr_matrix) at k for a non-collinear system

__abs__()

__add__(other)

__and__(other)

__array_ufunc__(ufunc, method, *inputs, **kwargs)

__contains__(key)

Check whether a sparse index is non-zero

__delattr__

Implement delattr(self, name).

__delitem__(key)

Delete elements of the sparse elements

__dir__

Default dir() implementation.

__divmod__(other)

__eq__(other)

__floordiv__(other)

__format__

Default object formatter.

__ge__(other)

__getattr__(attr)

Overload attributes from the hosting geometry

__getattribute__

Return getattr(self, name).

__getitem__(key)

Elements for the index(s)

__getstate__()

Return dictionary with the current state

__gt__(other)

__iadd__(other)

__iand__(other)

__ifloordiv__(other)

__ilshift__(other)

__imatmul__(other)

__imod__(other)

__imul__(other)

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

Create sparse object with element between orbitals

__init_subclass__

This method is called when a class is subclassed.

__invert__()

__ior__(other)

__ipow__(other)

__irshift__(other)

__isub__(other)

__iter__()

Iterations of the non-zero elements

__itruediv__(other)

__ixor__(other)

__le__(other)

__len__()

Returns number of rows in the basis (if non-collinear or spin-orbit, twice the number of orbitals)

__lshift__(other)

__lt__(other)

__matmul__(other)

__mod__(other)

__mul__(other)

__ne__(other)

__neg__()

__new__

Create and return a new object.

__or__(other)

__pos__()

__pow__(other)

__radd__(other)

__rand__(other)

__rdivmod__(other)

__reduce__

Helper for pickle.

__reduce_ex__

Helper for pickle.

__repr__()

Return repr(self).

__rfloordiv__(other)

__rlshift__(other)

__rmatmul__(other)

__rmod__(other)

__rmul__(other)

__ror__(other)

__rpow__(other)

__rrshift__(other)

__rshift__(other)

__rsub__(other)

__rtruediv__(other)

__rxor__(other)

__setattr__

Implement setattr(self, name, value).

__setitem__(key, val)

Set or create elements in the sparse data

__setstate__(state)

Return dictionary with the current state

__sizeof__

Size of object in memory, in bytes.

__str__()

Representation of the model

__sub__(other)

__subclasshook__

Abstract classes can override this to customize issubclass().

__truediv__(other)

__xor__(other)

_cls_kwargs()

Custom keyword arguments when creating a new instance

_dPk([k, dtype, gauge, format, _dim])

Sparse matrix (scipy.sparse.csr_matrix) at k differentiated with respect to k for a polarized system

_dPk_non_colinear([k, dtype, gauge, format])

Tuple of sparse matrix (scipy.sparse.csr_matrix) at k for a non-collinear system, differentiated with respect to k

_dPk_polarized([k, spin, dtype, gauge, format])

Tuple of sparse matrix (scipy.sparse.csr_matrix) at k, differentiated with respect to k

_dPk_spin_orbit([k, dtype, gauge, format])

Tuple of sparse matrix (scipy.sparse.csr_matrix) at k for a non-collinear system, differentiated with respect to k

_dPk_unpolarized([k, dtype, gauge, format])

Tuple of sparse matrix (scipy.sparse.csr_matrix) at k, differentiated with respect to k

_dSk([k, dtype, gauge, format])

Overlap matrix in a scipy.sparse.csr_matrix at k differentiated with respect to k

_dSk_non_colinear([k, dtype, gauge, format])

Overlap matrix (scipy.sparse.csr_matrix) at k for a non-collinear system

_ddPk([k, dtype, gauge, format, _dim])

Sparse matrix (scipy.sparse.csr_matrix) at k double differentiated with respect to k for a polarized system

_ddPk_non_colinear([k, dtype, gauge, format])

Tuple of sparse matrix (scipy.sparse.csr_matrix) at k for a non-collinear system, differentiated with respect to k twice

_ddPk_spin_orbit([k, dtype, gauge, format])

Tuple of sparse matrix (scipy.sparse.csr_matrix) at k for a non-collinear system, differentiated with respect to k

_ddSk([k, dtype, gauge, format])

Overlap matrix in a scipy.sparse.csr_matrix at k double differentiated with respect to k

_ddSk_non_colinear([k, dtype, gauge, format])

Overlap matrix in a scipy.sparse.csr_matrix at k for non-collinear spin, differentiated with respect to k

_reset()

Reset object according to the options, please refer to SparseOrbital.reset for details

_translate_cells(old, new)

Translates all columns in the old cell indices to the new cell indices

add(other[, axis, offset])

Add two sparse matrices by adding the parameters to one set.

append(other, axis[, eps, 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, param)

Create a simple function for passing to the construct function.

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

Cuts the sparse orbital model into different parts.

dSk([k, dtype, gauge, format])

Setup the \(k\)-derivatie of the overlap matrix for a given k-point

ddSk([k, dtype, gauge, format])

Setup the double \(k\)-derivatie of the overlap matrix for a given k-point

edges([atoms, exclude, orbitals])

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)

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(*args, **kwargs)

Removes all zero elements from the sparse matrix

empty([keep_nnz])

See empty for details

finalize()

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_col])

Indices row and column indices where non-zero elements exists

prepend(other, axis[, eps])

See append for details

remove(atoms)

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

remove_orbital(atoms, orbitals)

Remove a subset of orbitals on atom according to orbital

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(atoms)

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

sub_orbital(atoms, orbitals)

Retain only a subset of the orbitals on atom according to orbital

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([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 dimension

transpose([hermitian])

A transpose copy of this object, possibly apply the Hermitian conjugate as well

trs()

Create a new matrix with applied time-reversal-symmetry

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.

property S

Access the overlap elements associated with the sparse matrix

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

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.

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’/’dense’/’matrix’).

See also

dSk

Overlap matrix derivative with respect to k

ddSk

Overlap matrix double derivative with respect to k

Returns

matrix – the overlap matrix at \(k\). The returned object depends on format.

Return type

numpy.ndarray or scipy.sparse.*_matrix

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 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 appended

  • axis (int or None, optional) – 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 ((3,), optional) – offset in geometry of other when adding the atoms.

See also

append

append two matrices by also adding overlap couplings

prepend

see append

append(other, axis, eps=0.01, 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:

  1. Create the new extended geometry

  2. 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 eps 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()) * 0.5

To retain couplings only from the right sparse matrix, do:

>>> sporb = left.append(right, 0, scale=(0, 2.))
>>> sporb = (sporb + sporb.transpose()) * 0.5

Notes

The current implentation does not preserve the hermiticity of the matrix. If you want to preserve hermiticity of the matrix you have to do the following:

>>> h = (h + h.transpose()) / 2
Parameters
  • other (object) – must be an object of the same type as self

  • axis (int) – axis to append the two sparse geometries along

  • eps (float, optional) – 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 eps is used as a baseline, see above.

  • scale (float or array_like, optional) – the scale used for the overlapping region. For scalar values it corresponds to passing: (scale, scale). For array-like input scale[0] refers to the scale of the matrix elements coupling from self, while scale[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

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

a new instance with two sparse matrices joined and appended together

Return type

object

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, 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 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, atoms, atoms_xyz=None):
...     idx = self.geometry.close(ia, R=R, atoms=atoms, atoms_xyz=atoms_xyz)
...     for ix, p in zip(idx, param):
...         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 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)

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

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

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.

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’/’dense’/’matrix’).

See also

Sk

Overlap matrix at k

ddSk

Overlap matrix double derivative at k

Returns

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

Return type

tuple

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

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.

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’/’dense’/’matrix’).

See also

Sk

Overlap matrix at k

dSk

Overlap matrix derivative at k

Returns

for each of the Cartesian directions

Return type

tuple of tuples

property dim

Number of components per element

property dkind

Data type of sparse elements (in str)

property dtype

Data type of sparse elements

edges(atoms=None, exclude=None, orbitals=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
  • atoms (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 or None, optional) – remove edges which are in the exclude list, this list refers to orbitals.

  • 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)[source]

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

Parameters

spin (int, optional) – the spin-component to calculate the eigenvalue spectrum of, note that this parameter is only valid for Spin.POLARIZED matrices.

eigh(k=0, 0, 0, gauge='R', eigvals_only=True, **kwargs)[source]

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

Parameters

spin (int, optional) – the spin-component to calculate the eigenvalue spectrum of, note that this parameter is only valid for Spin.POLARIZED matrices.

eigsh(k=0, 0, 0, n=10, gauge='R', eigvals_only=True, **kwargs)[source]

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

Parameters

spin (int, optional) – the spin-component to calculate the eigenvalue spectrum of, note that this parameter is only valid for Spin.POLARIZED matrices.

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

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.

property finalized

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

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.sparse formats.

Parameters
  • geometry (Geometry) – geometry to describe the new sparse geometry

  • P (list of scipy.sparse or scipy.sparse) – the new sparse matrices that are to be populated in the sparse matrix

  • S (scipy.sparse, optional) – if provided this refers to the overlap matrix and will force the returned sparse matrix to be non-orthogonal

  • **kwargs (optional) – 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 S is not none

Return type

SparseGeometry

property geom

deprecated geometry

property geometry

Associated geometry

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 (int or array_like) – only loop on the non-zero elements coinciding with the orbitals on these atoms (not compatible with the orbital keyword)

  • orbitals (int or array_like) – only loop on the non-zero elements coinciding with the orbital (not compatible with the atom 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 io always belongs to atom ia and ia may be repeated according to the number of orbitals associated with the atom ia.

Parameters
  • atoms (int or array_like, optional) – only loop on the given atoms, default to all atoms

  • local (bool, optional) – whether the orbital index is the global index, or the local index relative to the atom it resides on.

Yields
  • ia – atomic index

  • io – orbital index

See also

Geometry.iter_orbitals

method used to iterate orbitals

property nnz

Number of non-zero elements

property non_orthogonal

True if the object is using a non-orthogonal basis

nonzero(atoms=None, only_col=False)

Indices row and column indices where non-zero elements exists

Parameters
  • atoms (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

property orthogonal

True if the object is using an orthogonal basis

prepend(other, axis, eps=0.01)

See append for details

This is currently equivalent to:

>>> other.append(self, axis, eps)
remove(atoms)

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

Parameters

atoms (array_like of int or Atom) – indices of removed atoms or Atom for direct removal of all 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

remove_orbital(atoms, orbitals)

Remove a subset of orbitals on atom according to orbital

Parameters
  • atoms (array_like of int or Atom) – indices of atoms or Atom that will be reduced in size according to orbital

  • orbitals (array_like of int or Orbital) – indices of the orbitals on atom that are removed from the sparse matrix.

Examples

>>> obj = SparseOrbital(...)
>>> # remove the second orbital on the 2nd atom
>>> # all other orbitals are retained
>>> obj.remove_orbital(1, 1)
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

property shape

Shape of sparse matrix

spalign(other)

See align for details

property spin

Associated spin class

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 atom

Negative indices are wrapped and thus works, supercell atoms are also wrapped to the unit-cell.

Parameters

atoms (array_like of int or Atom) – 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 routine

remove

the negative of sub, i.e. remove a subset of atoms

sub_orbital(atoms, orbitals)

Retain only a subset of the orbitals on atom according to orbital

This allows one to retain only a given subset of the sparse matrix elements.

Parameters
  • atoms (array_like of int or Atom) – indices of atoms or Atom that will be reduced in size according to orbital

  • orbitals (array_like of int or Orbital) – indices of the orbitals on atom 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-arange orbitals using this method.

Examples

>>> obj = SparseOrbital(...)
>>> # only retain the second orbital on the 2nd atom
>>> # all other orbitals are retained
>>> obj.sub_orbital(1, 1)
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(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, optional) – 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_matrix for the specified dimension

Parameters
  • dim (int, optional) – 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 (if isc=None)

transpose(hermitian=False)[source]

A transpose copy of this object, possibly apply the Hermitian conjugate as well

Parameters

hermitian (bool, optional) – if true, also emply a spin-box Hermitian operator to ensure TRS, otherwise only return the transpose values.

trs()[source]

Create a new matrix with applied time-reversal-symmetry

Time reversal symmetry is applied using the following equality:

\[2\mathbf M^{\mathrm{TRS}} = \mathbf M + \boldsymbol\sigma_y \mathbf M^* \boldsymbol\sigma_y\]

where \(*\) is the conjugation operator.