API


Functions

Crystalline.collect_compatibleMethod
collect_compatible(ptbm::ParameterizedTightBindingModel{D}; multiplicities_kws...)

Determine a decomposition of the bands associated with ptbm into a set of SymmetryVectors, with each symmetry vector corresponding to a set of compatibility-respecting (i.e., energy separable along high-symmetry k-lines) bands.

Keyword arguments

  • multiplicities_kws...: keyword arguments passed to Crystalline.collect_compatible used in determining the multiplicities of irreps across high-symmetry k-points.

Example

julia> using Crystalline, SymmetricTightBinding

julia> brs = calc_bandreps(221);

julia> cbr = @composite brs[1] + brs[2]
40-irrep CompositeBandRep{3}:
 (3d|A₁g) + (3d|A₁ᵤ) (6 bands)

julia> tbm = tb_hamiltonian(cbr); # a 4-term, 6-band model

julia> ptbm = tbm([1.0, 0.1, -1.0, 0.1]); # fix free coefficients

julia> collect_compatible(ptbm)
2-element Vector{SymmetryVector{3}}:
 [M₅⁺+M₁⁻, X₃⁺+X₁⁻+X₂⁻, Γ₁⁻+Γ₃⁻, R₄⁺] (3 bands)
 [M₁⁺+M₅⁻, X₁⁺+X₂⁺+X₃⁻, Γ₁⁺+Γ₃⁺, R₄⁻] (3 bands)

In the above example, the bands separate into two symmetry vectors, one for each of the original EBRs in cbr.

source
Crystalline.collect_irrep_annotationsMethod
collect_irrep_annotations(ptbm::ParameterizedTightBindingModel; kws...)

Collect the irrep labels across the high-symmetry k-points referenced by the underlying composite band representation of ptbm, across the bands of the model.

Useful for annotating irrep labels in band structure plots (via the Makie extension call plot(ks, energies; annotations=collect_irrep_annotations(ptbm)))

source
SymmetricTightBinding.energy_gradient_wrt_hoppingMethod
energy_gradient_wrt_hopping(
    ptbm::ParameterizedTightBindingModel{D},
    k::ReciprocalPointLike{D}
    (Es, us) = solve(ptbm, k; bloch_phase=Val(false));
    degen_rtol::Float64 = 1e-12,
    degen_atol::Float64 = 1e-12
) where D

Return the hopping gradient of the energy of each band in ptbm evaluated at momentum k.

The gradient is computed using the Feynman-Hellmann theorem. For degenerate bands (assessed energetically using relative and absolute tolerances degen_rtol and degen_atol), a degenerate variant is used, equivalent to degenerate perturbation theory.

The gradient is returned as column vectors, one for each band, with each column containing the gradient of the corresponding energy with respect to the hopping coefficients of ptbm.

source
SymmetricTightBinding.gradient_wrt_hoppingMethod
gradient_wrt_hopping(tbm :: TightBindingModel)
gradient_wrt_hopping(ptbm :: ParameterizedTightBindingModel)

Return a structure that encodes the gradient of a tight-binding model tbm or ptbm with respect to the hopping coefficients.

To evaluate the gradient at a particular momentum k, use the returned structure as a functor at k. I.e., gradient(tbm)(k) returns the gradient of the tight-binding Hamiltonian with respect to all hoppping coefficients at momentum k. This gradient is a vector of matrices.

source
SymmetricTightBinding.obtain_symmetry_related_hoppingsMethod
obtain_symmetry_related_hoppings(
    Rs::AbstractVector{V}, 
    brₐ::NewBandRep{D}, 
    brᵦ::NewBandRep{D};
    hermiticity::Bool = true,
) --> Vector{HoppingOrbit{D}}

Compute the symmetry related hopping terms from the points in WP of brₐ to the WP of brᵦ displaced a set of primitive lattice vectors representatives Rs.

Implementation

  1. Take a point a in the WP of brₐ and a point b in the WP of brᵦ. We

compute the displacement vector δ = b + R - a, where R ∈ Rs.

  1. If δ ∈ representatives then we add δ => (a, b, R) to the list of hoppings of that representative and continue. If not then, we search inside of all the representatives for the one that δ => (a, b, R) in the list of hoppings. If not found, then we add δ as a new representative and add δ => (a, b, R) to its list of hoppings.
  2. Take g ∈ generators and compute δ' = g δ and (a', b', R') = (g a, g b, g R), and repeat step 2.
  3. Repeat all steps 1 to 3 for all pair of points in the WPs of brₐ and brᵦ.

