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.73518594142 -0.88 11.0 1.46s
2 -36.65849067417 + -1.12 -1.47 1.0 256ms
3 +23.30962431441 + 1.78 -0.17 8.0 292ms
4 -36.50799677103 1.78 -1.06 6.0 255ms
5 -36.66864547901 -0.79 -1.46 3.0 156ms
6 -35.82030519069 + -0.07 -1.08 4.0 162ms
7 -36.69345320805 -0.06 -1.66 3.0 142ms
8 -36.73837188946 -1.35 -2.07 2.0 107ms
9 -36.74076206206 -2.62 -2.14 2.0 133ms
10 -36.74052227789 + -3.62 -2.05 2.0 132ms
11 -36.74238507250 -2.73 -2.55 2.0 113ms
12 -36.74244686592 -4.21 -2.67 2.0 118ms
13 -36.74238217509 + -4.19 -2.79 1.0 97.1ms
14 -36.74242559379 -4.36 -2.98 1.0 103ms
15 -36.74079006503 + -2.79 -2.45 3.0 153ms
16 -36.74238940768 -2.80 -2.98 3.0 148ms
17 -36.74229647558 + -4.03 -2.84 2.0 126ms
18 -36.74162035793 + -3.17 -2.59 3.0 149ms
19 -36.74247447858 -3.07 -3.58 3.0 148ms
20 -36.74247479291 -6.50 -3.59 3.0 138ms
21 -36.74248046975 -5.25 -4.20 2.0 120ms
22 -36.74247959991 + -6.06 -3.75 2.0 134ms
23 -36.74248044640 -6.07 -4.31 3.0 149ms
24 -36.74248065511 -6.68 -4.80 2.0 115ms
25 -36.74248066784 -7.90 -5.04 3.0 149ms
26 -36.74248066805 -9.68 -4.88 2.0 134ms
27 -36.74248067239 -8.36 -5.62 2.0 116ms
28 -36.74248066586 + -8.18 -5.15 3.0 158ms
29 -36.74248067178 -8.23 -5.55 3.0 144ms
30 -36.74248067261 -9.08 -5.86 2.0 119ms
31 -36.74248067252 + -10.08 -5.92 2.0 116ms
32 -36.74248067266 -9.85 -6.32 2.0 120ms
33 -36.74248067264 + -10.68 -6.14 3.0 157ms
34 -36.74248067267 -10.63 -6.27 2.0 125ms
35 -36.74248067268 -10.81 -6.81 2.0 118ms
36 -36.74248067268 + -12.42 -6.97 1.0 97.3ms
37 -36.74248067268 -12.29 -7.14 3.0 147ms
38 -36.74248067268 -12.38 -7.26 2.0 115ms
39 -36.74248067268 + -12.43 -7.17 2.0 222ms
40 -36.74248067268 -13.45 -7.23 3.0 139ms
41 -36.74248067268 -12.42 -7.67 2.0 1.16s
42 -36.74248067268 -13.55 -7.80 2.0 134ms
43 -36.74248067268 -13.85 -8.00 2.0 117ms
44 -36.74248067268 + -Inf -8.36 2.0 110ms
45 -36.74248067268 -13.85 -8.57 2.0 134ms
46 -36.74248067268 + -Inf -8.74 2.0 123ms
47 -36.74248067268 + -13.55 -9.07 1.0 112ms
48 -36.74248067268 -13.85 -8.80 3.0 169ms
49 -36.74248067268 + -Inf -9.57 3.0 152ms
50 -36.74248067268 -14.15 -9.59 3.0 169ms
51 -36.74248067268 -14.15 -9.78 2.0 131ms
52 -36.74248067268 + -Inf -9.82 2.0 136ms
53 -36.74248067268 + -13.85 -10.10 1.0 111ms
54 -36.74248067268 + -Inf -9.75 3.0 169ms
55 -36.74248067268 + -Inf -9.91 2.0 142ms
56 -36.74248067268 + -Inf -10.54 2.0 137ms
57 -36.74248067268 -13.85 -11.00 2.0 129ms
58 -36.74248067268 -13.85 -11.18 3.0 132ms
59 -36.74248067268 + -13.55 -10.84 3.0 145ms
60 -36.74248067268 + -Inf -11.15 3.0 137ms
61 -36.74248067268 -14.15 -11.67 2.0 110ms
62 -36.74248067268 + -14.15 -11.60 2.0 127ms
63 -36.74248067268 -13.85 -12.14 1.0 93.9ms
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.73057429241 -0.88 11.0 923ms
2 -36.73803469764 -2.13 -1.36 1.0 656ms
3 -36.73635346047 + -2.77 -1.57 3.0 131ms
4 -36.74190987571 -2.26 -2.15 1.0 91.0ms
5 -36.74156566425 + -3.46 -2.47 4.0 121ms
6 -36.74209601641 -3.28 -2.42 3.0 152ms
7 -36.74232009667 -3.65 -2.81 1.0 104ms
8 -36.74246878511 -3.83 -3.25 3.0 114ms
9 -36.74248034326 -4.94 -3.57 2.0 130ms
10 -36.74248009472 + -6.60 -3.62 1.0 92.5ms
11 -36.74248054674 -6.34 -3.91 1.0 92.9ms
12 -36.74248063974 -7.03 -4.60 5.0 120ms
13 -36.74248066890 -7.54 -4.71 3.0 137ms
14 -36.74248067010 -8.92 -4.80 1.0 100ms
15 -36.74248067231 -8.66 -5.21 1.0 98.8ms
16 -36.74248067239 -10.10 -5.53 2.0 110ms
17 -36.74248067262 -9.63 -5.74 3.0 122ms
18 -36.74248067265 -10.53 -5.94 2.0 100ms
19 -36.74248067268 -10.54 -6.23 3.0 119ms
20 -36.74248067268 -11.91 -6.70 3.0 108ms
21 -36.74248067268 -12.22 -6.70 3.0 149ms
22 -36.74248067268 -13.37 -7.05 1.0 100ms
23 -36.74248067268 -12.67 -7.25 3.0 105ms
24 -36.74248067268 -13.15 -7.59 2.0 134ms
25 -36.74248067268 + -Inf -7.76 2.0 100ms
26 -36.74248067268 -14.15 -8.06 2.0 115ms
27 -36.74248067268 + -Inf -8.26 3.0 123ms
28 -36.74248067268 -14.15 -8.72 2.0 128ms
29 -36.74248067268 + -14.15 -9.01 2.0 105ms
30 -36.74248067268 -13.85 -9.36 3.0 116ms
31 -36.74248067268 + -Inf -9.45 3.0 138ms
32 -36.74248067268 + -14.15 -10.03 1.0 97.9ms
33 -36.74248067268 + -Inf -10.28 4.0 141ms
34 -36.74248067268 + -Inf -10.47 3.0 139ms
35 -36.74248067268 -13.85 -10.93 2.0 110ms
36 -36.74248067268 + -14.15 -11.19 2.0 133ms
37 -36.74248067268 + -Inf -11.50 1.0 100ms
38 -36.74248067268 + -13.85 -11.72 2.0 104ms
39 -36.74248067268 -14.15 -12.07 5.0 130ms
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.02448898018825The 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.244211113643125This 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.723583649655636Since 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).