Tight Binding

The user tight-binding model contains all relevant information regarding the orbital basis, the model Hamiltonian (in addition to eigenvalues/eigenvectors), as well as the momentum domain of interest. In addition, the tight-binding class contains all relevant methods required to extract this information.

The user has a reasonable amount of flexibility in the format which they would use to generate the model Hamiltonian. This is intended to accommodate various programs which can be used to generate a tight-binding Hamiltonian (for example Wannier90), as well as the different presentations used in publications (Slater-Koster \(V_{lmk}\) , \(t_{ij}\) ), as well as alternative Hamiltonian types, such as low-energy effective models which do not adhere to the full translational symmetry required by Bloch’s theorem. These latter models do however provide a highly useful and physically transparent parameterization of the system Hamiltonian for narrow regions of momentum space near high-symmetry points. For these reasons, there are 4 general categories of Hamiltonian inputs we accept. The first three are described as: Slater-Koster, list, and text input. The last is described generically as executable.

While in principal a tight-binding Hamiltonian can be passed in the acceptable form for any of the above, the last option also supports these low-energy theories described above.

Slater-Koster Models:

In their 1954 PRB paper, Slater and Koster presented a fairly simple framework for defining the tight-binding model associated with hoppings between localized atomic-like orbitals on neighbouring lattice sites. To define the overlap for a pair of basis states with orbital angular momentum \(l_1\) and \(l_2\) for an arbitrary lattice geometry, we require only \((min( l_1, l_2 ) + 1)\) parameters. For example, for \(l_1=l_2=1\) we define \(V_{pp\sigma},\ V_{pp\pi}\). Intuitively, these parameters correspond to overlap integrals between the two orbitals when the lobes of the ‘p’ orbitals are aligned parallel to the connecting vector ( \(\sigma\)) and aligned perpendicular to the connecting vector ( \(\pi\) ). One can often use the frequently published table of Slater-Koster parameters to then define general rules for how these two parameters should be combined for a specific lattice geometry to arrive at the hopping amplitude between each pair of lattice sites.

This table is however restrictive as it provides rules for only hoppings between non-distorted cubic harmonics. One can alternatively take advantage of the representation of the orbital states in the basis of spherical harmonics \(Y_{l}^{m}(\Omega)\) to rotate an arbitrary pair of basis states into a common representation, and then rotate the frame of reference to align the bond-direction with a designated quantization axis: here the \(\hat z\) vector. A diagonal Hamiltonian matrix filled with the associated \(V_{l_1l_2\gamma}\) can then be applied. The rotation and basis transformation can be undone, ultimately producing a matrix of Hamiltonian parameters for the entire orbital shell along the designated bond vector. Mathematically, the procedure is represented by

\[\left<\psi|H|\phi\right> = \left<\psi| U^{\dagger} R^{-1}(\theta,\phi) V_{SK} R(\theta,\phi) U |\phi\right>\]

This formalism then allows for fast and user-friendly generation of a Hamiltonian over an arbitrary basis and geometry. Given only the \(V_{l_1,l_2,\gamma}\) parameters and the lattice geometry, a full tight-binding Hamiltonian can be built.

The Slater-Koster parameters are passed in the form of a dictionary, with the keys taking the form of \(a_1 a_2 n_1 n_2 l_1 l_2 \gamma\). For example, if I have two distinct atoms in my basis, where I have the Carbon 2p and Iron 3d \(e_g\) orbitals, the dictionary may look like shown here

V_SK = {'021':3.0,'13XY':0.0,'13ZR':-0.5,
'002211S':1.0,'002211P':-0.5,
'012312S':0.1,'012312P':0.6,
'113322S':0.05,'113322P':0.7,'113322D':-0.03}

In the first line I have written the on-site terms. While in a simple case I may have a single onsite energy for each orbital shell, here I have distinguished the onsite energy of the \(d_{x^2-y^2}\) and \(d_{3z^2-r^2}\) states. In the second-fourth lines, I have written the p-p, p-d, and d-d hopping amplitudes.

In many models, hopping beyond nearest neighbours are relevant, and will typically not have the same strength as the nearest neighbour terms. In these cases, we can pass a list of dictionaries to H_dict. For example

H_dict = {'type':'SK',
                'V':[SK_dict1,SK_dict2],
                'cutoff':[3.5,5.0],
                ...}

To include these next-nearest neighbour terms, I specify a list of hopping dictionaries, in addition to a list of cutoff distances, indicating the range of connecting vectors each of the indicated dictionaries should apply to. For this case, for connecting vectors where \(|R_{ij}|<3.5\) we use SK_dict1, whereas for \(3.5\leq |R_{ij}|<5.0\) we use SK_dict2.

Hamiltonian Construction:

