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.73356980262 -0.88 11.0 351ms
2 -36.62109064469 + -0.95 -1.42 1.0 86.5ms
3 +38.32251163025 + 1.87 -0.12 8.0 255ms
4 -36.44475922495 1.87 -1.07 7.0 230ms
5 -33.71008934292 + 0.44 -0.83 4.0 146ms
6 -36.62029989054 0.46 -1.43 4.0 150ms
7 -36.72140988354 -1.00 -1.71 2.0 112ms
8 -36.73720171591 -1.80 -2.07 2.0 102ms
9 -36.73870769914 -2.82 -1.97 2.0 129ms
10 -36.74066667598 -2.71 -2.10 2.0 124ms
11 -36.74035384004 + -3.50 -2.23 2.0 115ms
12 -36.74202584294 -2.78 -2.53 1.0 90.5ms
13 -36.74247329049 -3.35 -2.89 2.0 109ms
14 -36.74246673367 + -5.18 -2.97 2.0 127ms
15 -36.74064524346 + -2.74 -2.41 3.0 133ms
16 -36.74244416761 -2.74 -2.89 4.0 145ms
17 -36.73690352773 + -2.26 -2.18 3.0 143ms
18 -36.74228995785 -2.27 -2.74 4.0 148ms
19 -36.74249274849 -3.69 -3.31 3.0 125ms
20 -36.74246768693 + -4.60 -3.20 3.0 128ms
21 -36.74251353764 -4.34 -3.76 2.0 109ms
22 -36.74250364570 + -5.00 -3.51 3.0 137ms
23 -36.74251411611 -4.98 -3.98 2.0 106ms
24 -36.74251463974 -6.28 -4.14 2.0 119ms
25 -36.74251466543 -7.59 -4.33 1.0 91.2ms
26 -36.74251476771 -6.99 -5.04 2.0 104ms
27 -36.74251470327 + -7.19 -4.57 3.0 162ms
28 -36.74251476824 -7.19 -4.99 3.0 139ms
29 -36.74251476791 + -9.48 -5.15 2.0 105ms
30 -36.74251477166 -8.43 -5.31 2.0 110ms
31 -36.74251476462 + -8.15 -5.07 2.0 119ms
32 -36.74251477204 -8.13 -5.37 3.0 126ms
33 -36.74251477291 -9.06 -5.91 1.0 94.7ms
34 -36.74251476867 + -8.37 -5.24 4.0 150ms
35 -36.74251477286 -8.38 -5.85 3.0 141ms
36 -36.74251477300 -9.87 -6.20 1.0 94.8ms
37 -36.74251477302 -10.58 -6.37 3.0 138ms
38 -36.74251477304 -10.87 -6.64 3.0 114ms
39 -36.74251477303 + -11.97 -6.74 2.0 127ms
40 -36.74251477303 + -11.09 -6.48 3.0 134ms
41 -36.74251477304 -10.98 -7.24 2.0 118ms
42 -36.74251477304 -13.37 -7.30 3.0 142ms
43 -36.74251477304 + -12.28 -7.01 2.0 121ms
44 -36.74251477304 + -12.81 -7.01 3.0 132ms
45 -36.74251477304 -12.25 -7.32 3.0 138ms
46 -36.74251477304 -12.60 -7.87 2.0 106ms
47 -36.74251477304 + -12.72 -7.33 3.0 143ms
48 -36.74251477304 -12.70 -8.13 3.0 132ms
49 -36.74251477304 -14.15 -8.25 2.0 120ms
50 -36.74251477304 -14.15 -8.62 2.0 110ms
51 -36.74251477304 + -Inf -8.59 2.0 124ms
52 -36.74251477304 + -13.85 -8.96 2.0 110ms
53 -36.74251477304 -13.85 -9.14 2.0 106ms
54 -36.74251477304 + -Inf -9.37 1.0 97.3ms
55 -36.74251477304 + -14.15 -9.50 2.0 129ms
56 -36.74251477304 -13.85 -9.39 2.0 106ms
57 -36.74251477304 + -14.15 -9.82 2.0 110ms
58 -36.74251477304 + -Inf -9.85 2.0 101ms
59 -36.74251477304 + -13.85 -9.96 2.0 116ms
60 -36.74251477304 -13.85 -10.24 2.0 104ms
61 -36.74251477304 + -14.15 -9.78 3.0 143ms
62 -36.74251477304 -14.15 -10.23 3.0 131ms
63 -36.74251477304 + -Inf -10.85 2.0 112ms
64 -36.74251477304 + -14.15 -10.94 3.0 130ms
65 -36.74251477304 + -14.15 -10.44 3.0 139ms
66 -36.74251477304 -13.85 -10.64 3.0 134ms
67 -36.74251477304 + -Inf -11.14 2.0 110ms
68 -36.74251477304 + -Inf -10.97 2.0 117ms
69 -36.74251477304 + -14.15 -11.17 2.0 102ms
70 -36.74251477304 -14.15 -11.53 1.0 91.2ms
71 -36.74251477304 + -Inf -11.68 2.0 116ms
72 -36.74251477304 + -Inf -11.58 3.0 157ms
73 -36.74251477304 + -Inf -11.82 2.0 106ms
74 -36.74251477304 + -Inf -12.21 2.0 104ms
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.73177713461 -0.88 11.0 332ms
2 -36.73958079622 -2.11 -1.36 1.0 91.5ms
3 -36.74036340368 -3.11 -1.75 3.0 113ms
4 -36.74211190967 -2.76 -2.14 5.0 109ms
5 -36.74234156758 -3.64 -2.54 2.0 137ms
6 -36.74244811320 -3.97 -2.50 2.0 126ms
7 -36.74247939206 -4.50 -2.83 1.0 89.8ms
8 -36.74250905318 -4.53 -3.17 1.0 94.2ms
9 -36.74251196477 -5.54 -3.31 1.0 90.9ms
10 -36.74251453254 -5.59 -4.01 2.0 100ms
11 -36.74251467417 -6.85 -4.26 4.0 143ms
12 -36.74251476871 -7.02 -4.49 2.0 105ms
13 -36.74251476473 + -8.40 -4.81 3.0 107ms
14 -36.74251476789 -8.50 -4.96 3.0 141ms
15 -36.74251477016 -8.64 -5.03 1.0 111ms
16 -36.74251477284 -8.57 -5.50 2.0 102ms
17 -36.74251477280 + -10.44 -5.57 3.0 122ms
18 -36.74251477301 -9.67 -6.00 1.0 96.7ms
19 -36.74251477304 -10.69 -6.55 3.0 139ms
20 -36.74251477304 -11.80 -6.92 5.0 120ms
21 -36.74251477304 -12.97 -7.11 2.0 129ms
22 -36.74251477304 -13.37 -7.32 1.0 96.5ms
23 -36.74251477304 -13.37 -7.60 2.0 97.8ms
24 -36.74251477304 -13.37 -8.03 3.0 140ms
25 -36.74251477304 + -14.15 -8.16 2.0 106ms
26 -36.74251477304 + -Inf -8.60 2.0 99.2ms
27 -36.74251477304 + -Inf -8.56 3.0 135ms
28 -36.74251477304 + -13.85 -8.92 1.0 96.5ms
29 -36.74251477304 -13.85 -9.38 3.0 130ms
30 -36.74251477304 + -Inf -9.45 2.0 126ms
31 -36.74251477304 + -14.15 -9.75 1.0 96.6ms
32 -36.74251477304 + -Inf -10.06 3.0 131ms
33 -36.74251477304 -14.15 -10.54 2.0 129ms
34 -36.74251477304 + -Inf -10.90 2.0 130ms
35 -36.74251477304 + -13.85 -11.31 1.0 96.8ms
36 -36.74251477304 + -Inf -11.29 3.0 139ms
37 -36.74251477304 + -Inf -11.57 1.0 92.8ms
38 -36.74251477304 + -Inf -11.93 2.0 111ms
39 -36.74251477304 -13.85 -11.97 3.0 140ms
40 -36.74251477304 + -Inf -12.43 1.0 94.1ms
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.02448804311329
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.244210129320685
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.723586142597708
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).