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.73367354029 -0.88 13.0 1.32s
2 -36.65210708663 + -1.09 -1.43 1.0 251ms
3 +30.02911426867 + 1.82 -0.15 7.0 229ms
4 -36.49123554160 1.82 -0.99 6.0 275ms
5 -35.93681240428 + -0.26 -1.08 3.0 154ms
6 -35.63323406704 + -0.52 -1.04 6.0 213ms
7 -36.72825020214 0.04 -1.77 3.0 141ms
8 -36.73676424744 -2.07 -2.04 2.0 120ms
9 -36.73431523456 + -2.61 -1.87 2.0 131ms
10 -36.74126728239 -2.16 -2.17 2.0 130ms
11 -36.74161904473 -3.45 -2.37 1.0 100ms
12 -36.74227887921 -3.18 -2.66 2.0 117ms
13 -36.74241521601 -3.87 -2.75 1.0 98.2ms
14 -36.74246264536 -4.32 -2.98 1.0 104ms
15 -36.74170061885 + -3.12 -2.53 3.0 145ms
16 -36.74245995823 -3.12 -3.12 3.0 144ms
17 -36.74245446477 + -5.26 -3.25 2.0 120ms
18 -36.74147502504 + -3.01 -2.57 3.0 160ms
19 -36.74246203315 -3.01 -3.40 4.0 159ms
20 -36.74247372594 -4.93 -3.60 2.0 130ms
21 -36.74247940031 -5.25 -3.93 2.0 142ms
22 -36.74248056521 -5.93 -4.32 2.0 113ms
23 -36.74248044928 + -6.94 -4.18 2.0 137ms
24 -36.74248063580 -6.73 -4.67 2.0 141ms
25 -36.74248066240 -7.58 -4.75 2.0 119ms
26 -36.74248067244 -8.00 -5.46 1.0 101ms
27 -36.74248066547 + -8.16 -5.05 3.0 166ms
28 -36.74248067246 -8.16 -5.50 3.0 142ms
29 -36.74248067261 -9.83 -5.84 1.0 100ms
30 -36.74248067266 -10.28 -6.09 2.0 131ms
31 -36.74248067165 + -8.99 -5.55 3.0 146ms
32 -36.74248067247 -9.09 -5.89 3.0 152ms
33 -36.74248067267 -9.68 -6.47 3.0 127ms
34 -36.74248067266 + -11.00 -6.38 3.0 150ms
35 -36.74248067268 -10.74 -7.14 2.0 116ms
36 -36.74248067268 + -12.75 -7.05 2.0 138ms
37 -36.74248067268 -12.47 -7.40 2.0 121ms
38 -36.74248067268 -13.67 -7.57 1.0 100ms
39 -36.74248067268 -13.85 -7.72 2.0 122ms
40 -36.74248067268 + -Inf -7.78 1.0 97.2ms
41 -36.74248067268 -13.67 -8.11 1.0 104ms
42 -36.74248067268 -13.85 -8.40 2.0 115ms
43 -36.74248067268 + -13.67 -8.07 3.0 149ms
44 -36.74248067268 -14.15 -8.88 3.0 147ms
45 -36.74248067268 -13.85 -8.45 3.0 148ms
46 -36.74248067268 + -13.85 -8.92 3.0 150ms
47 -36.74248067268 + -Inf -8.98 3.0 135ms
48 -36.74248067268 + -Inf -9.13 2.0 114ms
49 -36.74248067268 + -Inf -9.85 1.0 236ms
50 -36.74248067268 + -Inf -9.80 3.0 156ms
51 -36.74248067268 -14.15 -10.11 1.0 1.21s
52 -36.74248067268 + -14.15 -9.95 3.0 146ms
53 -36.74248067268 + -Inf -9.96 3.0 146ms
54 -36.74248067268 + -Inf -10.60 2.0 114ms
55 -36.74248067268 + -Inf -10.37 3.0 141ms
56 -36.74248067268 + -Inf -10.61 3.0 143ms
57 -36.74248067268 + -Inf -10.67 2.0 114ms
58 -36.74248067268 + -13.85 -10.46 3.0 138ms
59 -36.74248067268 -13.85 -11.17 3.0 140ms
60 -36.74248067268 + -Inf -11.49 2.0 150ms
61 -36.74248067268 -14.15 -11.90 2.0 125ms
62 -36.74248067268 + -Inf -11.68 3.0 164ms
63 -36.74248067268 + -14.15 -11.83 2.0 146ms
64 -36.74248067268 + -Inf -12.16 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.73484780105 -0.88 11.0 937ms
2 -36.74038233356 -2.26 -1.36 1.0 451ms
3 -36.74163076856 -2.90 -2.00 4.0 131ms
4 -36.74219856957 -3.25 -2.10 2.0 127ms
5 -36.74238468273 -3.73 -2.52 1.0 95.1ms
6 -36.74243994931 -4.26 -2.53 2.0 133ms
7 -36.74245896366 -4.72 -3.15 2.0 103ms
8 -36.74247991006 -4.68 -3.50 6.0 132ms
9 -36.74248028850 -6.42 -3.64 3.0 123ms
10 -36.74248053256 -6.61 -4.04 1.0 107ms
11 -36.74248063192 -7.00 -4.26 2.0 153ms
12 -36.74248066936 -7.43 -4.66 1.0 97.0ms
13 -36.74248067102 -8.78 -4.83 3.0 130ms
14 -36.74248067245 -8.84 -5.23 3.0 110ms
15 -36.74248067264 -9.72 -5.67 3.0 143ms
16 -36.74248067261 + -10.45 -5.68 3.0 113ms
17 -36.74248067267 -10.19 -6.27 2.0 110ms
18 -36.74248067268 -11.48 -6.30 4.0 148ms
19 -36.74248067268 -11.34 -6.60 1.0 97.6ms
20 -36.74248067268 -11.92 -7.16 3.0 129ms
21 -36.74248067268 -12.73 -7.23 4.0 149ms
22 -36.74248067268 -13.67 -7.54 2.0 104ms
23 -36.74248067268 -13.85 -8.03 3.0 126ms
24 -36.74248067268 + -13.85 -8.03 3.0 142ms
25 -36.74248067268 -14.15 -8.42 1.0 97.8ms
26 -36.74248067268 + -Inf -9.09 3.0 126ms
27 -36.74248067268 + -Inf -9.03 4.0 162ms
28 -36.74248067268 -14.15 -9.37 1.0 97.1ms
29 -36.74248067268 + -Inf -9.74 3.0 127ms
30 -36.74248067268 + -14.15 -10.13 3.0 142ms
31 -36.74248067268 + -Inf -10.42 2.0 103ms
32 -36.74248067268 + -Inf -10.60 3.0 145ms
33 -36.74248067268 + -14.15 -11.08 1.0 98.9ms
34 -36.74248067268 + -Inf -11.28 3.0 146ms
35 -36.74248067268 -14.15 -11.44 2.0 124ms
36 -36.74248067268 + -Inf -11.87 1.0 98.0ms
37 -36.74248067268 + -13.85 -12.25 3.0 148ms
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.024489071947855The 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.24421121002926This 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.723583155238512Since 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).