From these matrices, we proceed to build a list of H_me objects. This class is the standard representation of the Hamiltonian in chinook. Each instance of H_me carries integers labelling the associated basis indices. For example \(\left<\phi_3|H|\phi_7\right>\) will be stored with H_me.i = 3, H_me.j = 7. We note here that the Hermiticity of the Hamiltonian allows one to explicitly define only the upper or lower diagonal of the Hamiltonian. Consequently, in chinook, we use only the upper diagonal, that is \(i\leq j\).

In addition to basis indices, the H_me object carries the functional information of the matrix element in its H_me.H attribute. This will be a list of connecting vectors and hopping strengths in the case of a standard tight-binding model, or a list of python executables otherwise.

In standard format then one might have

TB.mat_els[9].i = 3
TB.mat_els[9].j = 3
TB.mat_els[9].H = [[0.0,0.0,0.0,2.0],
                [1.5,0.0,0.0,3.7],
                [-1.5,0.0,0.0,3.7]]

This is the 10-th element in our model’s list TB.mat_els of H_me objects. This H_me instance contains an on-site term of strength 2.0, and a cosine-like hopping term along the \(x\) -axis of strength 3.7 eV. A closer consideration of the H attribute reveals the essential information. Each element of the list is a length-4 list of float, containing \(\vec{R_{ij}^n}\) and \(t_{ij}^n\). Ultimately, the full matrix element will be expressed as

\[H_{ij}(\vec{k}) = \sum_{n} t_{ij}^n e^{i\vec{k}\cdot\vec{R_{ij}^n}}\]

For this reason, if one does not have access to a suitable Slater-Koster representation of the tight-binding model, we can bypass the methods described above, passing a list of prepared matrix elements directly to the generator of the H_me objects. To accommodate this, we then also accept Hamiltonian input in the form of list, where each element of the list is written as for example

Hlist[10] = [3,3,1.5,0.0,0.0,3.7]

The elements then correspond to basis index 1, basis index 2, connecting vector x, y, and z components, and finally the \(t_{ij}\) value. Similarly, a textfile of this form is also acceptable, with each line in the textfile being a comma-separated list of values

3,3,1.5,0.0,0.0,3.7
...

A list of H_me objects can then be built as above.

Executable Hamiltonians:

The executable type Hamiltonian will ultimately require a slightly modified form of the input, as we do not intend to express the matrix elements as a Fourier transform of the real-space Hamiltonian in the same way as above.

In this form, we should have

\[H_{ij}(\vec{k}) = \sum_{n} f_n(\vec{k})\]

This then requires a modified form of H_me.H. By contrast with the above,

TB.mat_els[9].H = [np.cos,my_hopping_func]

where my_hopping_func is a user-defined executable. It could be for example a 2nd order polynomial, or any other function of momentum. The essential point is that the executables must take an Nx3 array of float as their only argument. This allows for these Hamiltonians to fit seamlessly into the chinook framework.

Warning

Automated tools for building these Hamiltonians are currently in the process of being built. Proceed with caution, and notify the developers of any unexpected performance.

For advanced users who would like to take advantage of this functionality now,

H_dict = {'type':'exec',
                'exec':exec_list,
                ...}

where the item

exec_list = [[(0,0),my_func1],
[(0,1),my_func2],
...]

The tuple of integers preceding the name of the executable corresponds to the basis indices. Here my_func1 applies to the diagonal elements of my first basis member, and my_func2 to the coupling of my first and second basis elements. As always, indexing in python is 0-based.

For the construction of the executable functions, we recommend the use of python’s lambda functions. For example, to define a function which evaluates to \(\alpha k_x^2 + \beta k_y^2\) , one may define

def my_func_generator(alpha,beta):
        return lambda kvec: alpha*k[:,0]**2 + beta*k[:,1]**2

#Define specific parameters for my executable functions
a1 = 2.3
b1 = 0.7

a2 = 5.6
b2 = 0.9

#Define my executables
my_func1 = my_func_generator(a1,b1)
my_func2 = my_func_generator(a2,b2)

#Populate a list of executables
exec_list = [[(0,0),my_func1],
[(0,1),my_func2]]

Please stay tuned for further developments which will facilitate more convenient construction of these Hamiltonians.

Below, we include the relevant documentation to the construction of tight-binding models.

Model Initialization

build_lib.gen_K(Kdic)

Generate k-path for TB model to be diagonalized along.

args:

  • Kdic: dictionary for generation of kpath with:

    • ‘type’: string ‘A’ (absolute) or ‘F’ (fractional) units
    • ‘avec’: numpy array of 3x3 float lattice vectors
    • ‘pts’: list of len3 array indicating the high-symmetry points

    along the path of interest

    • ‘grain’: int, number of points between each element of ‘pts’

    optional:

    • ‘labels’:list of strings with same length as ‘pts’, giving

    plotting labels for the kpath

