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.73387871262 -0.88 11.0 390ms
2 -36.61844390135 + -0.94 -1.42 1.0 94.9ms
3 +40.70055564237 + 1.89 -0.12 7.0 258ms
4 -36.42761476761 1.89 -1.09 7.0 270ms
5 -34.46647276878 + 0.29 -0.89 3.0 150ms
6 -36.41556085078 0.29 -1.26 4.0 172ms
7 -36.71851441242 -0.52 -1.76 3.0 139ms
8 -36.73869856578 -1.69 -2.04 1.0 102ms
9 -36.73996654366 -2.90 -2.02 2.0 134ms
10 -36.74065054476 -3.16 -2.22 2.0 121ms
11 -36.74118519154 -3.27 -2.13 2.0 124ms
12 -36.74198752246 -3.10 -2.58 2.0 111ms
13 -36.74223204586 -3.61 -2.69 1.0 103ms
14 -36.74219023883 + -4.38 -2.75 1.0 98.2ms
15 -36.74215998009 + -4.52 -2.71 1.0 104ms
16 -36.74234079120 -3.74 -2.93 3.0 127ms
17 -36.74241499210 -4.13 -3.07 3.0 130ms
18 -36.74201726850 + -3.40 -2.70 3.0 134ms
19 -36.74241649565 -3.40 -3.09 3.0 149ms
20 -36.74247343498 -4.24 -3.56 2.0 116ms
21 -36.74247966624 -5.21 -3.65 2.0 140ms
22 -36.74247955068 + -6.94 -3.89 2.0 120ms
23 -36.74248055606 -6.00 -4.12 1.0 100ms
24 -36.74247868678 + -5.73 -3.89 2.0 139ms
25 -36.74248065549 -5.71 -4.46 3.0 124ms
26 -36.74248037677 + -6.55 -4.30 3.0 150ms
27 -36.74248066576 -6.54 -4.85 2.0 115ms
28 -36.74248061207 + -7.27 -4.58 3.0 152ms
29 -36.74248065762 -7.34 -4.95 3.0 137ms
30 -36.74248066364 -8.22 -5.05 2.0 115ms
31 -36.74248066876 -8.29 -5.25 3.0 124ms
32 -36.74248067242 -8.44 -5.78 1.0 103ms
33 -36.74248067260 -9.74 -5.76 3.0 147ms
34 -36.74248067236 + -9.62 -5.71 2.0 123ms
35 -36.74248067249 -9.91 -5.84 3.0 135ms
36 -36.74248067257 -10.10 -5.97 2.0 121ms
37 -36.74248067264 -10.15 -6.19 1.0 97.2ms
38 -36.74248067266 -10.67 -6.30 3.0 143ms
39 -36.74248067267 -10.87 -6.42 2.0 134ms
40 -36.74248067268 -10.97 -6.70 2.0 118ms
41 -36.74248067268 + -11.78 -6.69 2.0 121ms
42 -36.74248067268 + -11.75 -6.74 3.0 143ms
43 -36.74248067268 -11.40 -7.48 2.0 107ms
44 -36.74248067268 + -12.08 -7.04 4.0 179ms
45 -36.74248067268 -12.05 -7.65 3.0 137ms
46 -36.74248067268 -13.67 -8.00 2.0 107ms
47 -36.74248067268 + -Inf -8.16 2.0 138ms
48 -36.74248067268 + -Inf -7.95 3.0 141ms
49 -36.74248067268 -14.15 -8.24 2.0 131ms
50 -36.74248067268 + -14.15 -8.80 2.0 110ms
51 -36.74248067268 -14.15 -8.68 3.0 157ms
52 -36.74248067268 + -Inf -9.06 2.0 107ms
53 -36.74248067268 + -Inf -9.00 3.0 143ms
54 -36.74248067268 -14.15 -9.42 2.0 115ms
55 -36.74248067268 + -14.15 -9.54 2.0 139ms
56 -36.74248067268 + -Inf -9.70 2.0 114ms
57 -36.74248067268 -14.15 -9.85 2.0 122ms
58 -36.74248067268 + -14.15 -9.86 1.0 97.9ms
59 -36.74248067268 -14.15 -10.33 1.0 98.1ms
60 -36.74248067268 + -Inf -10.20 3.0 150ms
61 -36.74248067268 + -14.15 -10.02 3.0 167ms
62 -36.74248067268 -14.15 -10.23 3.0 143ms
63 -36.74248067268 + -Inf -10.12 2.0 115ms
64 -36.74248067268 + -Inf -10.65 2.0 132ms
65 -36.74248067268 + -14.15 -10.94 2.0 106ms
66 -36.74248067268 + -Inf -10.65 3.0 153ms
67 -36.74248067268 + -Inf -11.20 3.0 143ms
68 -36.74248067268 + -13.85 -11.69 1.0 104ms
69 -36.74248067268 -13.67 -11.74 3.0 147ms
70 -36.74248067268 + -14.15 -12.06 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.73352333634 -0.88 13.0 383ms
2 -36.73993026933 -2.19 -1.37 1.0 93.7ms
3 -36.74099540517 -2.97 -1.95 3.0 140ms
4 -36.74208051075 -2.96 -2.06 1.0 92.7ms
5 -36.74226249831 -3.74 -2.63 4.0 111ms
6 -36.74243733365 -3.76 -2.68 3.0 150ms
7 -36.74246224554 -4.60 -3.04 1.0 96.8ms
8 -36.74247977237 -4.76 -3.52 4.0 120ms
9 -36.74247945913 + -6.50 -3.74 3.0 144ms
10 -36.74248057097 -5.95 -4.07 2.0 113ms
11 -36.74248062642 -7.26 -4.21 3.0 115ms
12 -36.74248050199 + -6.91 -4.36 3.0 129ms
13 -36.74248066224 -6.80 -4.76 2.0 112ms
14 -36.74248067210 -8.01 -5.19 3.0 147ms
15 -36.74248067045 + -8.78 -5.23 3.0 120ms
16 -36.74248067263 -8.66 -5.65 1.0 105ms
17 -36.74248067243 + -9.71 -5.81 3.0 185ms
18 -36.74248067263 -9.69 -6.11 1.0 106ms
19 -36.74248067268 -10.33 -6.69 3.0 142ms
20 -36.74248067268 -11.55 -6.98 3.0 150ms
21 -36.74248067268 -12.53 -7.15 2.0 110ms
22 -36.74248067268 -12.66 -7.39 2.0 120ms
23 -36.74248067268 -13.00 -7.84 3.0 127ms
24 -36.74248067268 + -13.85 -8.05 3.0 151ms
25 -36.74248067268 + -Inf -8.46 2.0 106ms
26 -36.74248067268 + -14.15 -8.62 3.0 151ms
27 -36.74248067268 -14.15 -8.97 1.0 100ms
28 -36.74248067268 -13.85 -9.54 3.0 128ms
29 -36.74248067268 + -13.85 -9.56 3.0 152ms
30 -36.74248067268 + -Inf -10.01 1.0 104ms
31 -36.74248067268 -13.85 -10.40 2.0 135ms
32 -36.74248067268 + -13.85 -10.52 4.0 132ms
33 -36.74248067268 + -Inf -10.92 2.0 111ms
34 -36.74248067268 + -Inf -11.08 3.0 146ms
35 -36.74248067268 -13.85 -11.56 1.0 100ms
36 -36.74248067268 + -13.55 -11.63 3.0 157ms
37 -36.74248067268 -13.85 -11.99 1.0 100ms
38 -36.74248067268 -14.15 -12.35 3.0 127msGiven 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.024489049660694The 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.24421118661838This 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.723584660079539Since 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).