Analysing SCF convergence
The goal of this example is to explain the differing convergence behaviour of SCF algorithms depending on the choice of the mixing. For this we look at the eigenpairs of the Jacobian governing the SCF convergence, that is
\[1 - α P^{-1} \varepsilon^\dagger \qquad \text{with} \qquad \varepsilon^\dagger = (1-\chi_0 K).\]
where $α$ is the damping $P^{-1}$ is the mixing preconditioner (e.g. KerkerMixing, LdosMixing) and $\varepsilon^\dagger$ is the dielectric operator.
We thus investigate the largest and smallest eigenvalues of $(P^{-1} \varepsilon^\dagger)$ and $\varepsilon^\dagger$. The ratio of largest to smallest eigenvalue of this operator is the condition number
\[\kappa = \frac{\lambda_\text{max}}{\lambda_\text{min}},\]
which can be related to the rate of convergence of the SCF. The smaller the condition number, the faster the convergence. For more details on SCF methods, see Self-consistent field methods.
For our investigation we consider a crude aluminium setup:
using AtomsBuilder
using DFTK
system_Al = bulk(:Al; cubic=true) * (4, 1, 1)FlexibleSystem(Al₁₆, periodicity = TTT):
cell_vectors : [ 16.2 0 0;
0 4.05 0;
0 0 4.05]u"Å"
and we discretise:
using PseudoPotentialData
pseudopotentials = PseudoFamily("dojo.nc.sr.lda.v0_4_1.standard.upf")
model_Al = model_DFT(system_Al; functionals=LDA(), temperature=1e-3,
symmetries=false, pseudopotentials)
basis_Al = PlaneWaveBasis(model_Al; Ecut=7, kgrid=[1, 1, 1]);On aluminium (a metal) already for moderate system sizes (like the 8 layers we consider here) the convergence without mixing / preconditioner is slow:
# Note: DFTK uses the self-adapting LdosMixing() by default, so to truly disable
# any preconditioning, we need to supply `mixing=SimpleMixing()` explicitly.
scfres_Al = self_consistent_field(basis_Al; tol=1e-12, mixing=SimpleMixing());n Energy log10(ΔE) log10(Δρ) Diag Δtime
--- --------------- --------- --------- ---- ------
1 -36.73442476688 -0.88 11.0 363ms
2 -36.67539210808 + -1.23 -1.51 1.0 91.8ms
3 +20.44752390907 + 1.76 -0.19 7.0 228ms
4 -36.69737662141 1.76 -1.33 6.0 231ms
5 -36.70665540002 -2.03 -1.60 2.0 123ms
6 -36.23154121737 + -0.32 -1.20 3.0 136ms
7 -36.65967129984 -0.37 -1.58 4.0 160ms
8 -36.74144445821 -1.09 -2.25 2.0 107ms
9 -36.74075217616 + -3.16 -2.11 3.0 162ms
10 -36.74181479063 -2.97 -2.31 2.0 115ms
11 -36.74187096043 -4.25 -2.39 2.0 130ms
12 -36.74228301357 -3.39 -2.77 2.0 111ms
13 -36.74245145249 -3.77 -3.08 1.0 103ms
14 -36.74246869452 -4.76 -3.06 3.0 141ms
15 -36.74216152314 + -3.51 -2.80 3.0 149ms
16 -36.74246010708 -3.52 -3.34 3.0 142ms
17 -36.74246058703 -6.32 -3.17 2.0 120ms
18 -36.74246004584 + -6.27 -3.38 2.0 121ms
19 -36.74247666405 -4.78 -3.69 2.0 114ms
20 -36.74247829766 -5.79 -3.85 3.0 142ms
21 -36.74248046191 -5.66 -4.02 2.0 132ms
22 -36.74248000542 + -6.34 -4.02 2.0 122ms
23 -36.74248061779 -6.21 -4.38 2.0 108ms
24 -36.74248065167 -7.47 -4.57 2.0 125ms
25 -36.74248066525 -7.87 -4.73 1.0 97.5ms
26 -36.74248066868 -8.46 -5.13 1.0 96.8ms
27 -36.74248066811 + -9.24 -5.16 2.0 137ms
28 -36.74248065558 + -7.90 -4.93 3.0 131ms
29 -36.74248067225 -7.78 -5.50 3.0 137ms
30 -36.74248066735 + -8.31 -5.16 3.0 151ms
31 -36.74248067254 -8.28 -5.66 3.0 147ms
32 -36.74248066126 + -7.95 -5.01 3.0 149ms
33 -36.74248067249 -7.95 -5.80 3.0 156ms
34 -36.74248067259 -9.98 -5.92 2.0 115ms
35 -36.74248067264 -10.27 -6.04 2.0 118ms
36 -36.74248067268 -10.48 -6.56 2.0 105ms
37 -36.74248067268 -11.50 -6.83 3.0 146ms
38 -36.74248067268 -12.03 -6.92 2.0 117ms
39 -36.74248067268 -12.46 -7.32 1.0 104ms
40 -36.74248067268 + -11.38 -6.73 3.0 149ms
41 -36.74248067268 -11.40 -7.29 3.0 146ms
42 -36.74248067268 -12.56 -7.99 2.0 109ms
43 -36.74248067268 + -Inf -7.94 3.0 149ms
44 -36.74248067268 + -Inf -7.94 1.0 104ms
45 -36.74248067268 -13.85 -8.20 2.0 114ms
46 -36.74248067268 -13.85 -8.21 3.0 142ms
47 -36.74248067268 + -13.85 -8.21 3.0 131ms
48 -36.74248067268 + -Inf -8.99 2.0 113ms
49 -36.74248067268 + -Inf -8.82 4.0 167ms
50 -36.74248067268 + -Inf -9.12 2.0 119ms
51 -36.74248067268 + -Inf -9.08 2.0 124ms
52 -36.74248067268 + -Inf -9.50 2.0 119ms
53 -36.74248067268 + -Inf -9.64 2.0 133ms
54 -36.74248067268 -14.15 -9.46 3.0 151ms
55 -36.74248067268 + -14.15 -9.86 2.0 114ms
56 -36.74248067268 -13.85 -10.12 2.0 115ms
57 -36.74248067268 + -13.85 -10.27 3.0 129ms
58 -36.74248067268 + -Inf -10.38 2.0 105ms
59 -36.74248067268 + -14.15 -10.54 2.0 114ms
60 -36.74248067268 -14.15 -10.58 3.0 133ms
61 -36.74248067268 + -Inf -11.24 2.0 120ms
62 -36.74248067268 + -14.15 -10.99 3.0 151ms
63 -36.74248067268 -14.15 -11.50 2.0 132ms
64 -36.74248067268 + -Inf -11.07 3.0 150ms
65 -36.74248067268 -14.15 -11.64 3.0 147ms
66 -36.74248067268 + -13.85 -11.34 3.0 134ms
67 -36.74248067268 + -Inf -11.28 3.0 146ms
68 -36.74248067268 -13.67 -11.77 3.0 141ms
69 -36.74248067268 + -13.85 -11.76 2.0 125ms
70 -36.74248067268 -14.15 -12.62 1.0 104mswhile when using the Kerker preconditioner it is much faster:
scfres_Al = self_consistent_field(basis_Al; tol=1e-12, mixing=KerkerMixing());n Energy log10(ΔE) log10(Δρ) Diag Δtime
--- --------------- --------- --------- ---- ------
1 -36.73182974311 -0.88 13.0 359ms
2 -36.73911866212 -2.14 -1.36 1.0 92.9ms
3 -36.73662422680 + -2.60 -1.52 3.0 130ms
4 -36.74222152464 -2.25 -2.29 1.0 93.5ms
5 -36.74216616487 + -4.26 -2.37 6.0 160ms
6 -36.74239391683 -3.64 -2.39 1.0 94.1ms
7 -36.74243321965 -4.41 -2.72 1.0 94.7ms
8 -36.74247526962 -4.38 -3.10 1.0 101ms
9 -36.74247889481 -5.44 -3.42 4.0 112ms
10 -36.74248046328 -5.80 -3.91 3.0 119ms
11 -36.74248047499 -7.93 -4.12 3.0 133ms
12 -36.74248066966 -6.71 -4.56 2.0 137ms
13 -36.74248063159 + -7.42 -4.69 1.0 98.4ms
14 -36.74248067209 -7.39 -5.32 1.0 103ms
15 -36.74248067216 -10.14 -5.55 4.0 153ms
16 -36.74248067262 -9.33 -5.75 1.0 103ms
17 -36.74248067266 -10.47 -6.18 2.0 113ms
18 -36.74248067268 -10.66 -6.52 2.0 137ms
19 -36.74248067268 -11.85 -6.77 4.0 115ms
20 -36.74248067268 + -12.38 -6.88 2.0 120ms
21 -36.74248067268 -12.33 -7.33 1.0 98.3ms
22 -36.74248067268 -13.07 -7.47 3.0 147ms
23 -36.74248067268 -13.85 -7.87 1.0 97.6ms
24 -36.74248067268 + -14.15 -7.87 3.0 133ms
25 -36.74248067268 -13.85 -8.22 1.0 98.0ms
26 -36.74248067268 + -14.15 -8.77 2.0 131ms
27 -36.74248067268 + -13.85 -9.19 2.0 119ms
28 -36.74248067268 + -Inf -9.33 3.0 138ms
29 -36.74248067268 + -Inf -9.80 5.0 125ms
30 -36.74248067268 -13.67 -9.87 2.0 132ms
31 -36.74248067268 + -13.67 -10.33 1.0 107ms
32 -36.74248067268 -14.15 -10.69 2.0 125ms
33 -36.74248067268 + -Inf -10.79 3.0 143ms
34 -36.74248067268 + -14.15 -11.13 1.0 97.1ms
35 -36.74248067268 -13.67 -11.47 2.0 108ms
36 -36.74248067268 + -13.85 -11.74 3.0 130ms
37 -36.74248067268 + -13.85 -12.04 2.0 136msGiven this scfres_Al we construct functions representing $\varepsilon^\dagger$ and $P^{-1}$:
# Function, which applies P^{-1} for the case of KerkerMixing
Pinv_Kerker(δρ) = DFTK.mix_density(KerkerMixing(), basis_Al, δρ)
# Function which applies ε† = 1 - χ0 K
function epsilon(δρ)
δV = apply_kernel(basis_Al, δρ; ρ=scfres_Al.ρ)
χ0δV = apply_χ0(scfres_Al, δV).δρ
δρ - χ0δV
endepsilon (generic function with 1 method)With these functions available we can now compute the desired eigenvalues. For simplicity we only consider the first few largest ones.
using KrylovKit
λ_Simple, X_Simple = eigsolve(epsilon, randn(size(scfres_Al.ρ)), 3, :LM;
tol=1e-3, eager=true, verbosity=2)
λ_Simple_max = maximum(real.(λ_Simple))44.024488980333956The smallest eigenvalue is a bit more tricky to obtain, so we will just assume
λ_Simple_min = 0.9520.952This makes the condition number around 30:
cond_Simple = λ_Simple_max / λ_Simple_min46.244211113796176This does not sound large compared to the condition numbers you might know from linear systems.
However, this is sufficient to cause a notable slowdown, which would be even more pronounced if we did not use Anderson, since we also would need to drastically reduce the damping (try it!).
Having computed the eigenvalues of the dielectric matrix we can now also look at the eigenmodes, which are responsible for the bad convergence behaviour. The largest eigenmode for example:
using Statistics
using Plots
mode_xy = mean(real.(X_Simple[1]), dims=3)[:, :, 1, 1] # Average along z axis
heatmap(mode_xy', c=:RdBu_11, aspect_ratio=1, grid=false,
legend=false, clim=(-0.006, 0.006))This mode can be physically interpreted as the reason why this SCF converges slowly. For example in this case it displays a displacement of electron density from the centre to the extremal parts of the unit cell. This phenomenon is called charge-sloshing.
We repeat the exercise for the Kerker-preconditioned dielectric operator:
λ_Kerker, X_Kerker = eigsolve(Pinv_Kerker ∘ epsilon,
randn(size(scfres_Al.ρ)), 3, :LM;
tol=1e-3, eager=true, verbosity=2)
mode_xy = mean(real.(X_Kerker[1]), dims=3)[:, :, 1, 1] # Average along z axis
heatmap(mode_xy', c=:RdBu_11, aspect_ratio=1, grid=false,
legend=false, clim=(-0.006, 0.006))Clearly the charge-sloshing mode is no longer dominating.
The largest eigenvalue is now
maximum(real.(λ_Kerker))4.7235844642762Since the smallest eigenvalue in this case remains of similar size (it is now around 0.8), this implies that the conditioning improves noticeably when KerkerMixing is used.
Note: Since LdosMixing requires solving a linear system at each application of $P^{-1}$, determining the eigenvalues of $P^{-1} \varepsilon^\dagger$ is slightly more expensive and thus not shown. The results are similar to KerkerMixing, however.
We could repeat the exercise for an insulating system (e.g. a Helium chain). In this case you would notice that the condition number without mixing is actually smaller than the condition number with Kerker mixing. In other words employing Kerker mixing makes the convergence worse. A closer investigation of the eigenvalues shows that Kerker mixing reduces the smallest eigenvalue of the dielectric operator this time, while keeping the largest value unchanged. Overall the conditioning thus workens.
Takeaways:
- For metals the conditioning of the dielectric matrix increases steeply with system size.
- The Kerker preconditioner tames this and makes SCFs on large metallic systems feasible by keeping the condition number of order 1.
- For insulating systems the best approach is to not use any mixing.
- The ideal mixing strongly depends on the dielectric properties of system which is studied (metal versus insulator versus semiconductor).