return:

Kobj: K-object including necessary attributes to be read by the TB_model
build_lib.gen_TB(basis_dict, hamiltonian_dict, Kobj=None, slab_dict=None)

Build a Tight-Binding Model using the user-input dictionaries

args:

  • basis_dict: dictionary, including the ‘bulk’ key value pair

generated by gen_basis

  • hamiltonian_dict: dictionary,

    • ‘spin’: same dictionary as passed to gen_basis
    • ‘type’: string, Hamiltonian type–‘list’ (list of matrix elements),

    ‘SK’ (Slater-Koster dictionaries, requires also a ‘V’ and ‘avec’ entry), ‘txt’ (textfile, requires a ‘filename’ key as well)

    • ‘cutoff’: float, cutoff hopping distance
    • ‘renorm’: optional float, renormalization factor default to 1.0
    • ‘offset’: optional float, offset of chemical potential, default to 0.0
    • ‘tol’: optional float, minimum matrix element tolerance, default to 1e-15
  • Kobj: optional, standard K-object, as generated by gen_K

  • slab_dict: dictionary for slab generation

    • ‘avec’: numpy array of 3x3 float, lattice vectors
    • ‘miller’: numpy array of 3 integers, indicating the Miller

    index of the surface normal in units of lattice vectors

    • ‘fine’: fine adjustment of the slab thickness, tuple of two

    numeric to get desired termination correct (for e.g. inversion symmetry)

    • ‘thick’: integer approximate number of unit cells in the

    slab (will not be exact, depending on the fine, and termination

    • ‘vac’: int size of the vacuum buffer – must be larger than

    the largest hopping length to ensure no coupling of slabs

    • ‘termination’: tuple of 2 integers: atom indices which

    terminate the top and bottom of the slab

return:

TB_model: tight-binding object, as defined in chinook.TB_lib.py
build_lib.gen_basis(basis)

Generate a list of orbital objects as the input basis for a tight-binding model. User passes a basis dictionary, function returns a modified version of this same dictionary, with the list of orbitals now appended as the ‘bulk’ entry

args:

  • basis–dictionary with keys:

    • ‘atoms’: list of integer, indices for distinct atoms,
    • ‘Z’: dictionary of integer: ‘atom’:element (integer) pairs
    • ‘orbs’: list of lists of string, for each atom containing the

    orbital labels (usually in conventional nlxx format)),

    • ‘pos’: list of numpy arrays of length 3 float indicating

    positions of the atoms in direct Angstrom units,

    • optional keys:

      • ‘orient’: list, one entry for each atom, indicating a

      local rotation of the indicated atom, various formats accepted; for more details, c.f. chinook.orbital.py

      • ‘spin’: dictionary of spin information:

        • ‘bool’: boolean, double basis into spinor basis,
        • ‘soc’: boolean, include spin-orbit coupling
        • ‘lam’: dictionary of SOC constants, integer:float

        pairs for atoms in ‘atoms’ list, and lambda_SOC in eV

return:

  • basis dictionary, modified to include the bulk list of orbital

objects

***

build_lib.recur_product(elements)

Utility function: Recursive evaluation of the product of all elements in a list

args:

  • elements: list of numeric type

return:

  • product of all elements of elements

***

Hamiltonian Library

H_library.AFM_order(basis, dS, p_up, p_dn)

Add antiferromagnetism to the tight-binding model, by adding a different on-site energy to orbitals of different spin character, on the designated sites.

args:

  • basis: list, orbital objects
  • dS: float, size of spin-splitting (eV)
  • p_up, p_dn: numpy array of float indicating the orbital positions

for the AFM order

return:

  • h_AF: list of matrix elements, as conventionally arranged [[o1,o2,0,0,0,H12],…]

***

H_library.FM_order(basis, dS)

Add ferromagnetism to the system. Take dS to assume that the splitting puts spin-up lower in energy by dS,and viceversa for spin-down. This directly modifies the TB_model’s mat_els attribute

args:

  • basis: list, of orbital objects in basis
  • dS: float, energy of the spin splitting (eV)

return:

  • list of matrix elements [[o1,o2,0,0,0,H12],…]

***

H_library.Lm(l)

L- operator in the l,m_l basis, organized with (0,0) = |l,l>… (2l,2l) = |l,-l>

The nonzero elements are on the upper diagonal

arg:

  • l: int orbital angular momentum

return:

  • M: numpy array (2l+1,2l+1) of real float

***

H_library.Lp(l)

L+ operator in the \(l\), \(m_l\) basis, organized with (0,0) = |l,l>… (2l,2l) = |l,-l> The nonzero elements are on the upper diagonal

