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.73380099725 -0.88 11.0 337ms
2 -36.71587417009 + -1.75 -1.59 1.0 90.2ms
3 -5.748644551249 + 1.49 -0.32 6.0 179ms
4 -36.58200345178 1.49 -1.12 6.0 187ms
5 -36.56995395716 + -1.92 -1.41 3.0 128ms
6 -36.23509437340 + -0.48 -1.21 4.0 143ms
7 -36.73308295903 -0.30 -1.83 3.0 127ms
8 -36.74192852480 -2.05 -2.26 1.0 92.4ms
9 -36.74165297800 + -3.56 -2.27 2.0 125ms
10 -36.74137117106 + -3.55 -2.26 2.0 115ms
11 -36.74226122486 -3.05 -2.64 2.0 107ms
12 -36.74245325080 -3.72 -2.95 2.0 100ms
13 -36.74247752858 -4.61 -3.30 2.0 123ms
14 -36.74247759615 -7.17 -3.24 2.0 107ms
15 -36.74248035922 -5.56 -3.87 1.0 93.8ms
16 -36.74240495630 + -4.12 -3.12 4.0 148ms
17 -36.74247462130 -4.16 -3.59 3.0 143ms
18 -36.74246472135 + -5.00 -3.43 3.0 134ms
19 -36.74247650518 -4.93 -3.72 3.0 133ms
20 -36.74248027089 -5.42 -4.17 3.0 124ms
21 -36.74248046928 -6.70 -4.23 2.0 110ms
22 -36.74248063062 -6.79 -4.48 2.0 108ms
23 -36.74248066528 -7.46 -4.99 2.0 126ms
24 -36.74248066684 -8.81 -4.99 2.0 107ms
25 -36.74248066906 -8.65 -5.04 2.0 100ms
26 -36.74248067019 -8.95 -5.16 1.0 90.7ms
27 -36.74248067133 -8.94 -5.47 1.0 93.9ms
28 -36.74248064354 + -7.56 -4.82 4.0 150ms
29 -36.74248067252 -7.54 -5.76 4.0 142ms
30 -36.74248067257 -10.26 -5.99 3.0 132ms
31 -36.74248067225 + -9.49 -5.72 3.0 134ms
32 -36.74248067259 -9.47 -6.07 3.0 214ms
33 -36.74248067268 -10.03 -6.66 2.0 755ms
34 -36.74248067268 + -11.35 -6.58 3.0 139ms
35 -36.74248067268 -11.33 -6.86 2.0 106ms
36 -36.74248067268 -12.51 -6.94 2.0 106ms
37 -36.74248067268 -12.05 -7.22 1.0 91.5ms
38 -36.74248067268 -12.50 -7.43 2.0 124ms
39 -36.74248067268 -13.55 -7.72 2.0 107ms
40 -36.74248067268 + -12.32 -7.19 3.0 137ms
41 -36.74248067268 -12.32 -7.68 3.0 131ms
42 -36.74248067268 + -14.15 -7.74 2.0 106ms
43 -36.74248067268 + -Inf -8.13 2.0 100ms
44 -36.74248067268 + -13.85 -7.83 3.0 134ms
45 -36.74248067268 -13.55 -7.91 3.0 132ms
46 -36.74248067268 + -13.45 -7.70 2.0 117ms
47 -36.74248067268 -13.37 -8.69 3.0 135ms
48 -36.74248067268 + -13.85 -8.46 3.0 125ms
49 -36.74248067268 -14.15 -8.32 3.0 136ms
50 -36.74248067268 -14.15 -9.08 3.0 123ms
51 -36.74248067268 + -14.15 -9.37 3.0 143ms
52 -36.74248067268 -14.15 -9.41 1.0 91.9ms
53 -36.74248067268 -14.15 -9.68 2.0 105ms
54 -36.74248067268 + -13.85 -9.64 2.0 109ms
55 -36.74248067268 + -Inf -10.20 1.0 95.4ms
56 -36.74248067268 -14.15 -9.88 3.0 138ms
57 -36.74248067268 + -Inf -9.50 4.0 158ms
58 -36.74248067268 + -14.15 -10.55 3.0 135ms
59 -36.74248067268 + -Inf -10.84 3.0 137ms
60 -36.74248067268 + -14.15 -10.93 2.0 110ms
61 -36.74248067268 -14.15 -10.58 3.0 124ms
62 -36.74248067268 + -Inf -11.01 3.0 134ms
63 -36.74248067268 -14.15 -11.38 2.0 110ms
64 -36.74248067268 + -14.15 -11.23 3.0 137ms
65 -36.74248067268 + -14.15 -11.69 1.0 91.1ms
66 -36.74248067268 -14.15 -12.06 4.0 133ms
while 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.73423390349 -0.88 11.0 338ms
2 -36.74024855333 -2.22 -1.36 1.0 91.8ms
3 -36.74072280249 -3.32 -1.81 3.0 127ms
4 -36.74217532766 -2.84 -2.14 1.0 92.7ms
5 -36.74231517824 -3.85 -2.68 2.0 97.6ms
6 -36.74243976912 -3.90 -2.55 6.0 146ms
7 -36.74246695074 -4.57 -3.05 1.0 90.9ms
8 -36.74247870596 -4.93 -3.22 2.0 103ms
9 -36.74247963592 -6.03 -3.41 2.0 106ms
10 -36.74248043285 -6.10 -4.00 1.0 91.5ms
11 -36.74248065724 -6.65 -4.20 4.0 140ms
12 -36.74248065143 + -8.24 -4.27 2.0 103ms
13 -36.74248067019 -7.73 -4.69 1.0 97.6ms
14 -36.74248067203 -8.74 -5.02 4.0 113ms
15 -36.74248067245 -9.38 -5.23 3.0 120ms
16 -36.74248067262 -9.76 -5.59 2.0 108ms
17 -36.74248067263 -11.07 -5.68 3.0 129ms
18 -36.74248067268 -10.33 -6.28 1.0 96.6ms
19 -36.74248067268 -12.53 -6.16 3.0 144ms
20 -36.74248067268 -11.52 -6.40 1.0 96.6ms
21 -36.74248067268 -12.27 -6.87 1.0 94.1ms
22 -36.74248067268 -12.81 -7.02 3.0 129ms
23 -36.74248067268 -13.07 -7.18 2.0 103ms
24 -36.74248067268 + -Inf -7.35 2.0 102ms
25 -36.74248067268 -14.15 -8.00 1.0 93.5ms
26 -36.74248067268 + -14.15 -7.87 4.0 149ms
27 -36.74248067268 -13.67 -8.28 2.0 103ms
28 -36.74248067268 + -Inf -8.54 2.0 102ms
29 -36.74248067268 + -14.15 -8.80 3.0 128ms
30 -36.74248067268 -14.15 -8.92 2.0 99.4ms
31 -36.74248067268 + -Inf -9.91 2.0 111ms
32 -36.74248067268 + -14.15 -9.76 5.0 164ms
33 -36.74248067268 -14.15 -10.02 1.0 97.6ms
34 -36.74248067268 + -Inf -10.46 2.0 121ms
35 -36.74248067268 + -Inf -10.61 2.0 102ms
36 -36.74248067268 + -Inf -11.08 2.0 107ms
37 -36.74248067268 + -13.85 -11.62 3.0 126ms
38 -36.74248067268 -13.85 -12.01 5.0 146ms
Given 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
end
epsilon (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.02448898072273
The smallest eigenvalue is a bit more tricky to obtain, so we will just assume
λ_Simple_min = 0.952
0.952
This makes the condition number around 30:
cond_Simple = λ_Simple_max / λ_Simple_min
46.24421111420455
This 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.723592816485761
Since 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).