Creating and modelling metallic supercells
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](@ref) 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.oncvpsp3.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/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/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/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/Al.upf")
ElementPsp(Al, "/home/runner/.julia/artifacts/1ea71a84cf375286564538a9cab789991f4bf1f4/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.
- This sounds like a trivial statement as per se the cost per SCF step increases as the system (and thus $N$) gets larger.
- But there is more to it: If one is not careful also the number of SCF iterations increases as the system gets larger.
- The aim of a proper computational treatment of such supercells is therefore to ensure that the number of SCF iterations remains constant when the system size increases.
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:
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.355452001316 -1.10 7.1 203ms
2 -9.356740378513 -2.89 -1.41 1.0 73.2ms
3 -9.356980608821 -3.62 -2.75 2.0 86.3ms
4 -9.357011588681 -4.51 -3.04 6.1 162ms
5 -9.357011986790 -6.40 -3.18 1.1 82.5ms
6 -9.357012274302 -6.54 -3.35 1.0 71.2ms
7 -9.357012436818 -6.79 -3.51 2.2 91.2ms
8 -9.357012496255 -7.23 -3.64 1.0 71.1ms
9 -9.357012540070 -7.36 -3.93 1.0 71.1ms
10 -9.357012547385 -8.14 -4.33 1.4 109ms
self_consistent_field(aluminium_setup(2); tol=1e-4);
n Energy log10(ΔE) log10(Δρ) Diag Δtime
--- --------------- --------- --------- ---- ------
1 -18.78032013252 -0.97 7.4 450ms
2 -18.79177798356 -1.94 -1.30 2.1 241ms
3 -18.79230767398 -3.28 -2.81 3.5 318ms
4 -18.79240731885 -4.00 -2.99 6.8 442ms
5 -18.79240843152 -5.95 -3.24 1.0 200ms
6 -18.79240856123 -6.89 -3.59 1.1 202ms
7 -18.79240859552 -7.46 -4.26 3.2 264ms
self_consistent_field(aluminium_setup(4); tol=1e-4);
n Energy log10(ΔE) log10(Δρ) Diag Δtime
--- --------------- --------- --------- ---- ------
1 -37.54975450160 -0.84 9.6 2.46s
2 -37.56025723248 -1.98 -1.19 1.0 731ms
3 -37.56180300605 -2.81 -2.53 5.5 1.15s
4 -37.56454772790 -2.56 -2.43 8.2 1.58s
5 -37.56455082803 -5.51 -2.49 6.8 1.06s
6 -37.56455373870 -5.54 -3.45 1.0 775ms
7 -37.56455399725 -6.59 -3.79 10.2 1.45s
8 -37.56455406771 -7.15 -4.42 4.4 1.04s
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.355273912164 -1.09 6.9 168ms
2 -9.356769936545 -2.83 -1.91 1.0 63.3ms
3 -9.356999604879 -3.64 -2.78 4.2 128ms
4 -9.356992359312 + -5.14 -2.61 3.5 117ms
5 -9.357012222878 -4.70 -3.56 1.1 72.5ms
6 -9.357012461129 -6.62 -3.85 3.8 109ms
7 -9.357012544119 -7.08 -4.50 1.5 73.6ms
self_consistent_field(aluminium_setup(4); tol=1e-4, mixing=SimpleMixing());
n Energy log10(ΔE) log10(Δρ) Diag Δtime
--- --------------- --------- --------- ---- ------
1 -37.55495287891 -0.84 9.5 1.76s
2 -37.55041392366 + -2.34 -1.59 1.0 644ms
3 -25.97570173438 + 1.06 -0.54 7.5 1.55s
4 -37.47910897366 1.06 -1.41 8.4 1.70s
5 -37.55061215931 -1.15 -1.93 2.5 888ms
6 -37.36757672793 + -0.74 -1.41 5.0 1.88s
7 -37.56096307569 -0.71 -2.24 5.2 1.13s
8 -37.56421231487 -2.49 -2.42 3.1 987ms
9 -37.56436562066 -3.81 -2.59 2.5 863ms
10 -37.56435540311 + -4.99 -2.83 1.8 790ms
11 -37.56454136223 -3.73 -3.42 1.8 787ms
12 -37.56455098137 -5.02 -3.59 4.5 1.23s
13 -37.56455365412 -5.57 -3.93 3.2 895ms
14 -37.56455400179 -6.46 -4.29 2.9 891ms
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.