arg:

  • l: int orbital angular momentum

return:

  • M: numpy array (2l+1,2l+1) of real float

***

H_library.Lz(l)

Lz operator in the l,:math:m_l basis

arg:

  • l: int orbital angular momentum

return:

  • numpy array (2*l+1,2*l+1)

***

H_library.SO(basis)

Generate L.S matrix-elements for a given basis. This is generic to all l, except the normal_order, which is defined here up to and including the f electrons. Otherwise, this method is generic to any orbital angular momentum.

In the factors dictionary defined here indicates the weight of the different \(L_iS_i\) terms. The keys are tuples of (L+/-/z,S+/-/z) in a bit of a cryptic way: for L, (0,1,2) ->(-1,0,1) and for S, (-1,0,1) = S1-S2 with S1,2 = +/- 1 here

L+,L-,Lz matrices are defined for each l shell in the basis, transformed into the basis of the tight-binding model. The nonzero terms will then just be used along with the spin and weighted by the factor value, and slotted into a len(basis)xlen(basis) matrix HSO

args:

  • basis: list of orbital objects

return:

  • HSO: list of matrix elements in standard format [o1,o2,0,0,0,H12]

***

H_library.Vlist_gen(V, pair)

Select the relevant hopping matrix elements to be used in defining the value of the Slater-Koster matrix elements for a given pair of orbitals. Handles situation where insufficient parameters have been passed to system.

args:

  • V: dictionary of Slater-Koster hopping terms
  • pair: tuple of int defining the orbitals to be paired, (a1,a2,n1,n2,l1,l2)

return:

  • Vvals: numpy array of Vllx related to a given pairing, e.g. for s-p np.array([Vsps,Vspp])

***

H_library.cluster_init(Vdict, cutoff, avec)

Generate a cluster of neighbouring lattice points to use in defining the hopping paths–ensuring that it extends sufficiently far enough to capture even the largest hopping vectors. Also reforms the SK dictionary and cutoff lengths to be in list format. Returns an array of lattice points which go safely to the edge of the cutoff range.

args:

  • Vdict: dictionary, or list of dictionaries of Slater Koster matrix elements
  • cutoff: float, or list of float
  • avec: numpy array of 3x3 float

return:

  • Vdict: list of length 1 if a single dictionary passed, else unmodified
  • cutoff: numpy array, append 0 to the beginning of the cutoff list,

else leave it alone.

  • pts: numpy array of lattice vector indices for a region of lattice points around

the origin.

***

H_library.index_ordering(basis)

We use an universal ordering convention for defining the Slater-Koster matrices which may (and most likely will) not match the ordering chosen by the user. To account for this, we define a dictionary which gives the ordering, relative to the normal order convention defined here, associated with a given a-n-l shell at each site in the lattice basis.

args:

  • basis: list of orbital objects

return:

  • indexing: dictionary of key-value pairs (a,n,l,x,y,z):numpy.array([…])

***

H_library.mat_els(Rij, SKmat, tol, i1, i2)

Extract the pertinent, and non-zero elements of the Slater-Koster matrix and transform to the conventional form of Hamiltonian list entries (o1,o2,Rij0,Rij1,Rij2,H12(Rij))

args:

  • Rij: numpy array of 3 float, relevant connecting vector
  • SKmat: numpy array of float, matrix of hopping elements

for the coupling of two orbital shells

  • tol: float, minimum hopping included in model
  • i1, i2: int,int, proper index ordering for the relevant

instance of the orbital shells involved in hopping

return:

  • out: list of Hamiltonian matrix elements, extracted from the

ordered SKmat, in form [[o1,o2,x12,y12,z12,H12],…]

***

H_library.mirror_SK(SK_in)

Generate a list of values which is the input appended with its mirror reflection. The mirror boundary condition suppresses the duplicate of the last value. e.g. [0,1,2,3,4] –> [0,1,2,3,4,3,2,1,0], [‘r’,’a’,’c’,’e’,’c’,’a’,’r’] –> [‘r’,’a’,’c’,’e’,’c’,’a’,’r’,’a’,’c’,’e’,’c’,’a’,’r’] Intended here to take an array of Slater-Koster hopping terms and reflect about its last entry i.e. [Vsps,Vspp] -> [Vsps,Vspp,Vsps]

args:

  • SK_in: iterable, of arbitrary length and data-type

return:

  • list of values with same data-type as input

***

H_library.on_site(basis, V, offset)

On-site matrix element calculation. Try both anl and alabel formats, if neither is defined, default the onsite energy to 0.0 eV

args:

  • basis: list of orbitals defining the tight-binding basis
  • V: dictionary, Slater Koster terms
  • offset: float, EF shift

