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.73216364321 -0.88 11.0 1.44s
2 -36.55079386063 + -0.74 -1.35 1.0 249ms
3 +50.24682644396 + 1.94 -0.09 8.0 402ms
4 -36.38894954380 1.94 -1.06 8.0 242ms
5 -31.85559180081 + 0.66 -0.71 4.0 154ms
6 -36.55114977436 0.67 -1.36 4.0 1.44s
7 -36.68828110012 -0.86 -1.58 3.0 121ms
8 -36.73119654197 -1.37 -1.90 2.0 105ms
9 -36.73339770511 -2.66 -1.87 2.0 123ms
10 -36.74081604482 -2.13 -2.16 2.0 106ms
11 -36.74114128035 -3.49 -2.24 1.0 90.6ms
12 -36.74215075874 -3.00 -2.50 1.0 91.3ms
13 -36.74232071480 -3.77 -2.61 2.0 100ms
14 -36.74242459499 -3.98 -2.84 2.0 98.1ms
15 -36.73478678560 + -2.12 -2.12 4.0 144ms
16 -36.74147077781 -2.17 -2.52 3.0 141ms
17 -36.71381983107 + -1.56 -1.82 4.0 170ms
18 -36.74242454486 -1.54 -3.06 4.0 169ms
19 -36.74243385370 -5.03 -3.14 3.0 130ms
20 -36.74231281835 + -3.92 -2.80 2.0 118ms
21 -36.74247647283 -3.79 -3.44 3.0 124ms
22 -36.74247758077 -5.96 -3.68 2.0 109ms
23 -36.74247974940 -5.66 -3.96 1.0 92.8ms
24 -36.74248049330 -6.13 -4.31 2.0 135ms
25 -36.74248063972 -6.83 -4.60 2.0 101ms
26 -36.74248066662 -7.57 -4.75 2.0 125ms
27 -36.74248062323 + -7.36 -4.69 3.0 115ms
28 -36.74248067194 -7.31 -5.42 2.0 107ms
29 -36.74248059625 + -7.12 -4.57 4.0 159ms
30 -36.74248066916 -7.14 -5.27 4.0 169ms
31 -36.74248067096 -8.75 -5.39 3.0 133ms
32 -36.74248067254 -8.80 -5.83 2.0 107ms
33 -36.74248067123 + -8.88 -5.40 3.0 141ms
34 -36.74248067265 -8.85 -6.22 3.0 133ms
35 -36.74248067260 + -10.29 -5.96 2.0 115ms
36 -36.74248067266 -10.21 -6.39 2.0 116ms
37 -36.74248067267 -11.05 -6.40 3.0 136ms
38 -36.74248067268 -11.09 -6.83 2.0 101ms
39 -36.74248067268 + -11.42 -6.60 2.0 124ms
40 -36.74248067268 -11.30 -7.31 3.0 122ms
41 -36.74248067268 + -11.96 -7.01 3.0 141ms
42 -36.74248067268 + -11.82 -6.82 4.0 146ms
43 -36.74248067268 -11.66 -7.16 3.0 148ms
44 -36.74248067268 -12.30 -7.87 2.0 107ms
45 -36.74248067268 + -13.07 -7.44 3.0 139ms
46 -36.74248067268 -13.15 -7.84 3.0 132ms
47 -36.74248067268 -13.85 -8.11 2.0 113ms
48 -36.74248067268 -13.67 -8.48 3.0 110ms
49 -36.74248067268 + -13.85 -8.61 2.0 131ms
50 -36.74248067268 + -Inf -8.88 2.0 100ms
51 -36.74248067268 + -Inf -8.74 2.0 122ms
52 -36.74248067268 + -14.15 -9.36 2.0 100ms
53 -36.74248067268 + -Inf -9.30 3.0 139ms
54 -36.74248067268 -13.85 -9.27 2.0 116ms
55 -36.74248067268 + -14.15 -9.53 3.0 117ms
56 -36.74248067268 -14.15 -9.88 3.0 112ms
57 -36.74248067268 + -14.15 -10.09 2.0 131ms
58 -36.74248067268 -14.15 -10.13 2.0 106ms
59 -36.74248067268 + -Inf -10.11 2.0 133ms
60 -36.74248067268 + -Inf -10.55 2.0 100ms
61 -36.74248067268 + -Inf -10.81 3.0 137ms
62 -36.74248067268 + -14.15 -10.47 3.0 132ms
63 -36.74248067268 + -Inf -11.04 2.0 112ms
64 -36.74248067268 + -Inf -11.41 2.0 107ms
65 -36.74248067268 -13.85 -11.01 3.0 142ms
66 -36.74248067268 + -13.85 -10.96 3.0 133ms
67 -36.74248067268 + -14.15 -11.69 3.0 128ms
68 -36.74248067268 -14.15 -11.35 3.0 142ms
69 -36.74248067268 + -13.85 -11.81 3.0 132ms
70 -36.74248067268 -13.85 -11.36 3.0 131ms
71 -36.74248067268 -14.15 -12.08 3.0 137ms
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.73171668377 -0.88 12.0 921ms
2 -36.73938154999 -2.12 -1.37 1.0 685ms
3 -36.74036846405 -3.01 -1.75 4.0 126ms
4 -36.74209506635 -2.76 -2.12 4.0 128ms
5 -36.74221996344 -3.90 -2.46 3.0 112ms
6 -36.74243157016 -3.67 -2.50 2.0 103ms
7 -36.74243759852 -5.22 -3.05 1.0 99.2ms
8 -36.74247783241 -4.40 -3.39 3.0 130ms
9 -36.74248033378 -5.60 -3.66 2.0 102ms
10 -36.74248029589 + -7.42 -3.61 5.0 112ms
11 -36.74248060735 -6.51 -4.33 1.0 97.8ms
12 -36.74248066775 -7.22 -4.72 3.0 134ms
13 -36.74248067129 -8.45 -4.83 2.0 109ms
14 -36.74248067198 -9.16 -5.26 1.0 92.7ms
15 -36.74248067262 -9.19 -5.72 5.0 136ms
16 -36.74248067267 -10.30 -6.01 3.0 149ms
17 -36.74248067267 + -11.46 -6.23 2.0 108ms
18 -36.74248067268 -10.83 -6.71 1.0 98.5ms
19 -36.74248067268 -12.03 -7.06 6.0 146ms
20 -36.74248067268 -12.81 -7.37 2.0 116ms
21 -36.74248067268 -14.15 -7.41 2.0 126ms
22 -36.74248067268 + -Inf -7.82 1.0 98.4ms
23 -36.74248067268 -13.85 -8.12 2.0 118ms
24 -36.74248067268 + -Inf -8.61 2.0 107ms
25 -36.74248067268 + -Inf -8.97 3.0 132ms
26 -36.74248067268 + -Inf -9.27 2.0 109ms
27 -36.74248067268 -14.15 -9.50 5.0 121ms
28 -36.74248067268 + -14.15 -10.05 2.0 105ms
29 -36.74248067268 + -14.15 -10.17 3.0 138ms
30 -36.74248067268 + -Inf -10.37 1.0 92.7ms
31 -36.74248067268 -13.85 -10.51 1.0 98.5ms
32 -36.74248067268 -14.15 -11.00 2.0 98.3ms
33 -36.74248067268 + -13.85 -11.31 4.0 143ms
34 -36.74248067268 + -13.85 -11.52 2.0 98.4ms
35 -36.74248067268 -13.85 -11.78 3.0 241ms
36 -36.74248067268 -13.85 -12.03 1.0 92.8ms
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.024488980307616The 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.24421111376851This 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.723583824108366Since 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).