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.73472967803 -0.88 11.0 332ms
2 -36.72289233881 + -1.93 -1.60 1.0 84.6ms
3 -13.99306335963 + 1.36 -0.38 6.0 186ms
4 -36.07200451854 1.34 -1.07 5.0 173ms
5 -36.48028038045 -0.39 -1.31 3.0 136ms
6 -36.07940414480 + -0.40 -1.15 3.0 132ms
7 -36.72438084631 -0.19 -1.75 3.0 124ms
8 -36.74062816528 -1.79 -2.06 2.0 103ms
9 -36.74151937055 -3.05 -2.28 2.0 122ms
10 -36.74106919403 + -3.35 -2.12 2.0 113ms
11 -36.74245010149 -2.86 -2.71 2.0 111ms
12 -36.74244651263 + -5.45 -2.70 2.0 114ms
13 -36.74245676283 -4.99 -2.98 1.0 89.8ms
14 -36.74249236677 -4.45 -3.27 1.0 93.2ms
15 -36.74069213646 + -2.74 -2.43 4.0 151ms
16 -36.74250299285 -2.74 -3.32 4.0 146ms
17 -36.74239213008 + -3.96 -2.98 3.0 132ms
18 -36.74222178294 + -3.77 -2.83 3.0 129ms
19 -36.74251303232 -3.54 -3.72 3.0 213ms
20 -36.74251407269 -5.98 -3.96 2.0 793ms
21 -36.74251450477 -6.36 -4.16 2.0 105ms
22 -36.74251469092 -6.73 -4.45 1.0 89.4ms
23 -36.74251476046 -7.16 -4.67 2.0 121ms
24 -36.74251476301 -8.59 -4.67 2.0 104ms
25 -36.74251477122 -8.09 -5.02 1.0 89.5ms
26 -36.74251477295 -8.76 -5.68 2.0 113ms
27 -36.74251477274 + -9.68 -5.70 3.0 133ms
28 -36.74251477069 + -8.69 -5.36 3.0 127ms
29 -36.74251477209 -8.86 -5.57 3.0 132ms
30 -36.74251477030 + -8.75 -5.34 3.0 128ms
31 -36.74251477303 -8.56 -6.38 3.0 128ms
32 -36.74251477302 + -11.36 -6.39 2.0 125ms
33 -36.74251477304 -10.82 -6.87 2.0 97.6ms
34 -36.74251477304 -12.11 -7.29 3.0 116ms
35 -36.74251477304 + -12.04 -6.96 3.0 136ms
36 -36.74251477304 -12.02 -7.28 3.0 129ms
37 -36.74251477304 -13.30 -7.60 1.0 93.0ms
38 -36.74251477304 -13.67 -7.96 1.0 93.4ms
39 -36.74251477304 + -13.85 -7.80 3.0 131ms
40 -36.74251477304 -13.85 -8.06 1.0 92.9ms
41 -36.74251477304 + -12.40 -7.26 3.0 141ms
42 -36.74251477304 -12.42 -7.78 4.0 147ms
43 -36.74251477304 -13.67 -8.42 2.0 118ms
44 -36.74251477304 + -14.15 -8.51 2.0 116ms
45 -36.74251477304 -13.85 -8.74 1.0 89.5ms
46 -36.74251477304 + -Inf -8.94 2.0 125ms
47 -36.74251477304 + -13.85 -9.39 1.0 92.3ms
48 -36.74251477304 -13.85 -9.02 3.0 134ms
49 -36.74251477304 + -14.15 -9.62 2.0 110ms
50 -36.74251477304 + -Inf -9.90 2.0 116ms
51 -36.74251477304 + -Inf -10.01 2.0 124ms
52 -36.74251477304 -14.15 -10.16 2.0 107ms
53 -36.74251477304 + -14.15 -9.98 3.0 125ms
54 -36.74251477304 -14.15 -10.08 3.0 127ms
55 -36.74251477304 + -13.85 -10.20 2.0 116ms
56 -36.74251477304 -13.85 -10.88 2.0 102ms
57 -36.74251477304 + -13.85 -10.94 3.0 138ms
58 -36.74251477304 + -14.15 -10.67 3.0 140ms
59 -36.74251477304 -13.85 -11.34 2.0 104ms
60 -36.74251477304 + -Inf -11.04 3.0 139ms
61 -36.74251477304 + -Inf -11.56 2.0 116ms
62 -36.74251477304 + -14.15 -11.56 3.0 125ms
63 -36.74251477304 -13.85 -11.82 2.0 107ms
64 -36.74251477304 + -13.85 -11.98 2.0 121ms
65 -36.74251477304 -14.15 -12.10 2.0 107ms
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.73335435553 -0.88 11.0 327ms
2 -36.74015206395 -2.17 -1.37 1.0 87.8ms
3 -36.74169413286 -2.81 -2.19 2.0 100ms
4 -36.74213596928 -3.35 -2.08 6.0 144ms
5 -36.74241790304 -3.55 -2.65 1.0 89.9ms
6 -36.74241677105 + -5.95 -2.56 5.0 114ms
7 -36.74247931263 -4.20 -2.96 1.0 87.8ms
8 -36.74250687898 -4.56 -3.44 2.0 100ms
9 -36.74251394167 -5.15 -3.63 2.0 123ms
10 -36.74251471950 -6.11 -4.02 2.0 96.4ms
11 -36.74251467633 + -7.36 -3.94 2.0 120ms
12 -36.74251476524 -7.05 -4.66 1.0 92.7ms
13 -36.74251476925 -8.40 -5.00 2.0 125ms
14 -36.74251477235 -8.51 -5.12 2.0 98.2ms
15 -36.74251477295 -9.22 -5.78 2.0 104ms
16 -36.74251477303 -10.12 -6.00 3.0 139ms
17 -36.74251477304 -11.26 -6.34 4.0 106ms
18 -36.74251477304 -11.94 -6.57 2.0 99.2ms
19 -36.74251477304 -13.55 -6.81 2.0 107ms
20 -36.74251477304 -12.36 -7.18 2.0 106ms
21 -36.74251477304 + -13.55 -7.41 4.0 136ms
22 -36.74251477304 -13.37 -7.61 2.0 108ms
23 -36.74251477304 -13.85 -7.99 2.0 96.2ms
24 -36.74251477304 + -13.85 -8.28 5.0 114ms
25 -36.74251477304 -13.85 -8.54 2.0 126ms
26 -36.74251477304 + -13.85 -8.86 1.0 94.2ms
27 -36.74251477304 -14.15 -9.23 1.0 94.4ms
28 -36.74251477304 + -Inf -9.48 2.0 123ms
29 -36.74251477304 -14.15 -9.90 2.0 99.1ms
30 -36.74251477304 + -13.85 -10.31 3.0 120ms
31 -36.74251477304 + -Inf -10.51 3.0 131ms
32 -36.74251477304 -14.15 -10.77 1.0 94.3ms
33 -36.74251477304 + -Inf -10.99 2.0 109ms
34 -36.74251477304 + -14.15 -11.49 2.0 96.3ms
35 -36.74251477304 -13.85 -11.73 4.0 130ms
36 -36.74251477304 + -Inf -12.00 3.0 132ms
37 -36.74251477304 + -13.85 -12.19 2.0 105ms
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.02448820590106
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.24421030031624
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.7235878926196655
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).