return:

  • Ho: list of Hamiltonian matrix elements

***

H_library.region(num)

Generate a symmetric grid of points in number of lattice vectors.

args:

  • num: int, grid will have size 2*num+1 in each direction

return:

  • numpy array of size ((2*num+1)**3,3) with centre value of first entry

of (-num,-num,-num),…,(0,0,0),…,(num,num,num)

***

H_library.sk_build(avec, basis, Vdict, cutoff, tol, renorm, offset)

Build SK model from using D-matrices, rather than a list of SK terms from table. This can handle orbitals of arbitrary orbital angular momentum in principal, but right now implemented for up to and including f-electrons. NOTE: f-hoppings require thorough testing

args:

  • avec: numpy array 3x3 float, lattice vectors
  • basis: list of orbital objects
  • Vdict: dictionary, or list of dictionaries, of Slater-Koster integrals/ on-site energies
  • cutoff: float or list of float, indicating range where Vdict is applicable
  • tol: float, threshold value below which hoppings are neglected
  • offset: float, offset for Fermi level

return:

  • H_raw: list of Hamiltonian matrix elements, in form [o1,o2,x12,y12,z12,t12]

***

H_library.spin_double(H, lb)

Duplicate the kinetic Hamiltonian terms to extend over the spin-duplicated orbitals, which are by construction in same order and appended to end of the original basis.

args:

  • H: list, Hamiltonian matrix elements [[o1,o2,x,y,z,H12],…]
  • lb: int, length of basis before spin duplication

return:

  • h2 modified copy of H, filled with kinetic terms for both

spin species

***

H_library.txt_build(filename, cutoff, renorm, offset, tol, Nonsite)

Build Hamiltonian from textfile, input is of form o1,o2,x12,y12,z12,t12, output in form [o1,o2,x12,y12,z12,t12]. To be explicit, each row of the textfile is used to generate a k-space Hamiltonian matrix element of the form:

\[H_{1,2}(k) = t_{1,2} e^{i (k_x x_{1,2} + k_y y_{1,2} + k_z z_{1,2})}\]

args:

  • filename: string, name of file
  • cutoff: float, maximum distance of hopping allowed, Angstrom
  • renorm: float, renormalization of the bandstructure
  • offset: float, energy offset of chemical potential, electron volts
  • tol: float, minimum Hamiltonian matrix element amplitude
  • Nonsite: int, number of basis states, use to apply offset

return:

  • Hlist: the list of Hamiltonian matrix elements

***

Momentum Library

klib.b_zone(a_vec, N, show=False)

Generate a mesh of points over the Brillouin zone. Each of the cardinal axes are divided by the same number of points (so points are not necessarily evenly spaced along each axis).

args:
  • a_vec: numpy array of size 3x3 float
  • N: int mesh density
kwargs:
  • show: boolean for optional plotting of the mesh points
return:
  • m_pts: numpy array of mesh points (float), shape (len(m_pts),3)
klib.bvectors(a_vec)

Define the reciprocal lattice vectors corresponding to the direct lattice in real space

args:
  • a_vec: numpy array of 3x3 float, lattice vectors
return:
  • b_vec: numpy array of 3x3 float, reciprocal lattice vectors

***

klib.kmesh(ang, X, Y, kz, Vo=None, hv=None, W=None)

Take a mesh of kx and ky with fixed kz and generate a Nx3 array of points which rotates the mesh about the z axis by ang. N is the flattened shape of X and Y.

args:
  • ang: float, angle of rotation
  • X: numpy array of float, one coordinate of meshgrid
  • Y: numpy array of float, second coordinate of meshgrid
  • kz: float, third dimension of momentum path, fixed
kwargs:
  • Vo: float, parameter necessary for inclusion of inner potential
  • hv: float, photon energy, to be used if Vo also included, for evaluating kz
  • W: float, work function
return:
  • k_arr: numpy array of shape Nx3 float, rotated kpoint array.
  • ph: numpy array of N float, angles of the in-plane momentum

points, before rotation.

***

klib.kmesh_hv(ang, x, y, hv, Vo=None)

Take a mesh of kx and ky with fixed kz and generate a Nx3 array of points which rotates the mesh about the z axis by ang. N is the flattened shape of X and Y.

args:
  • ang: float, angle of rotation
  • X: numpy array of float, one coordinate of meshgrid
  • Y: numpy array of float, second coordinate of meshgrid
  • kz: float, third dimension of momentum path, fixed
kwargs:
  • Vo: float, parameter necessary for inclusion of inner potential
  • hv: float, photon energy, to be used if Vo also included, for evaluating kz
  • W: float, work function
return:
  • k_arr: numpy array of shape Nx3 float, rotated kpoint array.
  • ph: numpy array of N float, angles of the in-plane momentum

