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.73083932089 -0.88 12.0 335ms
2 -35.95701822885 + -0.11 -1.10 2.0 102ms
3 +149.2093593747 + 2.27 0.08 28.0 414ms
4 -24.65953719718 2.24 -0.47 11.0 297ms
5 -35.38270439849 1.03 -0.95 4.0 158ms
6 -26.68864623024 + 0.94 -0.55 5.0 161ms
7 -36.67866560802 1.00 -1.43 4.0 153ms
8 -36.73260835838 -1.27 -1.73 2.0 108ms
9 -36.74070454634 -2.09 -1.97 2.0 100ms
10 -36.73961230392 + -2.96 -2.03 2.0 133ms
11 -36.73548292417 + -2.38 -1.94 2.0 101ms
12 -36.73802470556 -2.59 -2.18 1.0 88.8ms
13 -36.74159071796 -2.45 -2.42 1.0 92.2ms
14 -36.74168329823 -4.03 -2.54 1.0 92.0ms
15 -36.74223329097 -3.26 -2.78 2.0 103ms
16 -36.73244011364 + -2.01 -2.07 3.0 140ms
17 -36.74245103537 -2.00 -2.78 4.0 144ms
18 -36.73649863955 + -2.23 -2.17 3.0 131ms
19 -36.74234355145 -2.23 -2.76 4.0 262ms
20 -36.74233585777 + -5.11 -2.79 2.0 119ms
21 -36.74250885770 -3.76 -3.32 1.0 688ms
22 -36.74228839863 + -3.66 -2.82 4.0 147ms
23 -36.74250672064 -3.66 -3.50 3.0 131ms
24 -36.74251441834 -5.11 -3.87 2.0 109ms
25 -36.74251445346 -7.45 -3.78 3.0 134ms
26 -36.74251454504 -7.04 -4.31 1.0 100ms
27 -36.74251475701 -6.67 -4.65 2.0 121ms
28 -36.74250519429 + -5.02 -3.58 4.0 152ms
29 -36.74251473355 -5.02 -4.67 4.0 155ms
30 -36.74251474260 -8.04 -4.73 2.0 114ms
31 -36.74251476181 -7.72 -4.92 2.0 106ms
32 -36.74251476438 -8.59 -5.02 2.0 126ms
33 -36.74251477259 -8.09 -5.46 2.0 104ms
34 -36.74251477227 + -9.49 -5.40 2.0 117ms
35 -36.74251476706 + -8.28 -5.11 3.0 127ms
36 -36.74251477230 -8.28 -5.57 2.0 107ms
37 -36.74251477303 -9.14 -6.21 2.0 107ms
38 -36.74251477299 + -10.46 -6.22 4.0 157ms
39 -36.74251477302 -10.56 -6.40 1.0 92.4ms
40 -36.74251477302 + -11.22 -6.37 2.0 113ms
41 -36.74251477301 + -12.32 -6.38 3.0 130ms
42 -36.74251477304 -10.67 -6.85 2.0 107ms
43 -36.74251477304 -12.60 -7.14 2.0 124ms
44 -36.74251477303 + -11.74 -6.83 3.0 131ms
45 -36.74251477304 -11.64 -7.64 3.0 122ms
46 -36.74251477304 + -13.19 -7.29 3.0 139ms
47 -36.74251477304 + -13.55 -7.42 3.0 124ms
48 -36.74251477304 -12.94 -7.96 2.0 109ms
49 -36.74251477304 -13.85 -8.13 3.0 126ms
50 -36.74251477304 + -14.15 -8.06 2.0 107ms
51 -36.74251477304 + -14.15 -8.56 1.0 90.2ms
52 -36.74251477304 -13.85 -8.27 3.0 139ms
53 -36.74251477304 + -14.15 -8.73 3.0 126ms
54 -36.74251477304 -14.15 -8.47 2.0 110ms
55 -36.74251477304 + -Inf -8.45 2.0 114ms
56 -36.74251477304 + -13.85 -8.41 3.0 131ms
57 -36.74251477304 -14.15 -9.66 3.0 126ms
58 -36.74251477304 + -14.15 -9.69 3.0 149ms
59 -36.74251477304 -14.15 -9.19 3.0 133ms
60 -36.74251477304 + -Inf -9.90 3.0 123ms
61 -36.74251477304 -14.15 -10.53 2.0 123ms
62 -36.74251477304 + -14.15 -10.39 3.0 135ms
63 -36.74251477304 -14.15 -10.87 2.0 103ms
64 -36.74251477304 + -Inf -10.87 3.0 127ms
65 -36.74251477304 -13.85 -11.31 2.0 104ms
66 -36.74251477304 + -13.55 -11.46 3.0 134ms
67 -36.74251477304 -14.15 -11.52 2.0 109ms
68 -36.74251477304 + -Inf -11.25 3.0 144ms
69 -36.74251477304 + -14.15 -11.73 3.0 133ms
70 -36.74251477304 + -Inf -12.09 2.0 105ms
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.73425804378 -0.88 11.0 333ms
2 -36.73998514526 -2.24 -1.36 1.0 89.7ms
3 -36.73699599240 + -2.52 -1.52 2.0 103ms
4 -36.74228794951 -2.28 -2.44 1.0 90.7ms
5 -36.74234243138 -4.26 -2.38 4.0 145ms
6 -36.74247414912 -3.88 -2.56 4.0 107ms
7 -36.74248750394 -4.87 -2.67 1.0 89.1ms
8 -36.74250874683 -4.67 -3.43 1.0 93.4ms
9 -36.74251432481 -5.25 -3.55 3.0 131ms
10 -36.74251468589 -6.44 -3.93 2.0 95.1ms
11 -36.74251468868 -8.55 -4.43 6.0 126ms
12 -36.74251473362 -7.35 -4.51 3.0 136ms
13 -36.74251477088 -7.43 -4.88 1.0 95.0ms
14 -36.74251477220 -8.88 -5.36 2.0 103ms
15 -36.74251477291 -9.15 -5.61 3.0 135ms
16 -36.74251477299 -10.05 -5.82 1.0 91.5ms
17 -36.74251477300 -11.07 -6.20 5.0 119ms
18 -36.74251477303 -10.49 -6.68 2.0 107ms
19 -36.74251477304 -11.56 -6.89 3.0 132ms
20 -36.74251477304 + -13.11 -7.17 1.0 94.8ms
21 -36.74251477304 -12.72 -7.41 3.0 134ms
22 -36.74251477304 + -Inf -7.69 1.0 92.3ms
23 -36.74251477304 + -Inf -8.22 3.0 106ms
24 -36.74251477304 -13.85 -8.12 3.0 134ms
25 -36.74251477304 + -Inf -8.44 1.0 95.3ms
26 -36.74251477304 + -Inf -8.95 2.0 98.8ms
27 -36.74251477304 + -14.15 -9.17 4.0 136ms
28 -36.74251477304 + -Inf -9.52 2.0 128ms
29 -36.74251477304 + -14.15 -9.79 2.0 101ms
30 -36.74251477304 -14.15 -10.02 2.0 105ms
31 -36.74251477304 -14.15 -10.20 2.0 109ms
32 -36.74251477304 + -Inf -10.65 1.0 91.4ms
33 -36.74251477304 + -14.15 -10.68 4.0 140ms
34 -36.74251477304 -14.15 -11.17 1.0 97.4ms
35 -36.74251477304 + -Inf -11.25 4.0 133ms
36 -36.74251477304 + -14.15 -11.39 1.0 95.7ms
37 -36.74251477304 + -14.15 -11.63 4.0 107ms
38 -36.74251477304 + -Inf -12.04 2.0 101ms
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
end
epsilon (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.02448777178898
The smallest eigenvalue is a bit more tricky to obtain, so we will just assume
λ_Simple_min = 0.952
0.952
This makes the condition number around 30:
cond_Simple = λ_Simple_max / λ_Simple_min
46.24420984431616
This 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.723597567775634
Since 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).