Internal API

This page lists unexported functionality from SymmetricTightBinding.jl that may be of interest to developers.


Unexported, internal functionality

SymmetricTightBinding.OrbitalOrderingMethod
OrbitalOrdering(br::NewBandRep{D}) --> OrbitalOrdering{D}

Establishes a canonical, local ordering for the orbitals associated to a band representation br. This is the default ordering used when associating row/column indices in a tight-binding Hamiltonian block to specific orbitals in the associated band representations.

The canonical orbital ordering is stored in .ordering. The ith elements of ordering, ordering[i], is a NamedTuple with two fields: wp and idx:

  • wp: stores a Wyckoff position in the orbit of the Wyckoff position associated to br.wp.
  • idx: stores the index of the partner function of the site-symmetry irrep associated to br at wp.

I.e., the ith orbital associated to br is located at wp and transforms as the idxth partner function of the site-symmetry irrep of br.siteir. The total number of orbitals associated to a band representation, and hence the length of ordering, is the product of the site-symmetry irrep dimensionality and the number of sites in the Wyckoff position orbit.

source
SymmetricTightBinding.SiteInducedSGRepElementType
SiteInducedSGRepElement{D}(
    ρ::AbstractMatrix,
    positions::Vector{DirectPoint{D}},
    op::SymOperation{D}
)

Represents a matrix-valued element of a site-induced representation of a space group, including a global momentum-dependent phase factor.

This structure behaves like a functor: calling it with a momentum k :: AbstractVector returns the matrix representation at k.

Fields (internal)

  • ρ :: Matrix{ComplexF64} : The momentum-independent matrix part of the representation.
  • positions :: Vector{DirectPoint{D}}: Real-space positions corresponding to the orbitals in the orbit of the associated site-symmetry group.
source
SymmetricTightBinding._maybe_add_hoppings!Method

Computes and adds the symmetry related partners of a hopping term δ to the δ_orbit.

Warning

This function is an internal helper function for maybe_add_hoppings! and is not part of the public API.

source
SymmetricTightBinding._permute_symmetry_related_hoppings_under_symmetry_operationMethod

Build the P matrix for a particular symmetry operation acting on k-space, which permutes the rows of the M matrix.

For obtaining the P matrix, we make use that the action is on exponential of the type: $𝐞xp(2πk⋅δ)$, to instead act on δ ∈ h_orbit.orbit instead of k, which is a symbolic variable. Because of that, we need to use the inverse of the rotation part of the symmetry operation.

Sketch of proof

Assume g={R|τ} and Crystalline implements gk=(R⁻¹)ᵀk. Then (gk)⋅δ = ((R⁻¹)ᵀk)⋅δ + τ = k⋅(R⁻¹)δ.

Info

It is assumed that the operation op is provided in a primitive setting.

Warning

This function is an internal helper function for reciprocal_constraints_matrices and is not part of the public API.

source
SymmetricTightBinding.add_timereversal_related_orbits!Method
add_timereversal_related_orbits!(h_orbits::Vector{HoppingOrbit{D}}) where {D}

Adds the time-reversed hopping terms to the hopping orbits in h_orbits. The time-reversed hopping terms are added to the orbit of the hopping term they are related to, and if they are already present in another orbit, the two orbits are merged.

source
SymmetricTightBinding.construct_M_matrixMethod
construct_M_matrix(
    h_orbit::HoppingOrbit{D}, br1::NewBandRep{D}, br2::NewBandRep{D},
    [ordering1, ordering2]) 
    --> Array{Int,4}

Construct a set of matrices that encodes a Hamiltonian's term which resembles the hopping from EBR br1 to EBR br2.

The encoding is stored as a 4D matrix. Its last two axes correspond to elements of the Bloch Hamiltonian H(k); its first axis corresponds to orbit(h_orbit) and the associated complex exponentials stored in v; and its second axis to the elements of the vector t. That is:

Hₛₜ(k) = vᵢ(k) Mᵢⱼₛₜ tⱼ

See devdocs.md for details.

source
SymmetricTightBinding.evaluate_tight_binding_term!Method
evaluate_tight_binding_term!(tbt::TightBindingTerm, k, [c, H])

Evaluate the tight-binding term tbt at momentum k, possibly multiplied by a scalar coefficient c (unity if omitted). The term is added into the scratch space matrix H; if H is not provided, it is initialized as a zero matrix of the appropriate size.

The function returns the modified H matrix.

Note

The two-argument form of the function, i.e., returning the value of tbt at k, can be more simply achieved via tbt(k).

source
SymmetricTightBinding.maybe_add_hoppings!Method
maybe_add_hoppings!(h_orbits, δ, qₐ, qᵦ, R, ops) --> Vector{HoppingOrbit{D}}

Checks if a hopping term δ is already in the list of representatives. If not, it adds it and its symmetry related partners. If it is, it only adds the symmetry related partners.

source
SymmetricTightBinding.obtain_basis_free_parametersMethod
obtain_basis_free_parameters(
    h_orbit::HoppingOrbit{D},
    brₐ::NewBandRep{D}, 
    brᵦ::NewBandRep{D}, 
    [orderingₐ = OrbitalOrdering(brₐ), orderingᵦ = OrbitalOrdering(brᵦ)]
    )                            --> Tuple{Array{Int,4}, Vector{Vector{ComplexF64}}}

