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.73317929003 -0.88 12.0 349ms
2 -36.60741085754 + -0.90 -1.41 1.0 83.9ms
3 +41.69054442945 + 1.89 -0.11 8.0 233ms
4 -36.20199275035 1.89 -1.06 8.0 230ms
5 -33.15660984138 + 0.48 -0.79 4.0 141ms
6 -36.64580930379 0.54 -1.48 6.0 181ms
7 -36.72943166048 -1.08 -1.74 2.0 101ms
8 -36.73436041705 -2.31 -2.01 2.0 103ms
9 -36.73319207138 + -2.93 -1.91 2.0 122ms
10 -36.73939336882 -2.21 -2.13 2.0 105ms
11 -36.74122960579 -2.74 -2.23 1.0 91.3ms
12 -36.74224747154 -2.99 -2.61 1.0 88.7ms
13 -36.74242166124 -3.76 -2.72 2.0 106ms
14 -36.74245261532 -4.51 -2.86 2.0 101ms
15 -36.73614363878 + -2.20 -2.16 3.0 140ms
16 -36.74240776785 -2.20 -2.97 4.0 141ms
17 -36.73381973768 + -2.07 -2.09 4.0 165ms
18 -36.74248168544 -2.06 -3.23 3.0 139ms
19 -36.74248975109 -5.09 -3.25 2.0 124ms
20 -36.74251013467 -4.69 -3.54 2.0 103ms
21 -36.74250642804 + -5.43 -3.52 3.0 123ms
22 -36.74251467633 -5.08 -4.14 2.0 185ms
23 -36.74251458914 + -7.06 -4.36 3.0 780ms
24 -36.74251471920 -6.89 -4.54 1.0 91.9ms
25 -36.74251474558 -7.58 -4.69 2.0 107ms
26 -36.74251476362 -7.74 -4.87 2.0 106ms
27 -36.74251448719 + -6.56 -4.30 3.0 137ms
28 -36.74251476060 -6.56 -4.96 3.0 129ms
29 -36.74251471389 + -7.33 -4.66 3.0 133ms
30 -36.74251476416 -7.30 -4.98 3.0 129ms
31 -36.74251477256 -8.08 -5.57 2.0 99.4ms
32 -36.74251477230 + -9.58 -5.59 3.0 141ms
33 -36.74251477236 -10.26 -5.53 2.0 114ms
34 -36.74251477301 -9.18 -6.20 2.0 99.1ms
35 -36.74251477283 + -9.74 -5.85 3.0 143ms
36 -36.74251477302 -9.72 -6.42 3.0 113ms
37 -36.74251477303 -11.10 -6.69 3.0 137ms
38 -36.74251477303 -11.66 -6.70 3.0 124ms
39 -36.74251477304 -12.18 -6.75 2.0 110ms
40 -36.74251477304 -11.73 -7.22 2.0 104ms
41 -36.74251477303 + -11.54 -6.81 3.0 143ms
42 -36.74251477304 -11.55 -7.37 3.0 125ms
43 -36.74251477304 -12.89 -7.69 2.0 104ms
44 -36.74251477304 + -Inf -7.79 2.0 119ms
45 -36.74251477304 + -13.15 -7.50 3.0 124ms
46 -36.74251477304 -13.07 -7.95 3.0 134ms
47 -36.74251477304 + -Inf -8.08 2.0 109ms
48 -36.74251477304 -14.15 -8.42 1.0 94.4ms
49 -36.74251477304 + -14.15 -7.98 3.0 145ms
50 -36.74251477304 -13.85 -8.40 3.0 128ms
51 -36.74251477304 + -13.85 -9.14 2.0 109ms
52 -36.74251477304 + -Inf -9.27 3.0 142ms
53 -36.74251477304 -14.15 -9.38 1.0 94.6ms
54 -36.74251477304 + -14.15 -9.65 1.0 94.5ms
55 -36.74251477304 -13.85 -9.30 3.0 131ms
56 -36.74251477304 + -14.15 -9.31 3.0 135ms
57 -36.74251477304 + -Inf -9.59 3.0 126ms
58 -36.74251477304 + -14.15 -10.02 2.0 120ms
59 -36.74251477304 -13.85 -10.43 2.0 128ms
60 -36.74251477304 + -Inf -10.05 3.0 132ms
61 -36.74251477304 + -13.85 -10.52 3.0 129ms
62 -36.74251477304 + -Inf -10.51 3.0 135ms
63 -36.74251477304 -13.85 -11.13 1.0 95.1ms
64 -36.74251477304 + -14.15 -11.17 3.0 140ms
65 -36.74251477304 -14.15 -11.25 2.0 111ms
66 -36.74251477304 + -14.15 -11.72 1.0 94.7ms
67 -36.74251477304 + -14.15 -11.94 3.0 138ms
68 -36.74251477304 -13.85 -11.77 2.0 128ms
69 -36.74251477304 + -13.85 -11.70 3.0 130ms
70 -36.74251477304 -13.85 -12.00 2.0 110ms
71 -36.74251477304 + -Inf -12.06 2.0 120ms
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.73083312271 -0.88 14.0 343ms
2 -36.73873861004 -2.10 -1.36 1.0 91.2ms
3 -36.73965477396 -3.04 -1.76 4.0 136ms
4 -36.74202167443 -2.63 -2.01 1.0 88.4ms
5 -36.74226947427 -3.61 -2.38 5.0 113ms
6 -36.74243979840 -3.77 -2.44 3.0 116ms
7 -36.74250511986 -4.18 -2.99 1.0 93.0ms
8 -36.74251316516 -5.09 -3.39 4.0 129ms
9 -36.74251268520 + -6.32 -3.31 2.0 127ms
10 -36.74251438079 -5.77 -3.69 1.0 94.7ms
11 -36.74251457669 -6.71 -4.19 1.0 95.1ms
12 -36.74251474138 -6.78 -4.63 3.0 115ms
13 -36.74251476988 -7.55 -4.78 3.0 134ms
14 -36.74251477205 -8.66 -5.26 1.0 96.4ms
15 -36.74251477269 -9.20 -5.26 3.0 135ms
16 -36.74251477286 -9.75 -5.48 1.0 96.7ms
17 -36.74251477303 -9.79 -6.08 2.0 101ms
18 -36.74251477302 + -11.22 -6.05 3.0 140ms
19 -36.74251477303 -10.94 -6.30 2.0 104ms
20 -36.74251477304 -11.32 -6.60 2.0 103ms
21 -36.74251477304 + -12.10 -6.58 3.0 131ms
22 -36.74251477304 -11.81 -7.01 1.0 92.6ms
23 -36.74251477304 -12.92 -7.49 3.0 121ms
24 -36.74251477304 + -Inf -7.89 3.0 137ms
25 -36.74251477304 -13.67 -7.96 2.0 101ms
26 -36.74251477304 + -14.15 -8.58 2.0 101ms
27 -36.74251477304 + -14.15 -8.56 4.0 144ms
28 -36.74251477304 -13.85 -8.77 1.0 94.1ms
29 -36.74251477304 + -Inf -9.24 3.0 108ms
30 -36.74251477304 + -13.85 -9.23 3.0 136ms
31 -36.74251477304 -13.85 -9.46 1.0 95.5ms
32 -36.74251477304 + -Inf -10.01 2.0 107ms
33 -36.74251477304 + -Inf -9.85 3.0 135ms
34 -36.74251477304 + -Inf -9.92 1.0 95.4ms
35 -36.74251477304 + -Inf -10.65 2.0 110ms
36 -36.74251477304 + -14.15 -10.94 3.0 134ms
37 -36.74251477304 -14.15 -11.06 2.0 100ms
38 -36.74251477304 + -Inf -11.48 1.0 95.6ms
39 -36.74251477304 + -Inf -11.51 3.0 133ms
40 -36.74251477304 -14.15 -12.00 1.0 95.4ms
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.0244888105471
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.244210935448635
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.7235883245209696
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).