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.15656881832 -0.60 6.0 2 -11.16015822545 -2.44 -1.30 1.0 219ms 3 -11.16039367848 -3.63 -2.32 2.0 230ms 4 -11.16041668817 -4.64 -3.29 3.0 280ms 5 -11.16041704566 -6.45 -3.48 3.1 277ms 6 -11.16041704960 -8.40 -3.65 1.3 209ms 7 -11.16041705106 -8.83 -3.94 1.9 206ms 8 -11.16041705135 -9.54 -4.37 2.0 222ms 9 -11.16041705140 -10.30 -4.79 2.1 223ms 10 -11.16041705143 -10.61 -5.05 2.1 242ms 11 -11.16041705144 -10.84 -5.31 2.6 238ms 12 -11.16041705145 -11.12 -5.72 2.7 261ms 13 -11.16041705145 -11.88 -6.27 2.7 256ms Computing bands along kpath: Γ -> M -> K -> Γ Diagonalising Hamiltonian kblocks: 100%|████████████████| Time: 0:00:04