Obtain the basis of free parameters for the hopping terms between brₐ and brᵦ associated with the hopping orbit h_orbit.

Note

The presence or absence of time-reversal symmetry is inferred implicitly from brₐ and brᵦ.

source
SymmetricTightBinding.obtain_basis_free_parameters_TRSMethod
obtain_basis_free_parameters_TRS(
    h_orbit::HoppingOrbit{D}, 
    brₐ::NewBandRep{D}, 
    brᵦ::NewBandRep{D}, 
    orderingₐ::OrbitalOrdering{D} = OrbitalOrdering(brₐ),
    orderingᵦ::OrbitalOrdering{D} = OrbitalOrdering(brᵦ),
    Mm::AbstractArray{4, Int} = construct_M_matrix(h_orbit, brₐ, brᵦ, orderingₐ, orderingᵦ)
    )                             --> Tuple{Array{Int,4}, Vector{Vector{ComplexF64}}}}

Obtain the basis of free parameters for the hopping terms between brₐ and brᵦ associated with the hopping orbit h_orbit under time-reversal symmetry.

Real and imaginary parts of the basis vectors are differentiated explicitly: internally, we consider only variables.

source
SymmetricTightBinding.obtain_basis_free_parameters_hermiticityMethod
obtain_basis_free_parameters_hermiticity(
    h_orbit::HoppingOrbit{D},
    brₐ::NewBandRep{D},
    brᵦ::NewBandRep{D},
    orderingₐ::OrbitalOrdering{D} = OrbitalOrdering(brₐ),
    orderingᵦ::OrbitalOrdering{D} = OrbitalOrdering(brᵦ),
    Mm::AbstractArray{Int, 4} = construct_M_matrix(h_orbit, brₐ, brᵦ, orderingₐ, orderingᵦ);
    antihermitian::Bool = false,
) where {D}

Constructs a basis for the coefficient vectors t⁽ⁿ⁾ that span the space of Hermitian (or antihermitian if true) TB Hamiltonians Hₛₜ(k) = vᵢ(k) Mᵢⱼₛₜ tⱼ = vᵀ(k) M⁽ˢᵗ⁾ t. We do this by assuming that each coefficient vector t is sorted into a vector of the form [tᴿ; itᴵ] so that we can take the complex conjugate by as t* = σ₃t, which can then be moved onto M⁽ˢᵗ⁾ instead of t. The constraint Hₛₜ(k) = (H†)ₛₜ(k) = Hₜₛ*(k) can then be expressed as vᵀ(k) M⁽ˢᵗ⁾ tⱼ = v*ᵀ(k) M⁽ᵗˢ⁾ σ₃ t = vᵀ(k) (Pᵀ M⁽ᵗˢ⁾σ₃) t, which requires that t be a solution to the nullspace M⁽ˢᵗ⁾ - Pᵀ M⁽ᵗˢ⁾ σ₃ = 0. We cast this as Z - Q = 0, with Z = M⁽ˢᵗ⁾ and Q = Pᵀ M⁽ᵗˢ⁾ σ₃.

Notes

For anti-Hermitian symmetry, we require Hₛₜ(k) = -Hₜₛ*(k), which translates to M⁽ˢᵗ⁾ + Pᵀ M⁽ᵗˢ⁾ σ₃ = 0; i.e., simply swaps the sign of Q

source
SymmetricTightBinding.pin_freeMethod
pin_free(br::NewBandRep{D}, αβγ::AbstractVector{<:Real}) where D

Pin the free parameters of the Wyckoff position associated with the band representation br to the values in αβγ.

Returns a new band representation with all other properties, apart from the Wyckoff position, identical to (and sharing memory with) br.

source
SymmetricTightBinding.primitivized_orbitMethod
primitivized_orbit(br::NewBandRep{D}) where D

Return the orbit of the Wyckoff position associated with the band representation br. The coordinates of positions in the orbit are given relative to the primitive unit cell.

Positions are returned as a Vector{DirectPoint{D}}.

The following checks are made, producing an error if violated:

  1. There are no free parameters associated with the Wyckoff position.
  2. For every position, its coordinates, referred to the primitive basis, is in the range [0,1); i.e., every position lies in the parallepiped primitive unit cell [0,1)ᴰ.
source
SymmetricTightBinding.reciprocal_constraints_matricesMethod
reciprocal_constraints_matrices(
                                Mm::AbstractArray{Int,4}, 
                                gens::AbstractVector{SymOperation{D}}, 
                                h_orbit::HoppingOrbit{D}
                                ) --> Vector{Array{Int,4}}

Compute the reciprocal constraints matrices for the generators of the SG. This is done by permuting the rows of the M matrix according to the symmetry operation acting on k-space. See more details in permute_symmetry_related_hoppings_under_symmetry_operation and devdocs.md.