Additionally, if we have time-reversal symmetry or/and hermiticity , we check if orbits that relate δ and are present; if not, we add them. The presence or absence of time-reversal symmetry is automatically inferred from brₐ and brᵦ (which must be identical).

source
SymmetricTightBinding.pin_free!Method
pin_free!(
    brs::Collection{NewBandRep{D}},
    idx2αβγ::Pair{Int, <:AbstractVector{<:Real}}
)

pin_free!(
    brs::Collection{NewBandRep{D}},
    idx2αβγs::AbstractVector{<:Pair{Int, <:AbstractVector{<:Real}}}
)

For idx2αβγ = idx => αβγ, update brs[idx] such that the free parameters of its associated Wyckoff positions are pinned to αβγ.

A vector of pairs idx2αβγs can also be provided, to pin multiple distinct band representations.

See also pin_free for non-mutated input.

source
SymmetricTightBinding.sgrep_induced_by_siteirMethod
sgrep_induced_by_siteir(
    br::Union{NewBandRep, CompositeBandRep},
    op::SymOperation, [positions::Vector{<:DirectPoint}]
)
sgrep_induced_by_siteir(
    tbm::Union{TightBindingModel,ParameterizedTightBindingModel}, op::SymOperation
)
    --> SiteInducedSGRepElement

Computes the representation matrix of a symmetry operation op induced by the site symmetry group associated with an elementary or composite band representation br , including the global momentum-dependent phase factor, returning a SiteInducedSGRepElement, which is a functor over momentum inputs.

A (possibly parameterized) tight-binding model tbm can be specified instead of a band representation, in which case the latter is inferred from the former.

source
SymmetricTightBinding.spectrumMethod
spectrum(ptbm::ParameterizedTightBindingModel, ks; transform = identity)

Evaluate the spectrum, i.e., energies, of the tight-binding model ptbm over an iterable of input momenta ks.

Energies are returned as a matrix, with rows running over momenta and columns over distinct bands.

Keyword arguments

  • transform: a function to apply to the resulting matrix of energies, defaulting to the identity function. This can be used to e.g., convert the energies to a different scaling.

Example

As an example, we evaluating the band structure of graphene. Below, we first construct and parameterize a tight-binding model for the the (2b|A₁) EBR in plane group 17, corresponding to the highest-lying orbitals in graphene. Next, we construct a path along high-symmetry directions of the Brillouin zone using Brillouin.jl, calculate the spectrum across this path; and finally, plot the band structure using Brillouin and GLMakie (or PlotlyJS):

julia> using Crystalline, SymmetricTightBinding

julia> brs = calc_bandreps(17, Val(2));

julia> cbr = @composite brs[5]
13-irrep CompositeBandRep{2}:
 (2b|A₁) (2 bands)

julia> ptbm = tb_hamiltonian(cbr, [zeros(Int, dim(cbr))])([0.0, 1.0]);

julia> using Brillouin, GLMakie

julia> kpi = interpolate(irrfbz_path(17, directbasis(17, Val(2))), 100);

julia> plot(kpi, spectrum(ptbm, kpi))
source
SymmetricTightBinding.spectrumMethod
spectrum(ptbm::ParameterizedTightBindingModel, k::AbstractVector{<:Real})

Evaluate the spectrum, i.e., energies, of the tight-binding model ptbm at a single momentum k, across all the bands of ptbm.

source
SymmetricTightBinding.subduced_complementMethod
subduced_complement(tbm::TightBindingModel{D}, sgnumᴴ::Int; timereversal)
                                                    --> TightBindingModel{D}

Given a model tbm associated with a space group $G$, determine the new, independent tight-binding terms (i.e., the the orthogonal complement of terms) that become symmetry-allowed when the model's space group is reduced to a subgroup $H ≤ G$ with space group number sgnumᴴ and time-reversal symmetry timereversal.

Practically, the function answers the question: which new tight-binding terms become allowed if the symmetry of the model is reduced from space group $G$ to subgroup $H$?

Implementation

The function computes a basis of allowed tight-binding terms in the subgroup setting $H$ by simply restricting the constraints in $G$ to generators in $H$. This gives a basis for the tight-binding terms in the subduced $G ↓ H$ setting. The space spanned by this basis is compared to the space spanned in the original model; in particular new terms are identified as the orthogonal complement of the spaces associated with $G ↓ H$ relative to $G$.