points, before rotation.

***

class klib.kpath(pts, grain=None, labels=None)

Momentum object, defining a path in reciprocal space, for use in defining the Hamiltonian at different points in the Brillouin zone. ***

points()

Use the endpoints of kpath defined in kpath.pts to create numpy array of len(3) float which cover the entire path, based on method by I.S. Elfimov.

return:
  • kpath.kpts: numpy array of float, len(kpath.pts)(1+**kpath.grain**) by 3

***

klib.kz_kpt(hv, kpt, W=0, V=0)

Extract the kz associated with a given in-plane momentum, photon energy, work function and inner potential

args:
  • hv: float, photon energ in eV
  • kpt: float, in plane momentum, inverse Angstrom
  • W: float, work function in eV
  • V: float, inner potential in eV
return:
  • kz: float, out of plane momentum, inverse Angstrom

***

klib.mesh_reduce(blatt, mesh, inds=False)

Determine and select only k-points corresponding to the first Brillouin zone, by simply classifying points on the basis of whether or not the closest lattice point is the origin. By construction, the origin is index 13 of the blatt. If it is not, return error. Option to take only the indices of the mesh which we want, rather than the actual array points –this is relevant for tetrahedral interpolation methods

args:
  • blatt: numpy array of len(27,3), nearest reciprocal lattice vector points
  • mesh: numpy array of (N,3) float, defining a mesh of k points, before

being reduced to contain only the Brillouin zone.

kwargs:
  • inds: option to pass a list of bool, indicating the

indices one wants to keep, instead of autogenerating the mesh

return:
  • bz_pts: numpy array of (M,3) float, Brillouin zone points
klib.plt_pts(pts)

Plot an array of points iin 3D

args:
  • pts: numpy array shape N x 3
klib.raw_mesh(blatt, N)

Define a mesh of points filling the region of k-space bounded by the set of reciprocal lattice points generated by bvectors. These will be further reduced by mesh_reduce to find points which are within the first-Brillouin zone

args:
  • blatt: numpy array of 27x3 float
  • N: int, or iterable of len 3, defines a coarse estimation

of number of k-points

return:
  • mesh: numpy array of mesh points, size set roughly by N

***

klib.region(num)

Generate a symmetric grid of points in number of lattice vectors.

args:
  • num: int, grid will have size 2 num+1 in each direction
return:
  • numpy array of size ((2 num+1)**3,3) with centre value of

first entry of (-num,-num,-num),…,(0,0,0),…,(num,num,num)

***

Orbital Objects

orbital.fact(n)

Recursive factorial function, works for any non-negative integer.

args:

  • n: int, or integer-float

return:

  • int, recursively evaluates the factorial of the initial input value.

***

class orbital.orbital(atom, index, label, pos, Z, orient=[0.0], spin=1, lam=0.0, sigma=1.0, slab_index=None)

The orbital object carries all essential details of the elements of the model Hamiltonian basis, for both generation of the tight-binding model, in addition to the evaluation of expectation values and ARPES intensity.

copy()

Copy by value method for orbital object

return:

  • orbital_copy: duplicate of orbital object

***

orbital.rot_projection(l, proj, rotation)

Go through a projection array, and apply the intended transformation to the Ylm projections in order. Define Euler angles in the z-y-z convention THIS WILL BE A COUNTERCLOCKWISE ROTATION ABOUT a vector BY angle gamma expressed in radians. Note that we always define spin in the lab-frame, so spin degrees of freedom are not rotated when we rotate the orbital degrees of freedom.

args:

  • l: int,orbital angular momentum
  • proj: numpy array of shape Nx4 of float, each element is

[Re(projection),Im(projection),l,m]

  • rotation: float, or list, defining rotation of orbital. If float,

assume rotation about z-axis. If list, first element is a numpy array of len 3, indicating rotation vector, and second element is float, angle.

return:

  • proj: numpy array of Mx4 float, as above, but modified, and may

now include additional, or fewer elements than input proj.

  • Dmat: numpy array of (2l+1)x(2l+1) complex float indicating the

Wigner Big-D matrix associated with the rotation of this orbital shell about the intended axis.

***

orbital.slab_basis_copy(basis, new_posns, new_inds)

Copy elements of a slab basis into a new list of orbitals, with modified positions and index ordering.

args:

  • basis: list or orbital objects
  • new_posns: numpy array of len(basis)x3 float, new positions for

orbital

  • new_inds: numpy array of len(basis) int, new indices for orbitals

return:

  • new_basis: list of duplicated orbitals following modification.

***

orbital.sort_basis(basis, slab)

Utility script for organizing an orbital basis that is out of sequence

args:

  • basis: list of orbital objects
  • slab: bool, True or False if this is for sorting a slab