source
SymmetricTightBinding.reciprocal_constraints_trsMethod
reciprocal_constraints_trs(Mm::AbstractArray{Int,4}, h_orbit::HoppingOrbit{D}) 
--> Array{ComplexF64,4}

Time reversal symmetry action on reciprocal space. It is given by the association k -> -k => H(k) -> H(-k).

source
SymmetricTightBinding.representation_constraint_matricesMethod
representation_constraints_matrices(
    Mm::AbstractArray{Int,4}, 
    brₐ::NewBandRep{D},
    brᵦ::NewBandRep{D}) --> Vector{Array{ComplexF64,4}}

Build the Q matrix for a particular symmetry operation (or, equivalently, a particular matrix from the site-symmetry representation), acting on the M matrix. Relative to our white-board notes, Q has swapped indices, in the sense we below give Q[i,j,r,l].

(ρₐₐ)ᵣₛ Hₛₜ (ρᵦᵦ⁻¹)ₜₗ = (ρₐₐ)ᵣₛ vᵢ Mᵢⱼₛₜ tⱼ (ρᵦᵦ⁻¹)ₜₗ = vᵢ (ρₐₐ)ᵣₛ Mᵢⱼₛₜ (ρᵦᵦ⁻¹)ₜₗ tⱼ,

then we can define: Qᵢⱼᵣₗ = (ρₐₐ)ᵣₛ Mᵢⱼₛₜ (ρᵦᵦ⁻¹)ₜₗ

source
SymmetricTightBinding.representation_constraint_trsMethod
representation_constraint_trs(Mm::AbstractArray{Int,4}, h_orbit::HoppingOrbit{D})
--> Array{ComplexF64,4}

Time reversal symmetry action on the Hamiltonian. It is given by the association δ -> -δ and the complex conjugation in the free-parameter part: $tⱼ -> tⱼ* ⇒ [tⱼᴿ, itⱼᴵ] -> [tⱼᴿ, -itⱼᴵ] ⇒ H(k) -> H*(k)$.

source
SymmetricTightBinding.sgrep_induced_by_siteir_excl_phaseMethod
sgrep_induced_by_siteir_excl_phase(br::NewBandRep, op::SymOperation)
sgrep_induced_by_siteir_excl_phase(cbr::CompositeBandRep, op::SymOperation)
    --> Matrix{ComplexF64}

Return the representation matrix of a symmetry operation op induced by the site symmetry group of a band representation br or composite band representation cbr, excluding the global momentum-dependent phase factor.

source
SymmetricTightBinding.split_complexMethod
split_complex(t::Vector{<:Number}) -> Matrix{Real}

Consider αt where α ∈ ℂ and t ∈ ℂⁿ and build from t a matrix representation T that allows access to the real and imaginary parts of the product αt without using complex numbers by splitting α into a real 2-vector of its real and imaginary parts.

In particular, let $α = αᴿ + iαᴵ$ and $t = tᴿ + itᴵ$ with αᴿ, αᴵ ∈ ℝ and $tᴿ, tᴵ ∈ ℝⁿ$, then $αt$ can be rewritten as

\[αt = (αᴿ + iαᴵ)(tᴿ + itᴵ) = (αᴿtᴿ - αᴵtᴵ) + i(αᴿtᴵ + αᴵtᴿ) = [tᴿ, tᴵ]ᵀ [αᴿ, αᴵ] + i [tᴵ, tᴿ]ᵀ [αᴿ, αᴵ]\]

Then, defining T = [tᴿ -tᴵ; tᴵ tᴿ], the above product can then be reexpressed as: $Re(αt) = αᴿtᴿ - αᴵtᴵ =$ (T * [αᴿ; αᴵ])[1:n] and $Im(αt) = αᴿtᴵ + αᴵtᴿ =$ (T * [αᴿ; αᴵ])[n+1:2n]. I.e., the "upper half" of the product T * [real(α), imag(α)] is real(α * t) and the "lower half" is imag(αt).

This functionality is used to avoid complex numbers in amplitude basis coefficients, which simplifies the application of time-reversal symmetry and hermiticity.

Examples

julia> using SymmetricTightBinding: split_complex

julia> t = [im,0]
2-element Vector{Complex{Int64}}:
 0 + 1im
 0 + 0im

julia> T = split_complex(t)
4×2 Matrix{Int64}:
 0  -1
 0   0
 1   0
 0   0

julia> α = 0.5+0.2im; αv = [real(α), imag(α)];

julia> (T * αv)[1:2] == real(α*t) && (T * αv)[3:4] == imag(α*t)
julia> t = [1,im]
2-element Vector{Complex{Int64}}:
 1 + 0im
 0 + 1im

julia> split_complex(t)
4×2 Matrix{Int64}:
 1   0
 0  -1
 0   1
 1   0
source
SymmetricTightBinding.zassenhaus_intersectionMethod
zassenhaus_intersection(U::AbstractArray{<:Number}, W::AbstractArray{<:Number}) 
    --> AbstractArray{<:Number}

Finds the intersection of two bases U and W using the Zassenhaus algorithm. It assumes that the basis are given by columns.

References

  • https://en.wikipedia.org/wiki/Zassenhaus_algorithm
source