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.73201553146 -0.88 10.0 356ms
2 -36.67411460823 + -1.24 -1.48 1.0 85.0ms
3 +22.82137557195 + 1.77 -0.18 7.0 213ms
4 -36.60535033295 1.77 -1.09 6.0 211ms
5 -35.71220297287 + -0.05 -1.04 4.0 159ms
6 -36.06856179978 -0.45 -1.14 5.0 164ms
7 -36.73031740694 -0.18 -1.74 3.0 133ms
8 -36.73437048521 -2.39 -1.95 2.0 106ms
9 -36.73694155877 -2.59 -2.06 1.0 89.3ms
10 -36.74119001445 -2.37 -2.13 2.0 115ms
11 -36.74212809451 -3.03 -2.55 1.0 93.4ms
12 -36.74245357689 -3.49 -2.58 2.0 125ms
13 -36.74242514484 + -4.55 -2.82 1.0 95.2ms
14 -36.74233662533 + -4.05 -2.80 1.0 95.1ms
15 -36.72324783206 + -1.72 -1.93 4.0 161ms
16 -36.74227440387 -1.72 -2.71 4.0 152ms
17 -36.74199741377 + -3.56 -2.64 2.0 116ms
18 -36.74137354223 + -3.20 -2.53 3.0 136ms
19 -36.74250724059 -2.95 -3.35 3.0 137ms
20 -36.74251278246 -5.26 -3.65 2.0 128ms
21 -36.74251450041 -5.76 -3.90 2.0 108ms
22 -36.74251413594 + -6.44 -4.05 2.0 111ms
23 -36.74251445690 -6.49 -4.19 1.0 95.5ms
24 -36.74251469713 -6.62 -4.49 1.0 122ms
25 -36.74251473430 -7.43 -4.51 2.0 126ms
26 -36.74251476953 -7.45 -4.94 1.0 95.2ms
27 -36.74251475761 + -7.92 -4.95 3.0 128ms
28 -36.74251477098 -7.87 -5.37 2.0 108ms
29 -36.74251477294 -8.71 -5.78 2.0 128ms
30 -36.74251475136 + -7.67 -4.88 4.0 160ms
31 -36.74251477267 -7.67 -5.74 4.0 164ms
32 -36.74251477297 -9.51 -5.88 2.0 110ms
33 -36.74251477300 -10.53 -6.19 2.0 124ms
34 -36.74251477304 -10.51 -6.65 1.0 94.8ms
35 -36.74251477304 + -12.42 -6.66 3.0 150ms
36 -36.74251477304 -11.91 -6.88 1.0 95.2ms
37 -36.74251477304 -12.45 -6.96 1.0 92.5ms
38 -36.74251477303 + -11.56 -6.61 2.0 121ms
39 -36.74251477304 -11.53 -7.23 2.0 115ms
40 -36.74251477304 + -12.18 -7.03 3.0 136ms
41 -36.74251477304 -12.15 -7.61 2.0 108ms
42 -36.74251477304 + -13.11 -7.57 3.0 154ms
43 -36.74251477304 -13.55 -7.67 2.0 121ms
44 -36.74251477304 -13.37 -8.09 2.0 111ms
45 -36.74251477304 + -Inf -8.20 2.0 125ms
46 -36.74251477304 -13.85 -8.20 2.0 111ms
47 -36.74251477304 + -13.85 -8.53 2.0 111ms
48 -36.74251477304 -13.85 -8.32 3.0 136ms
49 -36.74251477304 + -14.15 -8.53 3.0 124ms
50 -36.74251477304 + -14.15 -8.87 2.0 109ms
51 -36.74251477304 -13.85 -9.28 2.0 215ms
52 -36.74251477304 + -14.15 -9.42 2.0 766ms
53 -36.74251477304 + -13.85 -9.15 3.0 127ms
54 -36.74251477304 -13.67 -9.52 3.0 146ms
55 -36.74251477304 + -Inf -9.70 2.0 109ms
56 -36.74251477304 + -13.85 -9.86 2.0 127ms
57 -36.74251477304 + -Inf -10.28 1.0 96.8ms
58 -36.74251477304 + -Inf -10.14 3.0 162ms
59 -36.74251477304 + -Inf -10.48 3.0 150ms
60 -36.74251477304 + -Inf -10.75 2.0 132ms
61 -36.74251477304 + -14.15 -10.48 3.0 140ms
62 -36.74251477304 -13.55 -10.97 3.0 137ms
63 -36.74251477304 + -14.15 -11.06 2.0 126ms
64 -36.74251477304 + -13.85 -11.40 2.0 108ms
65 -36.74251477304 -13.85 -11.86 2.0 105ms
66 -36.74251477304 + -14.15 -11.75 2.0 126ms
67 -36.74251477304 + -14.15 -11.89 2.0 107ms
68 -36.74251477304 -13.85 -11.96 1.0 91.7ms
69 -36.74251477304 + -Inf -12.53 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.73447095621 -0.88 13.0 364ms
2 -36.74022574722 -2.24 -1.36 1.0 87.9ms
3 -36.74068632701 -3.34 -1.71 4.0 119ms
4 -36.74224054288 -2.81 -2.28 1.0 91.8ms
5 -36.74228812447 -4.32 -2.49 6.0 141ms
6 -36.74243698616 -3.83 -2.43 2.0 97.1ms
7 -36.74246133958 -4.61 -2.98 1.0 92.7ms
8 -36.74251380694 -4.28 -3.61 2.0 98.7ms
9 -36.74251401047 -6.69 -3.49 4.0 165ms
10 -36.74251454220 -6.27 -3.79 1.0 91.3ms
11 -36.74251465346 -6.95 -4.26 2.0 100ms
12 -36.74251469408 -7.39 -4.42 6.0 133ms
13 -36.74251476879 -7.13 -4.74 2.0 128ms
14 -36.74251477175 -8.53 -5.03 1.0 96.5ms
15 -36.74251477269 -9.02 -5.28 2.0 119ms
16 -36.74251477298 -9.54 -5.60 1.0 95.3ms
17 -36.74251477302 -10.37 -5.85 2.0 101ms
18 -36.74251477302 + -11.83 -6.01 3.0 124ms
19 -36.74251477303 -10.90 -6.29 1.0 117ms
20 -36.74251477304 -11.37 -6.88 3.0 123ms
21 -36.74251477304 + -13.00 -6.81 4.0 140ms
22 -36.74251477304 -12.64 -7.19 1.0 95.6ms
23 -36.74251477304 -13.07 -7.57 3.0 124ms
24 -36.74251477304 + -14.15 -7.88 2.0 102ms
25 -36.74251477304 -13.85 -8.47 3.0 134ms
26 -36.74251477304 + -13.85 -8.58 3.0 134ms
27 -36.74251477304 + -Inf -8.94 2.0 111ms
28 -36.74251477304 -14.15 -9.42 2.0 122ms
29 -36.74251477304 -14.15 -9.65 4.0 141ms
30 -36.74251477304 + -Inf -9.94 1.0 95.8ms
31 -36.74251477304 + -Inf -10.11 5.0 116ms
32 -36.74251477304 + -Inf -10.25 2.0 118ms
33 -36.74251477304 + -13.85 -10.54 1.0 95.4ms
34 -36.74251477304 + -Inf -11.06 2.0 111ms
35 -36.74251477304 -13.85 -11.01 3.0 135ms
36 -36.74251477304 + -14.15 -11.48 1.0 92.4ms
37 -36.74251477304 -14.15 -11.72 2.0 171ms
38 -36.74251477304 + -14.15 -12.05 2.0 129ms
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.02448777129269
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.24420984379485
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.723584173550597
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).