This example plots the band structure of graphene, a 2D material. 2D band structures are not supported natively (yet), so we manually build a custom path in reciprocal space.
using DFTK
using Unitful
using UnitfulAtomic
using LinearAlgebra
# Define the convergence parameters (these should be increased in production)
L = 20 # height of the simulation box
kgrid = [6, 6, 1]
Ecut = 15
temperature = 1e-3
# Define the geometry and pseudopotential
a = 4.66 # lattice constant
a1 = a*[1/2,-sqrt(3)/2, 0]
a2 = a*[1/2, sqrt(3)/2, 0]
a3 = L*[0 , 0 , 1]
lattice = [a1 a2 a3]
C1 = [1/3,-1/3,0.0] # in reduced coordinates
C2 = -C1
positions = [C1, C2]
C = ElementPsp(:C, psp=load_psp("hgh/pbe/c-q4"))
atoms = [C, C]
# Run SCF
model = model_PBE(lattice, atoms, positions; temperature)
basis = PlaneWaveBasis(model; Ecut, kgrid)
scfres = self_consistent_field(basis)
# Construct 2D path through Brillouin zone
sgnum = 13 # Graphene space group number
kpath = irrfbz_path(model; dim=2, sgnum)
plot_bandstructure(scfres, kpath; kline_density=20)
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -11.15667314117 -0.60 5.9 2 -11.16017717058 -2.46 -1.30 1.0 233ms 3 -11.16039477848 -3.66 -2.33 2.1 321ms 4 -11.16041671082 -4.66 -3.29 2.7 340ms 5 -11.16041704729 -6.47 -3.47 3.1 375ms 6 -11.16041704985 -8.59 -3.63 1.0 238ms 7 -11.16041705076 -9.04 -3.87 1.3 270ms 8 -11.16041705117 -9.39 -4.29 2.0 273ms 9 -11.16041705134 -9.77 -4.60 2.0 283ms 10 -11.16041705140 -10.25 -4.86 2.0 301ms 11 -11.16041705144 -10.43 -5.21 2.7 319ms 12 -11.16041705145 -10.86 -5.79 3.0 387ms 13 -11.16041705145 -12.13 -6.36 2.6 328ms Computing bands along kpath: Γ -> M -> K -> Γ Diagonalising Hamiltonian kblocks: 100%|████████████████| Time: 0:00:06