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.73212182259 -0.88 12.0 1.36s
2 -36.69263937078 + -1.40 -1.54 1.0 245ms
3 +10.70980692788 + 1.68 -0.23 6.0 256ms
4 -36.64713926961 1.68 -1.18 6.0 205ms
5 -36.72276260934 -1.12 -1.63 2.0 115ms
6 -36.72242492148 + -3.47 -1.75 2.0 126ms
7 -36.28562789471 + -0.36 -1.24 4.0 140ms
8 -36.73745291634 -0.35 -2.06 3.0 261ms
9 -36.74190860507 -2.35 -2.25 2.0 123ms
10 -36.74074486860 + -2.93 -2.10 2.0 116ms
11 -36.74212388047 -2.86 -2.52 2.0 1.30s
12 -36.74229626498 -3.76 -2.76 2.0 113ms
13 -36.74247301352 -3.75 -3.20 1.0 91.8ms
14 -36.74246153739 + -4.94 -3.17 3.0 139ms
15 -36.74239026164 + -4.15 -3.06 3.0 139ms
16 -36.74241049345 -4.69 -3.12 3.0 135ms
17 -36.74247567383 -4.19 -3.62 2.0 109ms
18 -36.74247539079 + -6.55 -3.51 3.0 127ms
19 -36.74246258978 + -4.89 -3.42 2.0 118ms
20 -36.74246731018 -5.33 -3.50 3.0 132ms
21 -36.74248063939 -4.88 -4.36 2.0 107ms
22 -36.74248028115 + -6.45 -4.22 4.0 164ms
23 -36.74248055047 -6.57 -4.29 2.0 136ms
24 -36.74248065512 -6.98 -4.82 2.0 109ms
25 -36.74248066997 -7.83 -5.19 2.0 130ms
26 -36.74248066971 + -9.59 -5.04 3.0 117ms
27 -36.74248066140 + -8.08 -4.93 2.0 104ms
28 -36.74248064232 + -7.72 -4.80 3.0 137ms
29 -36.74248067229 -7.52 -5.55 3.0 147ms
30 -36.74248067265 -9.45 -5.91 2.0 104ms
31 -36.74248067253 + -9.92 -5.85 2.0 129ms
32 -36.74248066901 + -8.45 -5.28 4.0 156ms
33 -36.74248067265 -8.44 -6.08 3.0 144ms
34 -36.74248067043 + -8.65 -5.34 3.0 148ms
35 -36.74248067262 -8.66 -5.90 4.0 157ms
36 -36.74248067268 -10.27 -6.28 2.0 122ms
37 -36.74248067268 -11.70 -6.65 2.0 117ms
38 -36.74248067268 -11.39 -6.98 2.0 130ms
39 -36.74248067268 -12.32 -6.91 1.0 93.5ms
40 -36.74248067268 -12.50 -7.26 2.0 110ms
41 -36.74248067268 -12.79 -7.74 2.0 102ms
42 -36.74248067268 + -12.62 -7.30 3.0 152ms
43 -36.74248067268 -12.55 -7.93 3.0 132ms
44 -36.74248067268 + -14.15 -7.97 2.0 135ms
45 -36.74248067268 + -14.15 -8.38 2.0 110ms
46 -36.74248067268 + -13.67 -7.79 3.0 156ms
47 -36.74248067268 -13.55 -8.32 3.0 150ms
48 -36.74248067268 -14.15 -8.32 3.0 138ms
49 -36.74248067268 + -14.15 -8.55 3.0 122ms
50 -36.74248067268 + -Inf -8.37 3.0 132ms
51 -36.74248067268 -14.15 -9.03 2.0 116ms
52 -36.74248067268 + -Inf -9.15 3.0 141ms
53 -36.74248067268 + -14.15 -9.69 1.0 99.2ms
54 -36.74248067268 + -Inf -9.75 3.0 141ms
55 -36.74248067268 + -13.85 -9.78 2.0 116ms
56 -36.74248067268 -13.85 -10.03 2.0 103ms
57 -36.74248067268 + -Inf -10.27 1.0 99.0ms
58 -36.74248067268 + -Inf -10.20 2.0 129ms
59 -36.74248067268 + -Inf -10.50 2.0 105ms
60 -36.74248067268 -14.15 -10.01 3.0 150ms
61 -36.74248067268 -14.15 -10.49 3.0 145ms
62 -36.74248067268 + -14.15 -10.43 3.0 138ms
63 -36.74248067268 + -14.15 -11.29 2.0 110ms
64 -36.74248067268 + -Inf -11.10 3.0 163ms
65 -36.74248067268 + -Inf -11.36 2.0 122ms
66 -36.74248067268 + -Inf -11.51 2.0 110ms
67 -36.74248067268 -14.15 -11.59 3.0 121ms
68 -36.74248067268 + -14.15 -11.81 2.0 116ms
69 -36.74248067268 + -Inf -11.99 2.0 129ms
70 -36.74248067268 + -Inf -11.83 2.0 123ms
71 -36.74248067268 + -Inf -12.20 2.0 116ms
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.73300782732 -0.87 10.0 958ms
2 -36.73977462328 -2.17 -1.37 1.0 636ms
3 -36.73848028102 + -2.89 -1.58 4.0 166ms
4 -36.74228678228 -2.42 -2.27 1.0 93.0ms
5 -36.74234046808 -4.27 -2.48 8.0 148ms
6 -36.74242807230 -4.06 -2.47 2.0 119ms
7 -36.74246722788 -4.41 -2.99 1.0 95.1ms
8 -36.74247215587 -5.31 -2.99 1.0 101ms
9 -36.74247546141 -5.48 -3.29 1.0 95.8ms
10 -36.74248017695 -5.33 -3.92 2.0 124ms
11 -36.74248057650 -6.40 -4.14 2.0 133ms
12 -36.74248065694 -7.09 -4.34 1.0 102ms
13 -36.74248067095 -7.85 -4.92 2.0 110ms
14 -36.74248067215 -8.92 -5.09 3.0 145ms
15 -36.74248067249 -9.46 -5.20 1.0 96.8ms
16 -36.74248067222 + -9.57 -5.46 2.0 120ms
17 -36.74248067262 -9.39 -5.75 1.0 96.5ms
18 -36.74248067268 -10.29 -6.19 5.0 136ms
19 -36.74248067268 -11.32 -6.45 3.0 136ms
20 -36.74248067268 -11.82 -6.85 2.0 108ms
21 -36.74248067268 -12.92 -6.97 3.0 139ms
22 -36.74248067268 -13.55 -7.30 1.0 102ms
23 -36.74248067268 -13.15 -7.80 3.0 130ms
24 -36.74248067268 + -Inf -7.97 3.0 145ms
25 -36.74248067268 + -13.85 -8.32 2.0 102ms
26 -36.74248067268 -13.85 -8.74 3.0 136ms
27 -36.74248067268 + -14.15 -8.98 3.0 139ms
28 -36.74248067268 -14.15 -9.29 2.0 108ms
29 -36.74248067268 + -14.15 -9.70 3.0 139ms
30 -36.74248067268 -13.85 -9.90 3.0 130ms
31 -36.74248067268 + -13.85 -10.24 2.0 133ms
32 -36.74248067268 + -Inf -10.67 2.0 109ms
33 -36.74248067268 + -Inf -10.88 3.0 145ms
34 -36.74248067268 + -Inf -11.33 1.0 96.8ms
35 -36.74248067268 + -Inf -11.53 3.0 145ms
36 -36.74248067268 + -Inf -11.97 1.0 96.5ms
37 -36.74248067268 + -Inf -12.15 4.0 156ms
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.02448927747068The 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.24421142591458This 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.72359109391239Since 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).