Keywords

  • timereversal::Bool: Specifies whether time-reversal symmetry is present in the subgroup $H$. By default, the presence or absence is inherited from the original model tbm. Note that timereversal must be "lower or equal to" the time-reversal of the original model.

Example

It is well-known that the Dirac point of graphene is gapped under mirror and time-reversal symmetry breaking. We can see this by constructing a tight-binding model first for a model of graphene (plane group ⋕17) and then subducing it to a setting without mirror symmetry (plane group ⋕16) and without time-reversal symmetry (timereversal = false). First, we construct the tight-binding model for graphene (via the (2a|A₁) band representation):

julia> using SymmetricTightBinding, Crystalline

julia> brs = calc_bandreps(17, Val(2); timereversal = true);

julia> cbr = @composite brs[5]

julia> tbm = tb_hamiltonian(cbr, [[0,0], [1,0]])

Each of the 4 terms in this model is proportional to an identity matrix at K = (1/3, 1/3). Using subduced_complement, we can find the new terms that appear if we imagine lowering the symmetry from plane group ⋕17 to ⋕16 (which has no mirror symmetry) while also removing time-reversal symmetry.

julia> Δtbm = subduced_complement(tbm, 16; timereversal = false)
2-term 2×2 TightBindingModel{2} over (2b|A₁):
┌─
1. ⎡ i𝕖(δ₁)+i𝕖(δ₂)+i𝕖(δ₃)-i𝕖(δ₄)-i𝕖(δ₅)-i𝕖(δ₆)  0                                          ⎤
│  ⎣ 0                                          -i𝕖(δ₁)-i𝕖(δ₂)-i𝕖(δ₃)+i𝕖(δ₄)+i𝕖(δ₅)+i𝕖(δ₆) ⎦
└─ (2b|A₁) self-term:  δ₁=[1,0], δ₂=[0,1], δ₃=[-1,-1], δ₄=-δ₁, δ₅=-δ₂, δ₆=-δ₃
┌─
2. ⎡ 0                                       𝕖(δ₁)+𝕖(δ₂)+𝕖(δ₃)-𝕖(δ₇)-𝕖(δ₈)-𝕖(δ₉) ⎤
│  ⎣ 𝕖(δ₄)+𝕖(δ₅)+𝕖(δ₆)-𝕖(δ₁₀)-𝕖(δ₁₁)-𝕖(δ₁₂)  0                                   ⎦
└─ (2b|A₁) self-term:  δ₁=[4/3,-1/3], δ₂=[1/3,5/3], δ₃=[-5/3,-4/3], δ₄=-δ₁, δ₅=-δ₂, δ₆=-δ₃, δ₇=[1/3,-4/3], δ₈=[-5/3,-1/3], δ₉=[4/3,5/3], δ₁₀=-δ₇, δ₁₁=-δ₈, δ₁₂=-δ₉

The first of the of these terms is not diagonal at K and so opens a gap at the Dirac point:

julia> Δtbm[1](ReciprocalPoint(1/3, 1/3))
2×2 Matrix{ComplexF64}:
 5.19615+1.73195e-14im       0.0+0.0im
     0.0+0.0im          -5.19615+1.43774e-14im

Adding symmetry-breaking terms to the original model

To build a "complete" model, with both the original and symmetry-breaking terms, use vcat:

julia> tbm′ = vcat(tbm, Δtbm); length(tbm′) == length(tbm) + length(Δtbm)
true

Limitations

