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
using Plots
# 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
kpath = irrfbz_path(model; dim=2, space_group_number=13) # graphene space group number
bands = compute_bands(scfres, kpath; kline_density=20)
plot_bandstructure(bands)
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -11.15666063268 -0.60 6.1 203ms 2 -11.16021124379 -2.45 -1.30 1.0 115ms 3 -11.16039978162 -3.72 -2.33 2.1 149ms 4 -11.16041636799 -4.78 -3.20 2.1 143ms 5 -11.16041703709 -6.17 -3.39 2.3 140ms 6 -11.16041704742 -7.99 -3.54 1.1 112ms 7 -11.16041704979 -8.63 -3.84 1.0 108ms 8 -11.16041705072 -9.03 -4.15 1.1 112ms 9 -11.16041705116 -9.35 -4.51 1.4 118ms 10 -11.16041705137 -9.69 -4.91 1.4 118ms 11 -11.16041705143 -10.17 -5.33 1.9 128ms 12 -11.16041705145 -10.84 -5.76 1.9 147ms 13 -11.16041705145 -11.72 -6.25 2.0 132ms