Hamiltonian

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

Sparse Hamiltonian matrix object

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

>>> ham = Hamiltonian(...)
>>> ham.H[1,2] = 0.1

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

Parameters
geometryGeometry

parent geometry to create a density matrix from. The density matrix will have size equivalent to the number of orbitals in the geometry

dimint or Spin, optional

number of components per element, may be a Spin object

dtypenp.dtype, optional

data type contained in the density matrix. See details of Spin for default values.

nnzprint, optional

number of initially allocated memory per orbital in the density matrix. For increased performance this should be larger than the actual number of entries per orbital.

spinSpin, optional

equivalent to dim argument. This keyword-only argument has precedence over dim.

orthogonalbool, optional

whether the density matrix corresponds to a non-orthogonal basis. In this case the dimensionality of the density matrix is one more than dim. This is a keyword-only argument.

Attributes

H

Access elements to the sparse Hamiltonian

S

Access elements to the sparse overlap

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

spin

Associated spin class

Methods

DOS(self, E[, k, distribution])

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

Hk(self[, k, dtype, gauge, format])

Setup the Hamiltonian for a given k-point

PDOS(self, E[, k, distribution])

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

Rij(self[, what, dtype])

Create a sparse matrix with the vectors between atoms/orbitals

Sk(self[, k, dtype, gauge, format])

Setup the overlap matrix for a given k-point

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

Initialize Hamiltonian

append(self, other, axis[, eps])

Append other along axis to construct a new connected sparse matrix

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 orbital model into different parts.

dHk(self[, k, dtype, gauge, format])

Setup the Hamiltonian derivative for a given k-point

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

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

ddHk(self[, k, dtype, gauge, format])

Setup the Hamiltonian double derivative for a given k-point

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

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

edges(self[, atom, exclude, orbital])

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

eig(self[, k, gauge, eigvals_only])

Returns the eigenvalues of the physical quantity (using the non-Hermitian solver)

eigenstate(self[, k, gauge])

Calculate the eigenstates at k and return an EigenstateElectron object containing all eigenstates

eigenvalue(self[, k, gauge])

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

eigh(self[, k, gauge, eigvals_only])

Returns the eigenvalues of the physical quantity

eigsh(self[, k, n, gauge, eigvals_only])

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

eliminate_zeros(self[, atol])

Removes all zero elements from the sparse matrix

empty(self[, keep_nnz])

See empty for details

fermi_level(self[, bz, q, distribution, q_tol])

Calculate the Fermi-level using a Brillouinzone sampling and a target charge

finalize(self)

Finalizes the model

fromsp(geometry, P[, S])

Read and return the object with possible overlap

iter(self[, local])

Iterations of the orbital space in the geometry, two indices from loop

iter_nnz(self[, atom, orbital])

Iterations of the non-zero elements

nonzero(self[, atom, only_col])

Indices row and column indices where non-zero elements exists

prepend(self, other, axis[, eps])

See append for details

read(sile, \*args, \*\*kwargs)

Reads Hamiltonian from Sile using read_hamiltonian.

remove(self, atom[, orb_index])

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

remove_orbital(self, atom, orbital)

Remove a subset of orbitals on atom according to orbital

repeat(self, reps, axis)

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

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

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

rij(self[, what, dtype])

Create a sparse matrix with the distance between atoms/orbitals

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

Reset the number of allowed supercells in the sparse orbital

shift(self, E)

Shift the electronic structure by a constant energy

spalign(self, other)

See align for details

spin_moment(self[, k])

Calculate the spin moment for the eigenstates for a given k point

spin_squared(self[, k, n_up, n_down])

Calculate spin-squared expectation value, see spin_squared 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 atoms corresponding to atom

sub_orbital(self, atom, orbital)

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

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 orbital object, equivalent to Geometry.tile

toSparseAtom(self[, dim, dtype])

Convert the sparse object (without data) to a new sparse object with equivalent but reduced sparse pattern

tocsr(self[, dim, isc])

Return a csr_matrix for the specified dimension

transpose(self)

Create the transposed sparse geometry by interchanging supercell indices

velocity(self[, k])

Calculate the velocity for the eigenstates for a given k point

write(self, sile, \*args, \*\*kwargs)

Writes a Hamiltonian to the Sile as implemented in the Sile.write_hamiltonian method

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

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

Parameters
Earray_like

energies to calculate the DOS at

karray_like, optional

k-point at which the DOS is calculated

distributionfunc or str, optional

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

**kwargsoptional

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

EigenvalueElectron.DOS

Underlying method used to calculate the DOS

property H

Access elements to the sparse Hamiltonian

Hk(self, k=(0, 0, 0), dtype=None, gauge='R', format='csr', *args, **kwargs)[source]

Setup the Hamiltonian for a given k-point

Creation and return of the Hamiltonian for a given k-point (default to Gamma).

Parameters
karray_like

the k-point to setup the Hamiltonian at

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

