sisl.quantity.hamiltonian module

Tight-binding class to create tight-binding models.

class sisl.quantity.hamiltonian.Hamiltonian(geom, nnzpr=None, orthogonal=True, spin=1, dtype=None, *args, **kwargs)[source]

Bases: object

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
dtype Return data type of Hamiltonian (and overlap matrix)
finalized Whether the contained data is finalized and non-used elements have been removed
geom Return the attached geometry
geometry Return the attached geometry
nnz Returns number of non-zero elements in the tight-binding model
no Returns number of orbitals as used when the object was created
orthogonal Return whether the Hamiltonian is orthogonal
spin Return number of spin-components in Hamiltonian

Methods

construct(func[, na_iR, method, eta]) Automatically construct the Hamiltonian model based on a function that does the setting up of the Hamiltonian
copy([dtype]) Return a copy of the Hamiltonian object
create_construct(dR, param) Returns a simple function for passing to the construct function.
cut(seps, axis, *args, **kwargs) Cuts the tight-binding model into different parts.
eigh([k, atoms, eigvals_only, overwrite_a, ...]) Returns the eigenvalues of the Hamiltonian
eigsh([k, n, atoms, eigvals_only]) Returns the eigenvalues of the Hamiltonian
empty([keep]) See SparseCSR.empty for details
finalize() Finalizes the tight-binding model
iter([local]) Iterations of the orbital space in the geometry, two indices from loop
iter_nnz([atom, orbital]) Iterations of the non-zero elements, returns a tuple of orbital and coupling orbital
read(sile, *args, **kwargs) Reads Hamiltonian from Sile using read_H.
repeat(reps, axis) Refer to tile instead
reset([nnzpr, orthogonal, spin, dtype]) The sparsity pattern is cleaned and every thing is reset.
sp2HS(geom, H[, S]) Returns a tight-binding model from a preset H, S and Geometry
tile(reps, axis) Returns a repeated tight-binding model for this, much like the Geometry
tocsr(index[, isc]) Return a scipy.sparse.csr_matrix from the specified index
write(sile, *args, **kwargs) Writes a tight-binding model to the Sile as implemented in the Sile.write_es method

Create tight-binding model from geometry

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

Attributes

H
S
dtype Return data type of Hamiltonian (and overlap matrix)
finalized Whether the contained data is finalized and non-used elements have been removed
geom Return the attached geometry
geometry Return the attached geometry
nnz Returns number of non-zero elements in the tight-binding model
no Returns number of orbitals as used when the object was created
orthogonal Return whether the Hamiltonian is orthogonal
spin Return number of spin-components in Hamiltonian

Methods

