sisl.physics.hamiltonian module

Tight-binding class to create tight-binding models.

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

Bases: sisl.physics.sparse_physics.SparseOrbitalBZSpin

Hamiltonian object containing the coupling constants between orbitals.

The Hamiltonian object contains information regarding the
  • geometry
  • coupling constants between orbitals

It contains an intrinsic sparse matrix of the Hamiltonian elements.

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

Attributes

H
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
orthogonal True if the object is using an orthogonal basis
shape Shape of sparse matrix
spin Spin class

Methods

Hk([k, dtype, gauge, format]) Setup the Hamiltonian for a given k-point
Sk([k, dtype, gauge, format]) Setup the overlap matrix for a given k-point
align(other) See SparseCSR.align for details
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.
eigh([k, atoms, gauge, eigvals_only, ...]) Returns the eigenvalues of the physical quantity
eigsh([k, n, atoms, gauge, eigvals_only]) Calculates a subset of eigenvalues of the physical quantity (default 10)
eliminate_zeros() Removes all zero elements from the sparse matrix
empty([keep]) See SparseCSR.empty for details
finalize() Finalizes the model
fromsp(geom, 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
read(sile, *args, **kwargs) Reads Hamiltonian from Sile using read_hamiltonian.
repeat(reps, axis) Create a repeated sparse orbital object, equivalent to Geometry.repeat
reset([dim, dtype, nnzpr]) The sparsity pattern is cleaned and every thing is reset.
shift(E) Shift the electronic structure by a constant energy
spsame(other) Compare two sparse objects and check whether they have the same entries.
sub(atom) Create a subset of this sparse matrix by only retaining the elements corresponding to the atom
tile(reps, axis) Create a tiled sparse orbital object, equivalent to Geometry.tile
tocsr(index[, isc]) Return a scipy.sparse.csr_matrix of the specified index
write(sile, *args, **kwargs) Writes a tight-binding model to the Sile as implemented in the Sile.write_hamiltonian method

Create Hamiltonian model from geometry

Initializes a Hamiltonian using the geom object as the underlying geometry for the tight-binding parameters.

Attributes

H
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
orthogonal True if the object is using an orthogonal basis
shape Shape of sparse matrix
spin Spin class

Methods

Hk([k, dtype, gauge, format]) Setup the Hamiltonian for a given k-point
Sk([k, dtype, gauge, format]) Setup the overlap matrix for a given k-point
align(other) See SparseCSR.align for details
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.
eigh([k, atoms, gauge, eigvals_only, ...]) Returns the eigenvalues of the physical quantity
eigsh([k, n, atoms, gauge, eigvals_only]) Calculates a subset of eigenvalues of the physical quantity (default 10)
eliminate_zeros() Removes all zero elements from the sparse matrix
empty([keep]) See SparseCSR.empty for details
finalize() Finalizes the model
fromsp(geom, 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
read(sile, *args, **kwargs) Reads Hamiltonian from Sile using read_hamiltonian.
repeat(reps, axis) Create a repeated sparse orbital object, equivalent to Geometry.repeat
reset([dim, dtype, nnzpr]) The sparsity pattern is cleaned and every thing is reset.
shift(E) Shift the electronic structure by a constant energy
spsame(other) Compare two sparse objects and check whether they have the same entries.
sub(atom) Create a subset of this sparse matrix by only retaining the elements corresponding to the atom
tile(reps, axis) Create a tiled sparse orbital object, equivalent to Geometry.tile
tocsr(index[, isc]) Return a scipy.sparse.csr_matrix of the specified index
write(sile, *args, **kwargs) Writes a tight-binding model to the Sile as implemented in the Sile.write_hamiltonian method
H
Hk(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:

k : array_like

the k-point to setup the Hamiltonian 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').

See also

Sk
Overlap matrix at k

Notes

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

\[H(k) = H_{ij} e^{i k R}\]

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

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

\[H(k) = H_{ij} e^{i k r}\]

where \(r\) is the distance between the orbitals \(i\) and \(j\). Currently the second gauge is not implemented (yet).

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

Reads Hamiltonian from Sile using read_hamiltonian.

Parameters:

sile : Sile, str

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

shift(E)[source]

Shift the electronic structure by a constant energy

Parameters:

E : float

the energy (in eV) to shift the electronic structure

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

Writes a tight-binding model to the Sile as implemented in the Sile.write_hamiltonian method

sisl.physics.hamiltonian.TightBinding

alias of Hamiltonian