spinint, optional

if the Hamiltonian is a spin polarized one can extract the specific spin direction matrix by passing an integer (0 or 1). If the Hamiltonian is not Spin.POLARIZED this keyword is ignored.

Returns
objectthe Hamiltonian matrix at \(k\). The returned object depends on format.

See also

dHk

Hamiltonian derivative with respect to k

ddHk

Hamiltonian double derivative with respect to k

Notes

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

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

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

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

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

Parameters
Earray_like

energies to calculate the projected DOS at

karray_like, optional

k-point at which the projected DOS is calculated

distributionfunc or str, optional

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

**kwargsoptional

additional parameters passed to the eigenstate routine

See also

sisl.physics.distribution

setup a distribution function, see details regarding the distribution argument

eigenstate

method used to calculate the eigenstates

DOS

Calculate total DOS

EigenstateElectron.PDOS

Underlying method used to calculate the projected DOS

Rij(self, 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.

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.

property S

Access elements to the sparse overlap

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

the k-point to setup the overlap at (default Gamma point)

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

Returns
objectthe 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.

append(self, other, axis, eps=0.01)

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

Parameters
otherobject

must be an object of the same type as self

axisint

axis to append the two sparse geometries along

epsfloat, 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.

Returns
object

a new instance with two sparse matrices joined and appended together

Raises
ValueError if atomic coordinates does not overlap within eps

See also

prepend

equivalent scheme as this method

transpose

ensure hermiticity by using this routine

Geometry.append
Geometry.prepend

Notes

This routine and how it is functioning may change in future releases. There are many design choices in how to assign the matrix elements when combining two models and it is not clear what is the best procedure.

The current implentation does not preserve the hermiticity of the matrix.

Examples

>>> sporb = SparseOrbital(....)
>>> forced_hermitian = (sporb + sporb.transpose()) * 0.5
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)

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
sepsint

number of times the structure will be cut

axisint

the axis that will be cut

dHk(self, k=(0, 0, 0), dtype=None, gauge='R', format='csr', *args, **kwargs)[source]

Setup the Hamiltonian derivative for a given k-point

Creation and return of the Hamiltonian derivative for a given k-point (default to Gamma).

Parameters
karray_like

the k-point to setup the Hamiltonian at

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

spinint, optional

if the Hamiltonian is a spin polarized one can extract the specific spin direction matrix by passing an integer (0 or 1). If the Hamiltonian is not Spin.POLARIZED this keyword is ignored.

Returns
tuplefor each of the Cartesian directions a \(\partial \mathbf H(k)/\partial k_\alpha\) is returned.

See also

Hk

Hamiltonian with respect to k

ddHk

Hamiltonian double derivative with respect to k

Notes

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