return:

  • orb_basis: list of sorted orbital objects (by orbital.index value)

***

orbital.spin_double(basis, lamdict)

Double the size of a basis to introduce spin to the problem. Go through the basis and create an identical copy with opposite spin and incremented index such that the orbital basis order puts spin down in the first N/2 orbitals, and spin up in the second N/2.

args:

  • basis: list of orbital objects
  • lamdict: dictionary of int:float pairs providing the

spin-orbit coupling strength for the different inequivalent atoms in basis.

return:

  • doubled basis carrying all required spin information

***

Slater Koster Library

SlaterKoster.SK_cub(Ymats, l1, l2)

In order to generate a set of independent Lambda functions for rapid generation of Hamiltonian matrix elements, one must nest the definition of the lambda functions within another function. In this way, we avoid cross-contamination of unrelated functions. The variables which are fixed for a given lambda function are the cubic -to- spherical harmonics (Ymat) transformations, and the orbital angular momentum of the relevant basis channels. The output lambda functions will be functions of the Euler-angles pertaining to the hopping path, as well as the potential matrix V, which will be passed as a numpy array (min(l1,l2)*2+1) long of float.

We follow the method described for rotated d-orbitals in the thesis of JM Carter from Toronto (HY Kee), where the Slater-Koster hopping matrix can be defined as the following operation:

  1. Transform local orbital basis into spherical harmonics
  2. Rotate the hopping path along the z-axis
  3. Product with the diagonal SK-matrix
  4. Rotate the path backwards
  5. Rotate back into basis of local orbitals

6. Output matrix of hopping elements between all orbitals in the shell to fill Hamiltonian

args:

  • Ymats: list of numpy arrays corresponding to the relevant

transformation from cubic to spherical harmonic basis

  • l1, l2: int orbital angular momentum channels relevant

to a given hopping pair

return:

  • lambda function for the SK-matrix between these orbital shells,

for arbitrary hopping strength and direction.

***

SlaterKoster.SK_full(basis)

Generate a dictionary of lambda functions which take as keys the atom,orbital for both first and second element. Formatting is a1a2n1n2l1l2, same as for SK dictionary entries

args:

  • basis: list of orbital objects composing the TB-basis

return:

  • SK_funcs: a dictionary of hopping matrix functions

(lambda functions with args EA,EB,Ey,V as Euler angles and potential (V)) which can be executed for various hopping paths and potential strengths The keys of the dictionary will be organized similar to the way the SK parameters are passed, labelled by a1a2n1n2l1l2, which completely defines a given orbital-orbital coupling

***

SlaterKoster.Vmat(l1, l2, V)

For Slater-Koster matrix element generation, a potential matrix is sandwiched in between the two bond-rotating Dmatrices. It should be of the shape 2*l1+1 x 2*l2+1, and have the V_l,l’,D terms along the ‘diagonal’– a concept that is only well defined for a square matrix. For mismatched angular momentum channels, this turns into a diagonal square matrix of dimension min(2*l1+1,2*l2+1) centred along the larger axis. For channels where the orbital angular momentum change involves a change in parity, the potential should change sign, as per Slater Koster’s original definition from 1954. This is taken care of automatically in the Wigner formalism I use here, no need to have exceptions

args:

  • l1, l2: int orbital angular momentum of initial and final states
  • V: numpy array of float – length should be min(l1 ,**l2**)*2+1

return:

  • Vm: numpy array of float, shape 2 l1 +1 x 2 l2 +1

***

Tight-Binding Library

class TB_lib.H_me(i, j, executable=False)

This class contains the relevant executables and data structure pertaining to generation of the Hamiltonian matrix elements for a single set of coupled basis orbitals. Its attributes include integer values i, j indicating the basis indices, and a list of hopping vectors/matrix element values for the Hamiltonian.

The method H2Hk provides an executable function of momentum to allow broadcasting of the Hamiltonian over a large array of momenta. Python’s flexible protocol for equivalency and passing variables by reference/value require definition of a copy operator which allows one to produce safely, a copy of the object rather than its coordinates in memory alone. ***

H2Hk()

Transform the list of hopping elements into a Fourier-series expansion of the Hamiltonian. This is run during diagonalization for each matrix element index. If running a low-energy Hamiltonian, executable functions are simply summed for each basis index i,j, rather than computing a Fourier series. x is implicitly a numpy array of Nx3: it is essential that the executable conform to this input type.

return:

  • lambda function of a numpy array of float of length 3

***

append_H(H, R0=0, R1=0, R2=0)

Add a new hopping path to the coupling of the parent orbitals.

args:

  • H: complex float, matrix element strength, or if self.exectype,

should be an executable

  • R0, R1, R2: float connecting vector in cartesian

