In this section we will be concerned with modelling supercells of aluminium. When dealing with periodic problems there is no unique definition of the lattice: Clearly any duplication of the lattice along an axis is also a valid repetitive unit to describe exactly the same system. This is exactly what a supercell is: An n-fold repetition along one (or multiple) axes of the original lattice.
The following code achieves this for aluminium:
using AtomsBuilder
using DFTK
using LinearAlgebra
using Unitful
using UnitfulAtomic
using PseudoPotentialData
function aluminium_setup(repeat=1; Ecut=7.0, kgrid=[2, 2, 2])
# Use AtomsBuilder to setup aluminium cubic unit cell (4 Al atoms)
# with provided lattice constant, see AtomsBase integration for details.
unit_cell = bulk(:Al; a=7.65339u"bohr", cubic=true)
supercell = unit_cell * (repeat, 1, 1) # Make a supercell
# Select standard pseudodojo pseudopotentials, construct an LDA model, discretize
# Note: We disable symmetries explicitly here. Otherwise the problem sizes
# we are able to run on the CI are too simple to observe the numerical
# instabilities we want to trigger here.
pseudopotentials = PseudoFamily("dojo.nc.sr.lda.v0_4_1.standard.upf")
model = model_DFT(supercell; pseudopotentials, functionals=LDA(),
temperature=1e-3, symmetries=false)
PlaneWaveBasis(model; Ecut, kgrid)
end;
As expected we obtain the unit cell for repeat=1
:
aluminium_setup(1)
PlaneWaveBasis discretization: architecture : DFTK.CPU() num. mpi processes : 1 num. julia threads : 1 num. DFTK threads : 1 num. blas threads : 2 num. fft threads : 1 Ecut : 7.0 Ha fft_size : (24, 24, 24), 13824 total points kgrid : MonkhorstPack([2, 2, 2]) num. red. kpoints : 8 num. irred. kpoints : 8 Discretized Model(lda_x+lda_c_pw, 3D): lattice (in Bohr) : [7.65339 , 0 , 0 ] [0 , 7.65339 , 0 ] [0 , 0 , 7.65339 ] unit cell volume : 448.29 Bohr³ atoms : Al₄ atom potentials : ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") num. electrons : 12 spin polarization : none temperature : 0.001 Ha smearing : DFTK.Smearing.FermiDirac() terms : Kinetic() AtomicLocal() AtomicNonlocal() Ewald(nothing) PspCorrection() Hartree() Xc(lda_x, lda_c_pw) Entropy()
and 5-fold as large supercell with repeat=5
:
aluminium_setup(5)
PlaneWaveBasis discretization: architecture : DFTK.CPU() num. mpi processes : 1 num. julia threads : 1 num. DFTK threads : 1 num. blas threads : 2 num. fft threads : 1 Ecut : 7.0 Ha fft_size : (96, 24, 24), 55296 total points kgrid : MonkhorstPack([2, 2, 2]) num. red. kpoints : 8 num. irred. kpoints : 8 Discretized Model(lda_x+lda_c_pw, 3D): lattice (in Bohr) : [38.267 , 0 , 0 ] [0 , 7.65339 , 0 ] [0 , 0 , 7.65339 ] unit cell volume : 2241.5 Bohr³ atoms : Al₂₀ atom potentials : ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") ElementPsp(Al, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Al.upf") num. electrons : 60 spin polarization : none temperature : 0.001 Ha smearing : DFTK.Smearing.FermiDirac() terms : Kinetic() AtomicLocal() AtomicNonlocal() Ewald(nothing) PspCorrection() Hartree() Xc(lda_x, lda_c_pw) Entropy()
As we will see in this notebook the modelling of a system generally becomes harder if the system becomes larger.
For achieving the latter DFTK by default employs the LdosMixing
preconditioner [^HL2021] during the SCF iterations. This mixing approach is
completely parameter free, but still automatically adapts to the treated
system in order to efficiently prevent charge sloshing. As a result,
modelling aluminium slabs indeed takes roughly the same number of SCF iterations
irrespective of the supercell size:
[^HL2021]: M. F. Herbst and A. Levitt. Black-box inhomogeneous preconditioning for self-consistent field iterations in density functional theory. J. Phys. Cond. Matt 33 085503 (2021). ArXiv:2009.01665
self_consistent_field(aluminium_setup(1); tol=1e-4);
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -9.355192087429 -1.10 6.1 428ms 2 -9.356782960534 -2.80 -1.43 1.0 78.5ms 3 -9.357072155302 -3.54 -2.78 2.2 95.0ms 4 -9.357119492113 -4.32 -3.03 6.6 196ms 5 -9.357119865033 -6.43 -3.17 1.1 77.3ms 6 -9.357120085864 -6.66 -3.32 6.2 171ms 7 -9.357120222843 -6.86 -3.48 1.0 75.4ms 8 -9.357120282676 -7.22 -3.65 1.1 77.3ms 9 -9.357120314566 -7.50 -3.93 2.2 95.6ms 10 -9.357120321735 -8.14 -4.18 1.1 78.6ms
self_consistent_field(aluminium_setup(2); tol=1e-4);
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -18.74782587059 -0.97 6.4 525ms ┌ Warning: Eigensolver not converged │ n_iter = │ 8-element Vector{Int64}: │ 1 │ 1 │ 1 │ 6 │ 1 │ 7 │ 1 │ 1 └ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/scf/self_consistent_field.jl:76 2 -18.75924443948 -1.94 -1.32 2.4 356ms 3 -18.79217479336 -1.48 -2.21 3.9 566ms 4 -18.79261586414 -3.36 -2.52 5.1 401ms 5 -18.79262082936 -5.30 -3.18 2.9 257ms 6 -18.79262376235 -5.53 -3.58 6.2 377ms 7 -18.79262416471 -6.40 -3.96 2.4 243ms 8 -18.79262417627 -7.94 -4.46 2.0 241ms
self_consistent_field(aluminium_setup(4); tol=1e-4);
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -37.54938441249 -0.84 9.9 2.13s 2 -37.55822040752 -2.05 -1.22 3.4 903ms ┌ Warning: Eigensolver not converged │ n_iter = │ 8-element Vector{Int64}: │ 12 │ 6 │ 17 │ 14 │ 24 │ 6 │ 12 │ 11 └ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/scf/self_consistent_field.jl:76 3 -37.56446791305 -2.20 -2.15 12.8 1.70s 4 -37.56484227746 -3.43 -2.37 10.5 2.10s 5 -37.56497385310 -3.88 -3.03 3.1 974ms 6 -37.56497720341 -5.47 -3.40 6.8 1.30s 7 -37.56498444115 -5.14 -3.80 2.9 902ms 8 -37.56498518704 -6.13 -4.48 3.9 1.12s
When switching off explicitly the LdosMixing
, by selecting mixing=SimpleMixing()
,
the performance of number of required SCF steps starts to increase as we increase
the size of the modelled problem:
self_consistent_field(aluminium_setup(1); tol=1e-4, mixing=SimpleMixing());
n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -9.355228685331 -1.09 6.4 156ms 2 -9.356843293351 -2.79 -1.91 1.0 61.2ms 3 -9.357094899209 -3.60 -2.63 4.2 124ms 4 -9.357095495032 -6.22 -2.55 4.4 126ms 5 -9.357119904885 -4.61 -3.56 1.0 62.1ms 6 -9.357120213909 -6.51 -3.80 5.8 133ms 7 -9.357120321069 -6.97 -4.58 1.5 70.0ms
self_consistent_field(aluminium_setup(4); tol=1e-4, mixing=SimpleMixing());
┌ Warning: Eigensolver not converged │ n_iter = │ 8-element Vector{Int64}: │ 11 │ 3 │ 13 │ 13 │ 17 │ 13 │ 16 │ 2 └ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/scf/self_consistent_field.jl:76 n Energy log10(ΔE) log10(Δρ) Diag Δtime --- --------------- --------- --------- ---- ------ 1 -37.53784283860 -0.84 11.0 2.11s 2 -32.44804571825 + 0.71 -0.71 5.5 1.05s 3 +146.5738015668 + 2.25 0.08 17.6 3.41s 4 -24.39053013586 2.23 -0.37 8.0 2.64s 5 -36.06006869784 1.07 -0.93 4.9 1.47s 6 -36.92340356736 -0.06 -1.13 6.0 1.23s 7 -34.71926226884 + 0.34 -0.84 5.9 1.24s 8 -37.34474077953 0.42 -1.33 5.4 1.48s 9 -37.34347176973 + -2.90 -1.35 1.9 760ms 10 -37.43371685856 -1.04 -1.47 1.0 711ms 11 -37.54085560392 -0.97 -1.81 2.4 856ms 12 -37.56094929909 -1.70 -2.21 3.8 1.05s 13 -37.56481750469 -2.41 -2.69 3.8 1.23s 14 -37.56479665961 + -4.68 -2.88 4.8 1.57s 15 -37.56353153815 + -2.90 -2.46 4.8 1.18s 16 -37.56494517512 -2.85 -3.04 4.8 1.25s 17 -37.56497285634 -4.56 -3.43 2.9 942ms 18 -37.56498390410 -4.96 -3.87 2.5 943ms 19 -37.56485160180 + -3.88 -3.00 5.5 1.79s 20 -37.56498436759 -3.88 -3.99 6.4 1.38s 21 -37.56498505923 -6.16 -4.34 2.5 828ms
For completion let us note that the more traditional mixing=KerkerMixing()
approach would also help in this particular setting to obtain a constant
number of SCF iterations for an increasing system size (try it!). In contrast
to LdosMixing
, however, KerkerMixing
is only suitable to model bulk metallic
system (like the case we are considering here). When modelling metallic surfaces
or mixtures of metals and insulators, KerkerMixing
fails, while LdosMixing
still works well. See the Modelling a gallium arsenide surface example
or [^HL2021] for details. Due to the general applicability of LdosMixing
this
method is the default mixing approach in DFTK.