k-space paths
Default paths
To generate a k-path for, say, the space group of diamond (space group 227; a cubic face-centered Bravais lattice), we can call irrfbz_path
, which will return a minimal path in the irreducible Brillouin zone:
using Brillouin
sgnum = 227
Rs = [[1,0,0], [0,1,0], [0,0,1]] # conventional direct basis
kp = irrfbz_path(sgnum, Rs)
KPath{3} (6 points, 2 paths, 8 points in paths):
points: :U => [0.625, 0.25, 0.625]
:W => [0.5, 0.25, 0.75]
:K => [0.375, 0.375, 0.75]
:Γ => [0.0, 0.0, 0.0]
:L => [0.5, 0.5, 0.5]
:X => [0.5, 0.0, 0.5]
paths: [:Γ, :X, :U]
[:K, :Γ, :L, :W, :X]
basis: [-6.283185, 6.283185, 6.283185]
[6.283185, -6.283185, 6.283185]
[6.283185, 6.283185, -6.283185]
The path data is sourced from the HPKOT paper (or, equivalently, the SeeK-path Python package).
The coordinates of the path are given with respect to the primitive reciprocal basis (here, [[-2π,2π,2π], [2π,-2π,2π], [2π,2π,-2π]]
). To convert to a Cartesian basis, we can use cartesianize
or cartesianize!
(in-place):
cartesianize(kp)
KPath{3} (6 points, 2 paths, 8 points in paths):
points: :U => [1.570796, 6.283185, 1.570796]
:W => [3.141593, 6.283185, 0.0]
:K => [4.712389, 4.712389, 0.0]
:Γ => [0.0, 0.0, 0.0]
:L => [3.141593, 3.141593, 3.141593]
:X => [0.0, 6.283185, 0.0]
paths: [:Γ, :X, :U]
[:K, :Γ, :L, :W, :X]
basis: [-6.283185, 6.283185, 6.283185]
[6.283185, -6.283185, 6.283185]
[6.283185, 6.283185, -6.283185]
We can visualize the k-path using PlotlyJS.jl (conversion to a Cartesian basis for plotting is automatic):
using PlotlyJS
Pᵏ = plot(kp)
Usually, it'll be more helpful to understand the path's geometry in the context of the associated Brillouin zone. To visualize this, we can plot the combination of a Cell
(created via wignerseitz
) and a KPath
:
pGs = basis(kp) # primitive reciprocal basis associated with k-path
c = wignerseitz(pGs) # associated Brillouin zone
Pᶜ⁺ᵏ = plot(c, kp)
Interpolation
Interpolation of a KPath
structure can be achieved using either interpolate(::KPath, ::Integer)
or splice(::KPath, ::Integer)
, returning a KPathInterpolant
. As an example, interpolate(kp, N)
returns an interpolation with a target of N
interpolation points, distributed as equidistantly as possible (with the distance metric evaluated in Cartesian space):
kpi = interpolate(kp, 100)
99-element KPathInterpolant{3}:
[0.0, 0.0, 0.0]
[0.022727272727272728, 0.0, 0.022727272727272728]
[0.045454545454545456, 0.0, 0.045454545454545456]
[0.06818181818181818, 0.0, 0.06818181818181818]
[0.09090909090909091, 0.0, 0.09090909090909091]
[0.11363636363636363, 0.0, 0.11363636363636363]
[0.13636363636363635, 0.0, 0.13636363636363635]
[0.1590909090909091, 0.0, 0.1590909090909091]
[0.18181818181818182, 0.0, 0.18181818181818182]
[0.20454545454545453, 0.0, 0.20454545454545453]
⋮
[0.5, 0.18181818181818182, 0.6818181818181817]
[0.5, 0.1590909090909091, 0.6590909090909091]
[0.5, 0.13636363636363635, 0.6363636363636364]
[0.5, 0.11363636363636365, 0.6136363636363636]
[0.5, 0.09090909090909093, 0.5909090909090908]
[0.5, 0.06818181818181818, 0.5681818181818181]
[0.5, 0.045454545454545456, 0.5454545454545454]
[0.5, 0.022727272727272735, 0.5227272727272727]
[0.5, 0.0, 0.5]
The returned KPathInterpolant
implements the AbstractVector
interface, with iterants returning SVector{D, Float64}
elements. To get a conventional "flat" vector, we can simply call collect(kpi)
.
Internally, KPathInterpolant
includes additional structure and information: namely, the high-symmetry points and associated labels along the path as well as a partitioning into connected vs. disconnected path segments.
Band structure
The additional structure of KPathInterpolation
enables convenient and clear visualizations of band structure diagrams in combination with PlotlyJS.
To illustrate this, suppose we are considering a tight-binding problem for an s-orbital situated at the 1a Wyckoff position. Such a problem has a single band with dispersion [1] (assuming a cubic side length $a = 1$):
\[\epsilon(\mathbf{k}) = 4\gamma\Bigl( \cos \tfrac{1}{2}k_x \cos \tfrac{1}{2}k_y + \cos \tfrac{1}{2}k_y \cos \tfrac{1}{2}k_z + \cos \tfrac{1}{2}k_z \cos \tfrac{1}{2}k_x \Bigr)\]
with $k_{x,y,z}$ denoting coordinates in a Cartesian basis (which are related to the coordinates $k_{1,2,3}$ in a primitive reciprocal basis by $k_x = 2π(-k_1+k_2+k_3)$, $k_x = 2π(k_1-k_2+k_3)$, and $k_z = 2π(k_1+k_2-k_3)$).
We can calculate the energy band along our k-path using the interpolation object kpi
. To do so, we define a function that implements $\epsilon(\mathbf{k})$ and broadcast it over the elements of kpi
:
function ϵ(k; γ::Real=1.0)
kx = 2π*(-k[1]+k[2]+k[3])
ky = 2π*(+k[1]-k[2]+k[3])
kz = 2π*(+k[1]+k[2]-k[3])
return 4γ * (cos(kx/2)*cos(ky/2) + cos(ky/2)*cos(kz/2) + cos(kz/2)*cos(kx/2))
end
band = ϵ.(kpi)
99-element Vector{Float64}:
12.0
11.918571535047462
11.67594378891598
11.277055962836148
10.73002826264945
10.045996594834067
9.238885871562282
8.325126539644781
7.3233201040150915
6.253860454731438
⋮
-4.0
-3.9999999999999996
-4.0
-3.9999999999999996
-4.0
-4.0
-4.0
-4.0
-4.0
Finally, we can visualize the associated band using a Brillouin-overloaded PlotlyJS plot
-call:
P = plot(kpi, [band])
If we have multiple bands, say $\epsilon_1(\mathbf{k}) = \epsilon(\mathbf{k})$ and $\epsilon_2(\mathbf{k}) = 20 - \tfrac{1}{2}\epsilon(\mathbf{k})$, we can easily plot that by collecting the bands in a single vector (or concatenating into a matrix):
band1 = ϵ.(kpi)
band2 = 20 .- (1/2).*band1
P¹² = plot(kpi, [band1, band2])
PlotlyJS.plot
— Functionplot(kpi::KPathInterpolant, bands, [layout]; kwargs...)
Plot a dispersion diagram for provided bands
and k-path interpolant kpi
.
bands
must be an iterable of iterables of <:Real
s (e.g., a Vector{Vector{Float64}}
), with the first iteration running over distinct energy bands, and the second running over distinct k-points in kpi
. Note that the length of each iterant of bands
must equal length(kpi)
.
Alternatively, bands
can be an AbstractMatrix{<:Real}
, with columns interpreted as distinct energy bands and rows as distinct k-points.
A layout
can be supplied to overwrite default layout choices (set by BrillouinPlotlyJSExt.DEFAULT_PLOTLY_LAYOUT_DISPERSION
). Alternatively, some simple settings can be set directly via keyword arguments (see below).
Keyword arguments kwargs
ylims
: y-axis limits (default: quasi-tight aroundbands
's range)ylabel
: y-axis label (default: "Energy")title
: plot title (default:nothing
); can be aString
or anattr
dictionary of PlotlyJS propertiesband_highlights
: dictionary of non-default styling for specified band ranges (default:nothing
, indicating all-default styling).Example: To color bands 2 and 3 black, set
band_highlights = Dict(2:3 => attr(color=:black, width=3))
. Unlisted band ranges are shown in default band styling.annotations
: dictionary of hover-text annotations for labeled high-symmetry points inkpi
(default:nothing
, indicating no annotations). Suitable for labeling of irreps.Example: Assume bands 1 and 2 touch at X, but not at Γ. To label this, we set:
annotations = Dict(:X => [1:2 => "touching!"], :Γ => [1 => "isolated", 2 => "isolated"])
. If a band-range is provided, a single annotation is placed at the mean of the energies at these band-ranges. Alternatively, if the first element of each pair is a non-Integer
Real
number, it is interpreted as referring to the frequency of the annotation.
Treating unit cells in non-standard settings
irrfbz_path(sgnum, Rs)
requires Rs
to be provided in a standard setting. Often, the setting of Rs
may not be standard and it can be a hassle to convert existing calculations to such a setting. To avoid this, we can provide the unit cell information to irrfbz
via Spglib
's Cell
format (this functionality depends on separately loading the Spglib.jl package). If the provided unit cell is a supercell of a smaller primitive cell, irrfbz_path
returns the standard k-path of the smaller primitive cell in the basis of the supercell reciprocal lattice vectors.
For example, to construct a k-path for a non-standard trigonal lattice in space group 166 (R-3m):
using Spglib
a = 1.0
c = 8.0
pRs_standard = [[a*√3/2, a/2, c/3],
[-a*√3/2, a/2, c/3],
[0, -a, c/3]]
pRs_nonstandard = [[a*√3/2, -a/2, c/3],
[0, a, c/3],
[-a*√3/2, -a/2, c/3]]
cell_standard = Spglib.Cell(pRs_standard, [[0, 0, 0]], [0])
cell_nonstandard = Spglib.Cell(pRs_nonstandard, [[0, 0, 0]], [0])
kp_standard = irrfbz_path(cell_standard)
kp_nonstandard = irrfbz_path(cell_nonstandard)
KPath{3} (8 points, 3 paths, 10 points in paths):
points: :T => [0.5, 0.5, 0.5]
:F => [0.5, -0.0, 0.5]
:S₀ => [0.345052, -0.345052, 0.0]
:S₂ => [0.654948, 0.0, 0.345052]
:H₀ => [0.5, -0.190104, 0.190104]
:Γ => [0.0, 0.0, 0.0]
:L => [0.5, -0.0, 0.0]
:H₂ => [0.809896, 0.190104, 0.5]
paths: [:Γ, :T, :H₂]
[:H₀, :L, :Γ, :S₀]
[:S₂, :F, :Γ]
basis: [3.627599, -2.094395, 0.785398]
[0.0, 4.18879, 0.785398]
[-3.627599, -2.094395, 0.785398]
Note that the space group symmetry is inferred by Spglib from the atomic positions and provided basis.
We can check that the generated k-paths for the standard and non-standard lattices are equivalent by plotting and comparing their k-paths and associated Wigner-Seitz cells:
using Bravais # for `reciprocalbasis`
Pˢ = plot(wignerseitz(reciprocalbasis(pRs_standard)), kp_standard, Layout(title="standard cell"))
Pⁿˢ = plot(wignerseitz(reciprocalbasis(pRs_nonstandard)), kp_nonstandard, Layout(title="non-standard cell"))