construct(func[, na_iR, method, eta]) Automatically construct the Hamiltonian model based on a function that does the setting up of the Hamiltonian
copy([dtype]) Return a copy of the Hamiltonian object
create_construct(dR, param) Returns a simple function for passing to the construct function.
cut(seps, axis, *args, **kwargs) Cuts the tight-binding model into different parts.
eigh([k, atoms, eigvals_only, overwrite_a, ...]) Returns the eigenvalues of the Hamiltonian
eigsh([k, n, atoms, eigvals_only]) Returns the eigenvalues of the Hamiltonian
empty([keep]) See SparseCSR.empty for details
finalize() Finalizes the tight-binding model
iter([local]) Iterations of the orbital space in the geometry, two indices from loop
iter_nnz([atom, orbital]) Iterations of the non-zero elements, returns a tuple of orbital and coupling orbital
read(sile, *args, **kwargs) Reads Hamiltonian from Sile using read_H.
repeat(reps, axis) Refer to tile instead
reset([nnzpr, orthogonal, spin, dtype]) The sparsity pattern is cleaned and every thing is reset.
sp2HS(geom, H[, S]) Returns a tight-binding model from a preset H, S and Geometry
tile(reps, axis) Returns a repeated tight-binding model for this, much like the Geometry
tocsr(index[, isc]) Return a scipy.sparse.csr_matrix from the specified index
write(sile, *args, **kwargs) Writes a tight-binding model to the Sile as implemented in the Sile.write_es method
H
S
construct(func, na_iR=1000, method='rand', eta=False)[source]

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

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 dR the radii parameters for the corresponding tight-binding parameters. The second is the tight-binding parameters corresponding to the dR[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 the Hamiltonian object it-self (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.geom.close(ia, dR=[0.1, 1.44], idx=idxs, idx_xyz=idxs_xyz)
>>>     self.H[ia, idx[0]] = 0.   # on-site
>>>     self.H[ia, idx[1]] = -2.7 # nearest-neighbour

na_iR : int, 1000

number of atoms within the sphere for speeding up the iter_block loop.

method : str, ‘rand’

method used in Geometry.iter_block, see there for details

eta: bool, False

whether an ETA will be printed

copy(dtype=None)[source]

Return a copy of the Hamiltonian object

create_construct(dR, param)[source]

Returns a simple function for passing to the construct function.

This is simply to leviate the creation of simplistic functions needed for setting up the Hamiltonian.

Basically this returns a function: >>> def func(self, ia, idxs, idxs_xyz=None): >>> idx = self.geom.close(ia, dR=dR, idx=idxs) >>> for ix, p in zip(idx, param): >>> self[ia, ix] = p

Parameters:

dR : array_like

radii parameters for tight-binding parameters. Must have same length as param or one less. If one less it will be extended with dR[0]/100

param : array_like

coupling constants corresponding to the dR ranges. param[0,:] are the tight-binding parameter for the all atoms within dR[0] of each atom.

cut(seps, axis, *args, **kwargs)[source]

Cuts the tight-binding model into different parts.

Creates a tight-binding model by retaining the parameters for the cut-out region, possibly creating a super-cell.

Parameters:

seps : integer, optional

number of times the structure will be cut.

axis : integer

the axis that will be cut

dtype

Return data type of Hamiltonian (and overlap matrix)

eigh(k=(0, 0, 0), atoms=None, eigvals_only=True, overwrite_a=True, overwrite_b=True, *args, **kwargs)[source]

Returns the eigenvalues of the Hamiltonian

Setup the Hamiltonian and overlap matrix with respect to the given k-point, then reduce the space to the specified atoms and calculate the eigenvalues.

All subsequent arguments gets passed directly to scipy.linalg.eigh

eigsh(k=(0, 0, 0), n=10, atoms=None, eigvals_only=True, *args, **kwargs)[source]

Returns the eigenvalues of the Hamiltonian

Setup the Hamiltonian and overlap matrix with respect to the given k-point, then reduce the space to the specified atoms and calculate the eigenvalues.

All subsequent arguments gets passed directly to scipy.linalg.eigh

empty(keep=False)[source]

See SparseCSR.empty for details

finalize()[source]

Finalizes the tight-binding model

Finalizes the tight-binding model so that no new sparse elements can be added.

Sparse elements can still be changed.

finalized

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

geom

Return the attached geometry

geometry

Return the attached geometry

iter(local=False)[source]

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

An iterator returning the current atomic index and the corresponding orbital index.

>>> for ia, io in self:

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

Parameters:

local : bool=False

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

iter_nnz(atom=None, orbital=None)[source]

Iterations of the non-zero elements, returns a tuple of orbital and coupling orbital

An iterator returning the current orbital index and the corresponding connected orbital where a non-zero is defined

>>> for io, jo in self.iter_nnz():

In the above case io and jo are orbitals such that:

>>> self.H[io,jo] 

returns the non-zero element of the Hamiltonian.

One may reduce the iterated space by either requesting a specific set of atoms, or orbitals, _not_ both simultaneously.

Parameters:

atom : int/array_like

iterate on couplings to the set of atoms (not compatible with orbital)

orbital : int/array_like

iterate on couplings to the set of orbitals (not compatible with atom)

Examples

Looping only on one or more atoms:

>>> for io, jo in self.iter_nnz(atom=[2, 3]):
>>>     # loop on all orbitals on atom 3 and 4 (0 indexing)
>>> for io, jo in self.iter_nnz(orbital=[2, 3]):
>>>     # loop on orbitals 3 and 4 (0 indexing)
nnz

Returns number of non-zero elements in the tight-binding model

no

Returns number of orbitals as used when the object was created

orthogonal

Return whether the Hamiltonian is orthogonal

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

Reads Hamiltonian from Sile using read_H.

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

repeat(reps, axis)[source]

Refer to tile instead

reset(nnzpr=None, orthogonal=True, spin=1, dtype=None)[source]

The sparsity pattern is cleaned and every thing is reset.

The object will be the same as if it had been initialized with the same geometry as it were created with.

Parameters:

nnzpr: int

number of non-zero elements per row

orthogonal: boolean, True

if there is an overlap matrix associated with the Hamiltonian

spin: int, 1

number of spin-components

dtype: ``numpy.dtype``, `numpy.float64`

the datatype of the Hamiltonian

classmethod sp2HS(geom, H, S=None)[source]

Returns a tight-binding model from a preset H, S and Geometry

spin

Return number of spin-components in Hamiltonian

tile(reps, axis)[source]

Returns a repeated tight-binding model for this, much like the Geometry

The already existing tight-binding parameters are extrapolated to the new supercell by repeating them in blocks like the coordinates.

Parameters:

reps : number of tiles (repetitions)

axis : direction of tiling

0, 1, 2 according to the cell-direction

tocsr(index, isc=None)[source]

Return a scipy.sparse.csr_matrix from the specified index

Parameters:

index : int

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

isc : int, None

the supercell index (or all)

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

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

sisl.quantity.hamiltonian.TightBinding

alias of Hamiltonian