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.73314578126 -0.88 11.0 350ms
2 -36.70263177648 + -1.52 -1.56 1.0 90.7ms
3 +3.439197840320 + 1.60 -0.27 7.0 202ms
4 -36.56809312311 1.60 -1.15 6.0 205ms
5 -36.70088969328 -0.88 -1.62 2.0 118ms
6 -35.94243030072 + -0.12 -1.12 3.0 125ms
7 -36.71243585342 -0.11 -1.72 3.0 142ms
8 -36.74132663903 -1.54 -2.27 2.0 109ms
9 -36.73742791590 + -2.41 -1.98 3.0 138ms
10 -36.74002458136 -2.59 -2.17 2.0 120ms
11 -36.74174554218 -2.76 -2.38 2.0 137ms
12 -36.74246477960 -3.14 -2.87 1.0 92.0ms
13 -36.74250605570 -4.38 -3.08 2.0 129ms
14 -36.74248029362 + -4.59 -2.85 3.0 126ms
15 -36.74235189801 + -3.89 -2.94 2.0 117ms
16 -36.74233832053 + -4.87 -2.82 3.0 135ms
17 -36.74222148505 + -3.93 -2.77 2.0 122ms
18 -36.74209222143 + -3.89 -2.75 3.0 131ms
19 -36.74250669275 -3.38 -3.32 2.0 120ms
20 -36.74251382468 -5.15 -3.99 2.0 114ms
21 -36.74251454770 -6.14 -4.21 3.0 142ms
22 -36.74251448955 + -7.24 -4.18 2.0 113ms
23 -36.74251429931 + -6.72 -4.03 2.0 125ms
24 -36.74251474788 -6.35 -4.54 2.0 111ms
25 -36.74251475203 -8.38 -4.53 2.0 117ms
26 -36.74251475756 -8.26 -4.77 2.0 103ms
27 -36.74251476711 -8.02 -4.94 1.0 91.4ms
28 -36.74251475941 + -8.11 -4.89 2.0 129ms
29 -36.74251476079 -8.86 -4.91 2.0 125ms
30 -36.74251471057 + -7.30 -4.67 3.0 130ms
31 -36.74251477289 -7.21 -5.88 3.0 136ms
32 -36.74251477205 + -9.08 -5.55 4.0 161ms
33 -36.74251477286 -9.09 -5.89 3.0 130ms
34 -36.74251477303 -9.78 -6.37 2.0 110ms
35 -36.74251477303 -11.80 -6.52 2.0 128ms
36 -36.74251477303 + -11.70 -6.36 2.0 107ms
37 -36.74251477304 -11.08 -6.96 2.0 105ms
38 -36.74251477304 -12.89 -7.32 3.0 109ms
39 -36.74251477304 -13.67 -7.37 3.0 144ms
40 -36.74251477304 + -13.85 -7.36 1.0 95.7ms
41 -36.74251477304 -13.19 -7.50 2.0 116ms
42 -36.74251477304 + -13.15 -7.42 2.0 117ms
43 -36.74251477304 + -13.19 -7.42 3.0 129ms
44 -36.74251477304 -13.07 -7.51 3.0 123ms
45 -36.74251477304 -13.03 -8.11 2.0 105ms
46 -36.74251477304 + -Inf -8.37 3.0 155ms
47 -36.74251477304 + -13.07 -7.57 4.0 147ms
48 -36.74251477304 -13.11 -8.54 4.0 153ms
49 -36.74251477304 -14.15 -8.93 2.0 104ms
50 -36.74251477304 + -Inf -8.74 3.0 140ms
51 -36.74251477304 + -14.15 -9.35 2.0 111ms
52 -36.74251477304 -14.15 -9.25 3.0 130ms
53 -36.74251477304 + -14.15 -9.32 2.0 116ms
54 -36.74251477304 + -Inf -9.50 2.0 106ms
55 -36.74251477304 -14.15 -9.37 2.0 110ms
56 -36.74251477304 + -Inf -9.86 2.0 101ms
57 -36.74251477304 + -Inf -10.06 2.0 129ms
58 -36.74251477304 + -Inf -10.41 2.0 110ms
59 -36.74251477304 + -Inf -10.59 2.0 124ms
60 -36.74251477304 + -Inf -10.44 3.0 134ms
61 -36.74251477304 + -Inf -11.11 2.0 178ms
62 -36.74251477304 + -14.15 -11.15 2.0 125ms
63 -36.74251477304 -14.15 -11.30 2.0 665ms
64 -36.74251477304 + -13.85 -11.78 2.0 107ms
65 -36.74251477304 + -Inf -11.83 3.0 136ms
66 -36.74251477304 -13.85 -11.92 2.0 107ms
67 -36.74251477304 + -13.85 -12.12 3.0 145ms
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.73223540392 -0.88 12.0 339ms
2 -36.73972355993 -2.13 -1.36 1.0 92.0ms
3 -36.74032644836 -3.22 -1.67 4.0 156ms
4 -36.74228920105 -2.71 -2.21 2.0 94.5ms
5 -36.74240459562 -3.94 -2.45 6.0 121ms
6 -36.74244086323 -4.44 -2.44 2.0 102ms
7 -36.74250438069 -4.20 -3.19 1.0 91.8ms
8 -36.74251200847 -5.12 -3.17 3.0 136ms
9 -36.74251367524 -5.78 -3.53 1.0 91.2ms
10 -36.74251446645 -6.10 -4.04 2.0 97.5ms
11 -36.74251474108 -6.56 -4.33 3.0 131ms
12 -36.74251476830 -7.57 -4.50 2.0 103ms
13 -36.74251476698 + -8.88 -4.83 2.0 100ms
14 -36.74251477278 -8.24 -5.25 5.0 125ms
15 -36.74251477292 -9.86 -5.43 2.0 126ms
16 -36.74251477297 -10.36 -5.50 1.0 117ms
17 -36.74251477294 + -10.56 -5.58 2.0 118ms
18 -36.74251477303 -10.02 -6.19 1.0 97.9ms
19 -36.74251477304 -11.38 -6.56 3.0 139ms
20 -36.74251477304 + -12.22 -6.77 2.0 163ms
21 -36.74251477304 -12.05 -7.39 2.0 105ms
22 -36.74251477304 + -12.94 -7.31 3.0 148ms
23 -36.74251477304 -13.03 -7.48 2.0 113ms
24 -36.74251477304 + -Inf -7.63 1.0 97.7ms
25 -36.74251477304 -13.85 -8.40 2.0 129ms
26 -36.74251477304 + -Inf -8.65 3.0 147ms
27 -36.74251477304 -14.15 -8.91 2.0 103ms
28 -36.74251477304 + -14.15 -8.98 3.0 123ms
29 -36.74251477304 + -13.85 -9.31 1.0 93.4ms
30 -36.74251477304 -14.15 -9.64 4.0 114ms
31 -36.74251477304 + -Inf -10.05 2.0 122ms
32 -36.74251477304 -14.15 -10.35 3.0 139ms
33 -36.74251477304 + -13.85 -10.53 2.0 104ms
34 -36.74251477304 -13.85 -11.16 2.0 127ms
35 -36.74251477304 + -13.85 -11.16 4.0 141ms
36 -36.74251477304 -13.85 -11.39 1.0 97.6ms
37 -36.74251477304 + -13.85 -11.84 2.0 112ms
38 -36.74251477304 -13.85 -11.77 3.0 137ms
39 -36.74251477304 + -Inf -12.15 1.0 94.3ms
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.024487775913514
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.24420984864865
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.7236239152456525
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).