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.73155056969 -0.88 12.0 325ms
2 -36.67338466883 + -1.24 -1.42 1.0 83.9ms
3 +16.26415212675 + 1.72 -0.20 7.0 208ms
4 -36.03120079721 1.72 -0.94 6.0 210ms
5 -35.86413343355 + -0.78 -1.00 3.0 143ms
6 -35.06615587721 + -0.10 -0.95 4.0 147ms
7 -36.70369754410 0.21 -1.62 3.0 138ms
8 -36.74017603184 -1.44 -2.01 2.0 103ms
9 -36.74071961406 -3.26 -2.18 2.0 107ms
10 -36.74105839545 -3.47 -2.06 2.0 108ms
11 -36.74183524317 -3.11 -2.43 1.0 92.3ms
12 -36.74196029721 -3.90 -2.46 2.0 102ms
13 -36.74233446949 -3.43 -2.77 1.0 94.0ms
14 -36.74235968602 -4.60 -2.70 2.0 123ms
15 -36.74127107102 + -2.96 -2.51 3.0 134ms
16 -36.74240727431 -2.94 -3.03 3.0 133ms
17 -36.74214912616 + -3.59 -2.72 2.0 118ms
18 -36.74235138290 -3.69 -2.99 2.0 115ms
19 -36.74247044387 -3.92 -3.50 2.0 111ms
20 -36.74247940518 -5.05 -3.83 2.0 145ms
21 -36.74247977542 -6.43 -4.03 2.0 101ms
22 -36.74248050440 -6.14 -4.17 2.0 123ms
23 -36.74247957746 + -6.03 -3.98 2.0 109ms
24 -36.74248061887 -5.98 -4.59 2.0 103ms
25 -36.74248065472 -7.45 -4.71 3.0 124ms
26 -36.74248066715 -7.91 -4.92 1.0 93.4ms
27 -36.74248064807 + -7.72 -4.79 2.0 126ms
28 -36.74248065572 -8.12 -4.73 2.0 120ms
29 -36.74248066694 -7.95 -5.15 2.0 107ms
30 -36.74248067209 -8.29 -5.42 1.0 99.0ms
31 -36.74248066922 + -8.54 -5.24 3.0 125ms
32 -36.74248067263 -8.47 -5.76 2.0 106ms
33 -36.74248067083 + -8.75 -5.41 3.0 135ms
34 -36.74248067259 -8.75 -5.77 3.0 133ms
35 -36.74248067264 -10.37 -6.12 1.0 93.5ms
36 -36.74248067265 -10.77 -6.21 2.0 105ms
37 -36.74248067267 -10.89 -6.25 2.0 109ms
38 -36.74248067268 -10.81 -6.70 1.0 93.1ms
39 -36.74248067267 + -11.11 -6.50 3.0 137ms
40 -36.74248067267 -11.93 -6.49 3.0 135ms
41 -36.74248067266 + -10.87 -6.38 3.0 135ms
42 -36.74248067268 -10.68 -7.31 2.0 109ms
43 -36.74248067268 + -12.57 -7.23 3.0 152ms
44 -36.74248067268 -12.38 -7.83 2.0 106ms
45 -36.74248067268 + -13.15 -7.60 3.0 142ms
46 -36.74248067268 -13.55 -7.73 3.0 134ms
47 -36.74248067268 -13.37 -8.27 2.0 124ms
48 -36.74248067268 -14.15 -8.22 2.0 119ms
49 -36.74248067268 + -14.15 -8.39 2.0 118ms
50 -36.74248067268 + -Inf -8.49 2.0 109ms
51 -36.74248067268 -14.15 -8.66 1.0 90.4ms
52 -36.74248067268 + -13.85 -8.93 2.0 126ms
53 -36.74248067268 + -14.15 -9.06 2.0 109ms
54 -36.74248067268 -14.15 -9.14 2.0 109ms
55 -36.74248067268 -13.85 -9.36 1.0 93.6ms
56 -36.74248067268 + -Inf -9.54 2.0 132ms
57 -36.74248067268 + -Inf -10.39 2.0 104ms
58 -36.74248067268 + -Inf -9.59 4.0 176ms
59 -36.74248067268 + -Inf -10.39 4.0 160ms
60 -36.74248067268 + -13.85 -10.64 2.0 110ms
61 -36.74248067268 + -Inf -9.94 4.0 170ms
62 -36.74248067268 + -Inf -10.56 3.0 141ms
63 -36.74248067268 -13.85 -10.77 2.0 110ms
64 -36.74248067268 + -Inf -11.00 1.0 111ms
65 -36.74248067268 + -14.15 -11.23 2.0 97.4ms
66 -36.74248067268 + -Inf -11.39 2.0 108ms
67 -36.74248067268 -14.15 -10.39 4.0 170ms
68 -36.74248067268 + -Inf -11.27 4.0 160ms
69 -36.74248067268 + -13.85 -11.94 2.0 107ms
70 -36.74248067268 -13.67 -11.67 3.0 258ms
71 -36.74248067268 + -14.15 -11.80 3.0 746ms
72 -36.74248067268 + -Inf -11.87 2.0 119ms
73 -36.74248067268 + -Inf -12.06 2.0 112ms
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.73364764364 -0.88 12.0 345ms
2 -36.74001577941 -2.20 -1.37 1.0 91.1ms
3 -36.74072148360 -3.15 -1.85 3.0 114ms
4 -36.74212138590 -2.85 -2.14 5.0 131ms
5 -36.74232711252 -3.69 -2.72 2.0 126ms
6 -36.74228107071 + -4.34 -2.62 3.0 149ms
7 -36.74244488198 -3.79 -2.91 1.0 89.5ms
8 -36.74245356864 -5.06 -3.24 2.0 105ms
9 -36.74247763595 -4.62 -3.59 2.0 124ms
10 -36.74248062749 -5.52 -4.03 2.0 96.5ms
11 -36.74248065075 -7.63 -4.47 4.0 124ms
12 -36.74248067081 -7.70 -4.83 3.0 125ms
13 -36.74248067200 -8.92 -5.32 2.0 124ms
14 -36.74248067200 -11.76 -5.54 3.0 111ms
15 -36.74248067258 -9.24 -5.89 2.0 125ms
16 -36.74248067267 -10.01 -6.32 2.0 101ms
17 -36.74248067268 -11.27 -6.58 2.0 128ms
18 -36.74248067268 -11.59 -6.95 2.0 126ms
19 -36.74248067268 -12.43 -7.32 1.0 96.4ms
20 -36.74248067268 -13.19 -7.58 4.0 135ms
21 -36.74248067268 -13.55 -7.90 3.0 109ms
22 -36.74248067268 + -Inf -8.19 3.0 134ms
23 -36.74248067268 -13.85 -8.32 1.0 92.7ms
24 -36.74248067268 + -Inf -8.73 2.0 104ms
25 -36.74248067268 + -Inf -9.16 2.0 129ms
26 -36.74248067268 -14.15 -9.45 4.0 111ms
27 -36.74248067268 + -14.15 -9.61 2.0 120ms
28 -36.74248067268 + -Inf -10.02 2.0 111ms
29 -36.74248067268 + -13.85 -10.14 3.0 133ms
30 -36.74248067268 -13.85 -10.42 1.0 92.4ms
31 -36.74248067268 + -13.85 -10.57 3.0 122ms
32 -36.74248067268 -13.85 -10.94 1.0 96.0ms
33 -36.74248067268 + -Inf -11.30 2.0 129ms
34 -36.74248067268 + -Inf -11.54 2.0 105ms
35 -36.74248067268 + -Inf -12.05 1.0 96.0ms
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.02448898037906
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.244211113843555
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.723580312349016
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).