coordinate frame–this is the TOTAL vector, not the relevant lattice vectors only

return:

  • directly modifies the Hamiltonian list for these matrix

coordinates

***

clean_H()

Remove all duplicate instances of hopping elements in the matrix element list. This function is run automatically during slab generation.

The Hamiltonian list is not itself directly modified.

return:

  • list of hopping vectors and associated Hamiltonian matrix

element strengths

***

copy()

Copy by value of the H_me object

return:

  • H_copy: duplicate H_me object

***

class TB_lib.TB_model(basis, H_args, Kobj=None)

The TB_model object carries the model basis as a list of orbital objects, as well as the model Hamiltonian, as a list of H_me. The orbital indices paired in each instance of H_me are stored in a dictionary under ijpairs

append_H(new_elements)

Add new terms to the Hamiltonian by hand. This directly modifies the list of Hamiltonian matrix element, self.mat_els of the TB object.

args:

  • new_elements: list of Hamiltonian matrix elements, either a single element [i,j,x_ij,y_ij,z_ij,H_ij(x,y,z)]

or as a list of such lists. Here i, j are the related orbital-indices.

***

build_ham(H_args)

Buld the Hamiltonian using functions from chinook.H_library.py

args:

  • H_args: dictionary, containing all relevant information for

defining the Hamiltonian list. For details, see TB_model.__init__.

return:

  • sorted list of matrix element objects. These objects have

i,j attributes referencing the orbital basis indices, and a list of form [R0,R1,R2,Re(H)+1.0jIm(H)]

***

copy()

Copy by value of the TB_model object

return:

  • TB_copy: duplicate of the TB_model object.
plot_unitcell(ax=None)

Utility script for visualizing the lattice and orbital basis. Distinct atoms are drawn in different colours

kwargs:

  • ax: matplotlib Axes, for plotting on existing Axes

return:

  • ax: matplotlib Axes, for further modifications to plot

***

plotting(win_min=None, win_max=None, ax=None)

Plotting routine for a tight-binding model evaluated over some path in k. If the model has not yet been diagonalized, it is done automatically before proceeding.

kwargs:

  • win_min, win_max: float, vertical axis limits for plotting

in units of eV. If not passed, a reasonable choice is made which covers the entire eigenspectrum.

  • ax: matplotlib Axes, for plotting on existing Axes

return:

  • ax: matplotlib axes object

***

print_basis_summary()

Very basic print function for printing a summary of the orbital basis, including their label, atomic species, position and spin character.

solve_H(Eonly=False)

This function diagonalizes the Hamiltonian over an array of momentum vectors. It uses the mat_el objects to quickly define lambda functions of momentum, which are then filled into the array and diagonalized. According to https://ntrs.nasa.gov/archive/nasa/casi.ntrs.nasa.gov/20050192421.pdf SVD algorithms require memory of 2*order*(4*order + 1) ~ 8*order^2. The matrices are complex float, so this should be 16 bytes per entry: so len(k)*(2048*order**2). If the diagonalization is requesting more than 85% of the available memory, then split up the k-path into sequential diagonalizations.

return:

  • self.Eband: numpy array of float, shape(len(self.Kobj.kpts),len(self.basis)),

eigenvalues

  • self.Evec: numpy array of complex float, shape(len(self.Kobj.kpts),len(self.basis),len(self.basis))

eigenvectors

***

unpack()

Reduce a Hamiltonian object down to a list of matrix elements. Include the Hermitian conjugate terms

return:

  • Hlist: list of Hamiltonian matrix elements

***

TB_lib.atom_coords(basis)

Define a dictionary organizing the distinct coordinates of instances of each atomic species in the basis

args:

  • basis: list of orbital objects

return:

  • **dictionary with integer keys, numpy array of float values. atom:locations are

encoded in this way

***.

TB_lib.cell_edges(avec)

Define set of line segments which enclose the unit cell.

args:

  • avec: numpy array of 3x3 float

return:

  • edges: numpy array of 12 x 6, endpoints of the 12 edges of the unit cell parallelepiped

***

TB_lib.gen_H_obj(htmp, executable=False)

Take a list of Hamiltonian matrix elements in list format: [i,j,Rij[0],Rij[1],Rij[2],Hij(R)] and generate a list of H_me objects instead. This collects all related matrix elements for a given orbital-pair for convenient generation of the matrix Hamiltonians over an input array of momentum

args:

  • htmp: list of numeric-type values (mixed integer[:2], float[2:5], complex-float[-1])

kwargs:

  • executable: boolean, if True, we don’t have a standard Fourier-type Hamiltonian,

but perhaps a low-energy expansion. In this case, the htmp elements are

return:

  • Hlist: list of Hamiltonian matrix element, H_me objects

***