DFTK is a Julia package for playing with plane-wave density-functional theory algorithms. In its basic formulation it solves periodic Kohn-Sham equations.
This document provides an overview of the structure of the code and how to access basic information about calculations. Basic familiarity with the concepts of plane-wave density functional theory is assumed throughout. Feel free to take a look at the Periodic problems or the Introductory resources chapters for some introductory material on the topic.
Convergence parameters in the documentation
We use rough parameters in order to be able to automatically generate this documentation very quickly. Therefore results are far from converged. Tighter thresholds and larger grids should be used for more realistic results.
For our discussion we will use the classic example of computing the LDA ground state of the silicon crystal. Performing such a calculation roughly proceeds in three steps.
using DFTK
using Plots
using Unitful
using UnitfulAtomic
using PseudoPotentialData
# 1. Define lattice and atomic positions
a = 5.431u"angstrom" # Silicon lattice constant
lattice = a / 2 * [[0 1 1.]; # Silicon lattice vectors
[1 0 1.]; # specified column by column
[1 1 0.]];
By default, all numbers passed as arguments are assumed to be in atomic units. Quantities such as temperature, energy cutoffs, lattice vectors, and the k-point grid spacing can optionally be annotated with Unitful units, which are automatically converted to the atomic units used internally. For more details, see the Unitful package documentation and the UnitfulAtomic.jl package.
We use a pseudodojo pseudopotential
(see PseudoPotentialData
for more details on PseudoFamily
):
pd_lda_family = PseudoFamily("dojo.nc.sr.lda.v0_4_1.standard.upf")
Si = ElementPsp(:Si, pd_lda_family)
# Specify type and positions of atoms
atoms = [Si, Si]
positions = [ones(3)/8, -ones(3)/8]
2-element Vector{Vector{Float64}}: [0.125, 0.125, 0.125] [-0.125, -0.125, -0.125]
Note that DFTK supports a few other ways to supply atomistic structures, see for example the sections on AtomsBase integration and Input and output formats for details.
# 2. Select model and basis
model = model_DFT(lattice, atoms, positions; functionals=LDA())
kgrid = [4, 4, 4] # k-point grid (Regular Monkhorst-Pack grid)
Ecut = 7 # kinetic energy cutoff
# Ecut = 190.5u"eV" # Could also use eV or other energy-compatible units
basis = PlaneWaveBasis(model; Ecut, kgrid)
# Note the implicit passing of keyword arguments here:
# this is equivalent to PlaneWaveBasis(model; Ecut=Ecut, kgrid=kgrid)
# 3. Run the SCF procedure to obtain the ground state
scfres = self_consistent_field(basis, tol=1e-5);
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -8.505304500322 -0.93 5.4 4.12s 2 -8.508288857562 -2.53 -1.79 1.0 4.48s 3 -8.508467218562 -3.75 -2.93 1.8 53.9ms 4 -8.508483060644 -4.80 -3.26 2.6 56.6ms 5 -8.508483189659 -6.89 -3.62 1.1 23.2ms 6 -8.508483221233 -7.50 -4.97 1.4 28.2ms 7 -8.508483223795 -8.59 -5.05 3.2 36.2ms
That's it! Now you can get various quantities from the result of the SCF. For instance, the different components of the energy:
scfres.energies
Energy breakdown (in Ha): Kinetic 3.0841920 AtomicLocal -2.3554760 AtomicNonlocal 1.3116657 Ewald -8.3979253 PspCorrection 0.3948681 Hartree 0.5559173 Xc -3.1017249 total -8.508483223795
Eigenvalues:
stack(scfres.eigenvalues)
7×8 Matrix{Float64}: -0.264798 -0.235351 -0.179073 … -0.189308 -0.112071 -0.106305 0.174324 0.0295324 -0.0825367 -0.0265547 -0.112071 -0.106305 0.174324 0.146148 0.129859 0.0347397 0.068656 0.0309541 0.174324 0.146148 0.129859 0.125211 0.068656 0.0309541 0.267152 0.247827 0.229557 0.263041 0.198144 0.329876 0.267152 0.302405 0.297036 … 0.349552 0.198144 0.329955 0.267152 0.302405 0.297036 0.362161 0.541866 0.356808
eigenvalues
is an array (indexed by k-points) of arrays (indexed by
eigenvalue number).
The resulting matrix is 7 (number of computed eigenvalues) by 8
(number of irreducible k-points). There are 7 eigenvalues per
k-point because there are 4 occupied states in the system (4 valence
electrons per silicon atom, two atoms per unit cell, and paired
spins), and the eigensolver gives itself some breathing room by
computing some extra states (see the bands
argument to
self_consistent_field
as well as the AdaptiveBands
documentation).
There are only 8 k-points (instead of 4x4x4) because symmetry has been used to reduce the
amount of computations to just the irreducible k-points (see
Crystal symmetries
for details).
We can check the occupations ...
stack(scfres.occupation)
7×8 Matrix{Float64}: 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 2.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
... and density, where we use that the density objects in DFTK are indexed as ρ[ix, iy, iz, iσ], i.e. first in the 3-dimensional real-space grid and then in the spin component.
rvecs = collect(r_vectors(basis))[:, 1, 1] # slice along the x axis
x = [r[1] for r in rvecs] # only keep the x coordinate
plot(x, scfres.ρ[1, :, 1, 1], label="", xlabel="x", ylabel="ρ", marker=2)
We can also perform various postprocessing steps: We can get the Cartesian forces (in Hartree / Bohr):
compute_forces_cart(scfres)
2-element Vector{StaticArraysCore.SVector{3, Float64}}: [-2.051359950112395e-14, -2.1339749420919406e-14, -2.081491976567319e-14] [2.078754985381231e-14, 2.1373452336197004e-14, 2.1073053722947937e-14]
As expected, they are numerically zero in this highly symmetric configuration. We could also compute a band structure,
plot_bandstructure(scfres; kline_density=10)
┌ Warning: Calling plot_bandstructure without first computing the band data is deprecated and will be removed in the next minor version bump. └ @ DFTKPlotsExt ~/work/DFTK.jl/DFTK.jl/ext/DFTKPlotsExt.jl:27
or plot a density of states, for which we increase the kgrid a bit to get smoother plots:
bands = compute_bands(scfres, MonkhorstPack(6, 6, 6))
plot_dos(bands; temperature=1e-3, smearing=Smearing.FermiDirac())
Note that directly employing the scfres
also works, but the results
are much cruder:
plot_dos(scfres; temperature=1e-3, smearing=Smearing.FermiDirac())
Where to go from here
- Background on DFT:
- Periodic problems,
- Introduction to density-functional theory,
- Self-consistent field methods
- Running calculations:
- Temperature and metallic systems
- Performing a convergence study
- Geometry optimization
- Tips and tricks:
- Using DFTK on compute clusters,
- Using DFTK on GPUs,
- Saving SCF results on disk and SCF checkpoints