The subgroup $H$ must be a volume-preserving subgroup of the original group $G$. I.e. $H$ must be a translationen-gleiche subgroup of $G$ (or $G$ itself), and there must exist a transformation from $G$ to $H$ that preserves volume (i.e., has det(t.P) == 1 for t denoting an element returned by Crystalline.jl's conjugacy_relations).

source
SymmetricTightBinding.symmetry_eigenvaluesMethod
symmetry_eigenvalues(
    ptbm::ParameterizedTightBindingModel{D},
    ops::AbstractVector{SymOperation{D}},
    k::ReciprocalPointLike{D},
    [sgreps::AbstractVector{SiteInducedSGRepElement{D}}]
)
symmetry_eigenvalues(
    ptbm::ParameterizedTightBindingModel{D},
    lg::LittleGroup{D},
    [sgreps::AbstractVector{SiteInducedSGRepElement{D}}]
)
    --> Matrix{ComplexF64}

Compute the symmetry eigenvalues of a coefficient-parameterized tight-binding model ptbm at the k-point k for the symmetry operations ops. A LittleGroup can also be provided instead of ops and k.

Representations of the symmetry operations ops as acting on the orbitals of the tight-binding setting can optionally be provided in sgreps (see sgrep_induced_by_siteir) and are otherwise initialized by the function.

The symmetry eigenvalues are returned as a matrix, with rows running over the elements of ops and columns running over the bands of ptbm.

source
SymmetricTightBinding.tb_hamiltonianMethod
tb_hamiltonian(cbr::CompositeBandRep{D}, Rs::AbstractVector{Vector{Int}}) 
    --> Vector{TightBindingTerm{D}}

Construct the TB Hamiltonian matrix from a given composite band representation cbr and a set of global translation-representatives Rs. The Hamiltonian is constructed block by block according to the symmetry-related hoppings between the band representations in cbr. Several models returned, each representing a term that is closed under the symmetry operations of the underlying space group.

source

Types

SymmetricTightBinding.HoppingOrbitType
HoppingOrbit{D}

A structure holding information about symmetry-related spatial hopping vectors.

An orbit includes all symmetry related vectors {δ} = {δ₁, δ₂, …} obtained by applying the symmetry operations of the underlying space group to a representative vector δ.

For each element of the orbit δᵢ there may be multiple hopping terms, from site a to site b, possibly related by a lattice translation R. Each such set of possible hopping term (a, b, R) is for each δᵢ is stored elements as the ith element of the vector hoppings.

Fields

  • representative :: RVec{D}: the representative hopping vector δ
  • orbit :: Vector{RVec{D}}: the δᵢ elements of the orbit generated by δ. Generally orbit[1] == representative, i.e., δ₁ = δ
  • hoppings :: Vector{Vector{NTuple{3,RVec{D}}}}: the ith element gives the possible physical hopping terms (a,b,R) associated to orbit[i]. Multiple physical hopping terms may correspond to each element of the orbit. For (a, b, R) = hoppings[i][j], we have orbit[i] = δᵢ = b + R - a
source
SymmetricTightBinding.ParameterizedTightBindingModelType
ParameterizedTightBindingModel{D}

A coefficient-parameterized tight-binding model, that can be used as a functor for evaluation at input momenta k.

Fields

  • tbm :: TightBindingModel{D}: A tight-binding model, consisting of a set of a list of TightBindingTerm{D}s.
  • cs :: Vector{Float64}: A vector of coefficients, each associated to a corresponding element of tbm.
  • scratch :: Matrix{ComplexF64}: A scratch space for evaluating the Hamiltonian matrix at at specific momenta. This is a N×N matrix, where N is the number of orbitals in tbm (i.e., tbm.N). The scratch space is instantiated automatically on construction.

Functor over momenta

A ParameterizedTightBindingModel ptbm can be be evaluated at any ´D-dimensional momentumkby usingptbmas a functor. That is,ptbm(k)returns a numerical representation of the Hamiltonian matrix forptbmevaluated at momentumk.

source
SymmetricTightBinding.TightBindingElementStringType
TightBindingElementString

A structure for pretty-printing tight-binding matrix elements.

Fields

  • s :: String: the string representing the tight-binding matrix element
  • active :: Bool: whether the element belongs to an "active" block - i.e., one we want to highlight (then shown in blue).
source
SymmetricTightBinding.TightBindingModelType
TightBindingModel{D}

A structure storing a list of TightBindingTerm{D}s. Each term is assumed to associated with an identical list of EBRs.

To associate a set of coefficients to each term, see ParameterizedTightBindingModel, which also allows evaluation at specific momenta.

Fields

  • terms :: Vector{TightBindingTerm{D}}: a vector of TightBindingTerm{D}s, each of which represents a block (or conjugated pairs of blocks) of the Hamiltonian matrix.
  • cbr :: CompositeBandRep{D}: the composite band representation associated to the model.
  • positions :: Vector{DirectPoint{D}}: a vector of positions, specified in the lattice basis, associated to each orbital of the model.
  • N :: Int: the total number of orbitals in the model, i.e., the size of the Hamiltonian matrix associated to each element of terms.
source