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.73225309124 -0.88 11.0 1.38s
2 -36.66297466507 + -1.16 -1.49 1.0 270ms
3 +28.94727239843 + 1.82 -0.15 7.0 276ms
4 -36.60049439522 1.82 -1.23 6.0 249ms
5 -36.07978972808 + -0.28 -1.14 3.0 145ms
6 -36.11316152652 -1.48 -1.15 5.0 190ms
7 -36.73001181128 -0.21 -1.84 3.0 142ms
8 -36.73992948706 -2.00 -2.12 1.0 104ms
9 -36.73953108960 + -3.40 -2.09 2.0 136ms
10 -36.74156838279 -2.69 -2.21 2.0 243ms
11 -36.74159893951 -4.51 -2.29 2.0 124ms
12 -36.74219382394 -3.23 -2.55 2.0 1.34s
13 -36.74240066800 -3.68 -2.86 2.0 112ms
14 -36.74245992769 -4.23 -2.96 1.0 101ms
15 -36.74016780233 + -2.64 -2.35 4.0 162ms
16 -36.74246466648 -2.64 -3.04 4.0 163ms
17 -36.74243414239 + -4.52 -3.07 2.0 115ms
18 -36.74201676544 + -3.38 -2.73 2.0 132ms
19 -36.74246688308 -3.35 -3.45 3.0 147ms
20 -36.74247819039 -4.95 -3.77 2.0 112ms
21 -36.74247768631 + -6.30 -3.77 3.0 146ms
22 -36.74248010832 -5.62 -4.08 1.0 100ms
23 -36.74248051807 -6.39 -4.14 2.0 126ms
24 -36.74248054371 -7.59 -4.38 2.0 124ms
25 -36.74248062102 -7.11 -4.64 2.0 158ms
26 -36.74248062387 -8.55 -4.65 2.0 112ms
27 -36.74248064815 -7.61 -4.79 2.0 128ms
28 -36.74248066841 -7.69 -5.18 2.0 117ms
29 -36.74248065409 + -7.84 -4.87 3.0 146ms
30 -36.74248067073 -7.78 -5.37 3.0 153ms
31 -36.74248067249 -8.75 -5.74 2.0 118ms
32 -36.74248067243 + -10.23 -5.84 3.0 127ms
33 -36.74248067262 -9.73 -6.11 2.0 137ms
34 -36.74248067267 -10.28 -6.51 1.0 99.4ms
35 -36.74248067268 -11.25 -6.66 3.0 161ms
36 -36.74248067268 -11.79 -6.56 2.0 125ms
37 -36.74248067268 -11.76 -7.14 2.0 113ms
38 -36.74248067268 + -12.19 -7.02 3.0 152ms
39 -36.74248067268 -12.11 -7.11 2.0 118ms
40 -36.74248067268 -13.37 -7.31 2.0 117ms
41 -36.74248067268 + -12.17 -7.02 2.0 134ms
42 -36.74248067268 -12.18 -7.38 3.0 145ms
43 -36.74248067268 -12.73 -7.68 2.0 123ms
44 -36.74248067268 -13.85 -7.95 1.0 100ms
45 -36.74248067268 + -13.25 -7.57 3.0 155ms
46 -36.74248067268 -13.25 -8.06 3.0 146ms
47 -36.74248067268 + -13.85 -8.68 2.0 123ms
48 -36.74248067268 -14.15 -8.54 3.0 153ms
49 -36.74248067268 + -Inf -8.82 2.0 117ms
50 -36.74248067268 -14.15 -9.17 2.0 123ms
51 -36.74248067268 -14.15 -8.95 3.0 149ms
52 -36.74248067268 -14.15 -9.15 3.0 151ms
53 -36.74248067268 + -13.85 -9.73 2.0 111ms
54 -36.74248067268 -14.15 -8.92 4.0 180ms
55 -36.74248067268 + -14.15 -9.72 4.0 194ms
56 -36.74248067268 + -Inf -9.91 2.0 117ms
57 -36.74248067268 -14.15 -10.34 2.0 113ms
58 -36.74248067268 + -14.15 -10.34 3.0 148ms
59 -36.74248067268 + -Inf -10.06 3.0 152ms
60 -36.74248067268 -13.85 -10.46 2.0 117ms
61 -36.74248067268 + -13.85 -11.06 2.0 135ms
62 -36.74248067268 + -14.15 -11.09 3.0 155ms
63 -36.74248067268 -14.15 -11.30 2.0 123ms
64 -36.74248067268 -14.15 -11.28 3.0 136ms
65 -36.74248067268 + -14.15 -11.52 2.0 128ms
66 -36.74248067268 + -Inf -11.60 2.0 120ms
67 -36.74248067268 + -Inf -11.49 3.0 146ms
68 -36.74248067268 + -Inf -12.01 2.0 124ms
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.73305526567 -0.88 12.0 1.02s
2 -36.73969266172 -2.18 -1.36 1.0 666ms
3 -36.74053717411 -3.07 -1.78 4.0 153ms
4 -36.74218178669 -2.78 -2.18 2.0 107ms
5 -36.74230601328 -3.91 -2.53 3.0 120ms
6 -36.74239997354 -4.03 -2.41 3.0 145ms
7 -36.74243883410 -4.41 -2.98 1.0 94.8ms
8 -36.74247611669 -4.43 -3.19 1.0 102ms
9 -36.74247607189 + -7.35 -3.44 2.0 102ms
10 -36.74248061691 -5.34 -4.25 1.0 103ms
11 -36.74248065495 -7.42 -4.37 9.0 185ms
12 -36.74248066100 -8.22 -4.55 1.0 98.2ms
13 -36.74248067034 -8.03 -4.86 2.0 113ms
14 -36.74248067241 -8.68 -5.22 1.0 98.4ms
15 -36.74248067262 -9.67 -5.65 3.0 132ms
16 -36.74248067266 -10.47 -5.83 3.0 145ms
17 -36.74248067268 -10.66 -6.45 1.0 97.8ms
18 -36.74248067268 + -12.94 -6.40 4.0 155ms
19 -36.74248067268 -11.55 -6.76 1.0 98.2ms
20 -36.74248067268 -12.69 -6.95 3.0 245ms
21 -36.74248067268 -12.77 -7.52 4.0 1.33s
22 -36.74248067268 + -Inf -7.47 2.0 135ms
23 -36.74248067268 + -13.85 -7.58 1.0 98.5ms
24 -36.74248067268 -13.55 -8.26 1.0 98.1ms
25 -36.74248067268 + -14.15 -8.32 6.0 166ms
26 -36.74248067268 + -Inf -8.72 1.0 98.0ms
27 -36.74248067268 + -Inf -8.94 3.0 143ms
28 -36.74248067268 -13.85 -9.23 1.0 97.6ms
29 -36.74248067268 + -Inf -9.58 1.0 98.4ms
30 -36.74248067268 + -13.67 -9.88 3.0 143ms
31 -36.74248067268 -14.15 -10.10 2.0 109ms
32 -36.74248067268 + -13.85 -10.38 1.0 100ms
33 -36.74248067268 -13.67 -10.63 3.0 156ms
34 -36.74248067268 + -14.15 -10.95 1.0 121ms
35 -36.74248067268 + -Inf -11.39 3.0 148ms
36 -36.74248067268 + -Inf -11.67 3.0 155ms
37 -36.74248067268 + -Inf -11.93 2.0 106ms
38 -36.74248067268 -14.15 -12.27 2.0 134ms
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
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.0244889801804The 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.244211113634876This 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.7235876536896875Since 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).