Tutorial

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.oncvpsp3.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.505325524955                   -0.93    5.1   44.8ms
  2   -8.508293429362       -2.53       -1.79    1.0   21.8ms
  3   -8.508463491511       -3.77       -2.92    1.5   24.6ms
  4   -8.508487321466       -4.62       -3.09    3.0   34.3ms
  5   -8.508487484493       -6.79       -3.28    1.0   22.0ms
  6   -8.508487585885       -6.99       -4.63    1.0   21.7ms
  7   -8.508487591546       -8.25       -4.71    3.1   35.0ms
  8   -8.508487591701       -9.81       -5.53    1.0   22.1ms

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.0841999 
    AtomicLocal         -2.3554951
    AtomicNonlocal      1.3116733 
    Ewald               -8.3979253
    PspCorrection       0.3948635 
    Hartree             0.5559240 
    Xc                  -3.1017280

    total               -8.508487591701

Eigenvalues:

stack(scfres.eigenvalues)
7×8 Matrix{Float64}:
 -0.264799  -0.235353   -0.179075   …  -0.18931    -0.112074  -0.106308
  0.174322   0.0295303  -0.0825383     -0.026557   -0.112074  -0.106308
  0.174322   0.146146    0.129856       0.0347383   0.068654   0.0309527
  0.174322   0.146146    0.129856       0.125209    0.068654   0.0309527
  0.26715    0.247825    0.229556       0.26304     0.198144   0.329875
  0.26715    0.302404    0.297035   …   0.34955     0.198144   0.330501
  0.267152   0.302404    0.297035       0.362168    0.541872   0.356807

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)
Example block output

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.0348950058570192e-14, -2.145461385302733e-14, -2.1390170699434477e-14]
 [2.0749756207563516e-14, 2.1368163814917752e-14, 2.14839497274024e-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)
Example block output

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())
Example block output

Note that directly employing the scfres also works, but the results are much cruder:

plot_dos(scfres; temperature=1e-3, smearing=Smearing.FermiDirac())
Example block output