\[\nabla_k \mathbf H_\alpha(k) = i R_\alpha \mathbf H_{\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 H_\alpha(k) = i r_\alpha \mathbf H_{\nu\mu} e^{i k r}\]

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

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

the k-point to setup the overlap at (default Gamma point)

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

Returns
tuplefor 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.

ddHk(self, k=(0, 0, 0), dtype=None, gauge='R', format='csr', *args, **kwargs)[source]

Setup the Hamiltonian double derivative for a given k-point

Creation and return of the Hamiltonian double derivative for a given k-point (default to Gamma).

Parameters
karray_like

the k-point to setup the Hamiltonian at

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

spinint, optional

if the Hamiltonian is a spin polarized one can extract the specific spin direction matrix by passing an integer (0 or 1). If the Hamiltonian is not Spin.POLARIZED this keyword is ignored.

Returns
tuple of tuplesfor each of the Cartesian directions

See also

Hk

Hamiltonian with respect to k

dHk

Hamiltonian derivative with respect to k

Notes

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

\[\nabla_k^2 \mathbf H_{\alpha\beta}(k) = - R_\alpha R_\beta \mathbf H_{\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 H_{\alpha\beta}(k) = - r_\alpha r_\beta \mathbf H_{\nu\mu} e^{i k r}\]

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

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

the k-point to setup the overlap at (default Gamma point)

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

Returns
tuple of tuplesfor 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.

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=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
atomint 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.

excludeint or list of int, optional

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

orbitalint 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(self, 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

Parameters
spinint, optional

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

eigenstate(self, k=(0, 0, 0), gauge='R', **kwargs)[source]

Calculate the eigenstates at k and return an EigenstateElectron object containing all eigenstates

Parameters
karray_like*3, optional

the k-point at which to evaluate the eigenstates at

gaugestr, optional

the gauge used for calculating the eigenstates

sparsebool, optional

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

**kwargsdict, optional

passed arguments to the eigh routine

Returns
EigenstateElectron

See also

eigh

eigenvalue routine

eigsh

eigenvalue routine

eigenvalue(self, k=(0, 0, 0), gauge='R', **kwargs)[source]

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

Parameters
karray_like*3, optional

the k-point at which to evaluate the eigenvalues at

gaugestr, optional

the gauge used for calculating the eigenvalues

sparsebool, optional

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

**kwargsdict, optional

passed arguments to the eigh routine

Returns
EigenvalueElectron

See also

eigh

eigenvalue routine

eigsh

eigenvalue routine

eigh(self, 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

Parameters
spinint, optional

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

eigsh(self, 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

Parameters
spinint, optional

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

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

fermi_level(self, bz=None, q=None, distribution='fermi_dirac', q_tol=1e-12)[source]

Calculate the Fermi-level using a Brillouinzone sampling and a target charge

The Fermi-level will be calculated using an iterative approach by first calculating all eigenvalues and subsequently fitting the Fermi level to the final charge (q).

Parameters
bzBrillouinzone, optional

sampled k-points and weights, the bz.parent will be equal to this object upon return default to Gamma-point

qfloat, list of float, optional

seeked charge, if not set will be equal to self.geometry.q0. If a list of two is passed there will be calculated a Fermi-level per spin-channel. If the Hamiltonian is not spin-polarized the sum of the list will be used and only a single fermi-level will be returned.

distributionstr, func, optional

used distribution, must accept the keyword mu as parameter for the Fermi-level

q_tolfloat, optional

tolerance of charge for finding the Fermi-level

Returns
fermi-levelthe Fermi-level of the system (or two if two different charges are passed)
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(geometry, P, S=None)

Read and return the object with possible overlap

property geom

Associated geometry

property geometry

Associated geometry

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

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

iter_nnz(self, atom=None, orbital=None)

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 orbitals on these atoms (not compatible with the orbital keyword)

orbitalint 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
property nnz

Number of non-zero elements

property non_orthogonal

True if the object is using a non-orthogonal basis

nonzero(self, atom=None, only_col=False)

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 But for all orbitals.

only_colbool, 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(self, other, axis, eps=0.01)

See append for details

This is currently equivalent to:

>>> other.append(self, axis, eps)
static read(sile, *args, **kwargs)[source]

Reads Hamiltonian from Sile using read_hamiltonian.

Parameters
sileSile, str or pathlib.Path

a Sile object which will be used to read the Hamiltonian and the overlap matrix (if any) if it is a string it will create a new sile using get_sile.

*args passed directly to read_hamiltonian(,**)
remove(self, atom, orb_index=None)

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

Parameters
atomarray_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(self, atom, orbital)

Remove a subset of orbitals on atom according to orbital

Parameters
atomarray_like of int or Atom

indices of atoms or Atom that will be reduced in size according to orbital

orbitalarray_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(self, 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
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, 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.

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)

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

shift(self, E)[source]

Shift the electronic structure by a constant energy

This is equal to performing this operation:

\[\mathbf H_\sigma = \mathbf H_\sigma + E \mathbf S\]

where \(\mathbf H_\sigma\) correspond to the spin diagonal components of the Hamiltonian.

Parameters
Efloat or (2,)

the energy (in eV) to shift the electronic structure, if two values are passed the two first spin-components get shifted individually.

spalign(self, other)

See align for details

property spin

Associated spin class

spin_moment(self, k=(0, 0, 0), **kwargs)[source]

Calculate the spin moment for the eigenstates for a given k point

Parameters
karray_like, optional

k-point at which the spin moments are calculated

**kwargsoptional

additional parameters passed to the eigenstate routine

See also

eigenstate

method used to calculate the eigenstates

EigenvalueElectron.spin_moment

Underlying method used to calculate the spin moment

spin_squared(self, k=(0, 0, 0), n_up=None, n_down=None, **kwargs)[source]

Calculate spin-squared expectation value, see spin_squared for details

Parameters
karray_like, optional

k-point at which the spin-squared expectation value is

n_upint, optional

number of states for spin up configuration, default to all. All states up to and including n_up.

n_downint, optional

same as n_up but for the spin-down configuration

**kwargsoptional

additional parameters passed to the eigenstate routine

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)

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
atomarray_like of int or Atom

indices of retained atoms or Atom for retaining only that atom

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
sub_orbital(self, atom, orbital)

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
atomarray_like of int or Atom

indices of atoms or Atom that will be reduced in size according to orbital

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

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

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

number of dimensions allocated in the SparseAtom object, default to the same

dtypenumpy.dtype, optional

used data-type for the sparse object. Defaults to the same.

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
velocity(self, k=(0, 0, 0), **kwargs)[source]

Calculate the velocity for the eigenstates for a given k point

Parameters
karray_like, optional

k-point at which the velocities are calculated

**kwargsoptional

additional parameters passed to the eigenstate routine

See also

eigenstate

method used to calculate the eigenstates

EigenvalueElectron.velocity

Underlying method used to calculate the velocity

write(self, sile, *args, **kwargs)[source]

Writes a Hamiltonian to the Sile as implemented in the Sile.write_hamiltonian method