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

# 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.

# Load HGH pseudopotential for Silicon
Si = ElementPsp(:Si, psp=load_psp("hgh/lda/Si-q4"))

# Specify type and positions of atoms
atoms     = [Si, Si]
positions = [ones(3)/8, -ones(3)/8]

# 2. Select model and basis
model = model_LDA(lattice, atoms, positions)
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   -7.900350763154                   -0.70    4.8
  2   -7.904988914081       -2.33       -1.52    1.0   33.4ms
  3   -7.905174805247       -3.73       -2.53    1.2   35.3ms
  4   -7.905210574503       -4.45       -2.83    2.6   46.9ms
  5   -7.905211041488       -6.33       -2.93    1.1   35.4ms
  6   -7.905211518849       -6.32       -4.65    1.0   34.2ms
  7   -7.905211531025       -7.91       -4.54    3.1   53.4ms
  8   -7.905211531381       -9.45       -5.12    1.0   34.5ms

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.1020952 
    AtomicLocal         -2.1987824
    AtomicNonlocal      1.7296089 
    Ewald               -8.3979253
    PspCorrection       -0.2946254
    Hartree             0.5530384 
    Xc                  -2.3986210

    total               -7.905211531381

Eigenvalues:

hcat(scfres.eigenvalues...)
7×8 Matrix{Float64}:
 -0.176942  -0.14744   -0.091169   …  -0.101219   -0.0239768  -0.0184078
  0.261073   0.116915   0.0048252      0.0611645  -0.0239768  -0.0184078
  0.261073   0.23299    0.216734       0.121636    0.155532    0.117748
  0.261073   0.23299    0.216734       0.212134    0.155532    0.117748
  0.354532   0.335109   0.317102       0.350436    0.285692    0.417258
  0.354532   0.389829   0.384601   …   0.436926    0.285692    0.417391
  0.354532   0.389829   0.384601       0.449232    0.627537    0.443806

eigenvalues is an array (indexed by k-points) of arrays (indexed by eigenvalue number). The "splatting" operation ... calls hcat with all the inner arrays as arguments, which collects them into a matrix.

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 ...

hcat(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 ρ[iσ, ix, iy, iz], i.e. first in the spin component and then in the 3-dimensional real-space grid.

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: for instance compute a band structure

plot_bandstructure(scfres; kline_density=10)

or get the cartesian forces (in Hartree / Bohr)

compute_forces_cart(scfres)
2-element Vector{StaticArraysCore.SVector{3, Float64}}:
 [-6.346135937309126e-16, -2.6637050621587342e-17, 4.1764709737706414e-16]
 [-5.639359245861691e-17, -4.0738782468307493e-16, -5.028681177537508e-16]

As expected, they are numerically zero in this highly symmetric configuration.