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.73568991388 -0.88 11.0 344ms
2 -36.73858453437 -2.54 -1.69 1.0 85.1ms
3 -32.87883802005 + 0.59 -0.77 4.0 157ms
4 -36.72131992200 0.58 -1.74 4.0 151ms
5 -36.67427425087 + -1.33 -1.60 2.0 107ms
6 -36.43121724782 + -0.61 -1.32 3.0 126ms
7 -36.73279732969 -0.52 -1.84 3.0 128ms
8 -36.73992821210 -2.15 -2.19 2.0 102ms
9 -36.74214538899 -2.65 -2.28 2.0 112ms
10 -36.74240559283 -3.58 -2.63 2.0 103ms
11 -36.74248494957 -4.10 -2.84 2.0 121ms
12 -36.74249400571 -5.04 -3.12 1.0 93.0ms
13 -36.74248629957 + -5.11 -3.28 1.0 93.1ms
14 -36.74249442500 -5.09 -3.34 2.0 102ms
15 -36.74243014245 + -4.19 -3.04 3.0 132ms
16 -36.74250250660 -4.14 -3.49 2.0 241ms
17 -36.74251177329 -5.03 -3.79 2.0 731ms
18 -36.74251450354 -5.56 -4.08 2.0 123ms
19 -36.74251370709 + -6.10 -3.89 3.0 127ms
20 -36.74251306053 + -6.19 -3.93 3.0 131ms
21 -36.74251323653 -6.75 -3.93 3.0 145ms
22 -36.74251470814 -5.83 -4.55 2.0 117ms
23 -36.74251467521 + -7.48 -4.50 3.0 131ms
24 -36.74251476058 -7.07 -4.94 1.0 89.5ms
25 -36.74251476963 -8.04 -5.13 3.0 126ms
26 -36.74251477273 -8.51 -5.54 2.0 124ms
27 -36.74251477288 -9.82 -5.85 1.0 90.0ms
28 -36.74251477288 + -11.25 -5.73 2.0 125ms
29 -36.74251477289 -10.87 -5.82 2.0 106ms
30 -36.74251477279 + -10.00 -5.83 2.0 119ms
31 -36.74251477301 -9.66 -6.15 2.0 102ms
32 -36.74251477303 -10.72 -6.37 2.0 104ms
33 -36.74251477274 + -9.55 -5.79 3.0 142ms
34 -36.74251477298 -9.62 -6.14 3.0 141ms
35 -36.74251477302 -10.36 -6.38 3.0 125ms
36 -36.74251477302 + -11.07 -6.37 3.0 126ms
37 -36.74251477303 -10.71 -6.73 2.0 107ms
38 -36.74251477303 + -11.98 -6.70 3.0 135ms
39 -36.74251477304 -11.50 -7.19 1.0 90.0ms
40 -36.74251477304 -12.87 -7.27 3.0 132ms
41 -36.74251477304 -12.56 -7.58 2.0 99.3ms
42 -36.74251477304 + -12.28 -7.15 3.0 140ms
43 -36.74251477304 -12.25 -8.19 3.0 132ms
44 -36.74251477304 + -Inf -8.21 3.0 131ms
45 -36.74251477304 + -14.15 -8.16 2.0 117ms
46 -36.74251477304 + -13.45 -7.76 3.0 131ms
47 -36.74251477304 -13.55 -8.67 3.0 132ms
48 -36.74251477304 -14.15 -8.69 3.0 126ms
49 -36.74251477304 + -14.15 -8.79 2.0 109ms
50 -36.74251477304 -13.85 -9.30 1.0 89.8ms
51 -36.74251477304 + -13.85 -8.98 3.0 140ms
52 -36.74251477304 + -Inf -9.44 3.0 127ms
53 -36.74251477304 -14.15 -9.62 2.0 103ms
54 -36.74251477304 + -14.15 -9.36 3.0 131ms
55 -36.74251477304 + -Inf -10.06 3.0 112ms
56 -36.74251477304 -13.85 -10.08 3.0 133ms
57 -36.74251477304 + -13.85 -10.49 1.0 92.3ms
58 -36.74251477304 + -Inf -10.73 3.0 131ms
59 -36.74251477304 -13.85 -10.68 2.0 114ms
60 -36.74251477304 + -Inf -10.99 2.0 100ms
61 -36.74251477304 + -Inf -11.38 2.0 108ms
62 -36.74251477304 + -Inf -11.28 2.0 125ms
63 -36.74251477304 + -14.15 -11.27 2.0 110ms
64 -36.74251477304 -14.15 -11.01 3.0 141ms
65 -36.74251477304 + -14.15 -11.51 3.0 123ms
66 -36.74251477304 + -14.15 -11.53 3.0 132ms
67 -36.74251477304 -14.15 -12.14 2.0 108ms
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.73268115704 -0.88 11.0 325ms
2 -36.73927073774 -2.18 -1.35 1.0 89.6ms
3 -36.73567475021 + -2.44 -1.40 3.0 128ms
4 -36.74232764964 -2.18 -2.24 1.0 89.9ms
5 -36.74240338296 -4.12 -2.40 5.0 125ms
6 -36.74240627539 -5.54 -2.36 2.0 112ms
7 -36.74250114539 -4.02 -2.87 1.0 91.6ms
8 -36.74250325650 -5.68 -2.87 1.0 89.3ms
9 -36.74251099499 -5.11 -3.25 1.0 93.2ms
10 -36.74251436258 -5.47 -3.99 2.0 104ms
11 -36.74251467368 -6.51 -4.13 4.0 144ms
12 -36.74251473347 -7.22 -4.30 2.0 100ms
13 -36.74251477089 -7.43 -4.96 2.0 100ms
14 -36.74251477202 -8.95 -5.09 4.0 146ms
15 -36.74251477263 -9.21 -5.16 2.0 96.3ms
16 -36.74251477217 + -9.33 -5.31 2.0 109ms
17 -36.74251477285 -9.16 -5.50 1.0 94.5ms
18 -36.74251477295 -10.00 -5.75 1.0 95.0ms
19 -36.74251477303 -10.14 -6.28 5.0 123ms
20 -36.74251477303 -11.08 -6.43 3.0 133ms
21 -36.74251477304 -11.47 -7.00 1.0 94.7ms
22 -36.74251477304 + -11.88 -6.66 3.0 141ms
23 -36.74251477304 -11.94 -7.04 2.0 109ms
24 -36.74251477304 -12.62 -7.68 2.0 118ms
25 -36.74251477304 -14.15 -7.87 2.0 127ms
26 -36.74251477304 + -14.15 -8.37 1.0 91.6ms
27 -36.74251477304 + -Inf -8.54 3.0 136ms
28 -36.74251477304 + -Inf -8.79 1.0 94.5ms
29 -36.74251477304 -14.15 -9.25 3.0 132ms
30 -36.74251477304 + -Inf -9.72 5.0 126ms
31 -36.74251477304 + -14.15 -9.93 3.0 124ms
32 -36.74251477304 -14.15 -10.01 2.0 110ms
33 -36.74251477304 -14.15 -10.46 1.0 91.5ms
34 -36.74251477304 + -Inf -10.97 3.0 133ms
35 -36.74251477304 + -13.85 -11.14 6.0 131ms
36 -36.74251477304 + -Inf -11.44 1.0 95.8ms
37 -36.74251477304 + -Inf -11.58 2.0 127ms
38 -36.74251477304 -13.85 -11.98 1.0 91.8ms
39 -36.74251477304 + -13.85 -12.17 3.0 137ms
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.02448777549218
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.244209848206076
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.723606478093838
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).