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.73130141087 -0.88 11.0 1.33s
2 -36.41936442271 + -0.51 -1.27 1.0 289ms
3 +124.2913605595 + 2.21 0.05 16.0 354ms
4 -33.09505071428 2.20 -0.69 10.0 379ms
5 -36.26104333231 0.50 -1.12 4.0 1.09s
6 -31.18696284816 + 0.71 -0.69 4.0 155ms
7 -36.67407550002 0.74 -1.43 4.0 158ms
8 -36.73094010264 -1.25 -1.88 2.0 106ms
9 -36.73859116652 -2.12 -1.89 2.0 105ms
10 -36.74205326854 -2.46 -2.16 1.0 91.5ms
11 -36.73975535123 + -2.64 -2.19 3.0 116ms
12 -36.74007365739 -3.50 -2.27 1.0 94.3ms
13 -36.74135077904 -2.89 -2.48 1.0 93.9ms
14 -36.74199951533 -3.19 -2.65 1.0 97.7ms
15 -36.74207782993 -4.11 -2.66 2.0 133ms
16 -36.74061493578 + -2.83 -2.43 2.0 134ms
17 -36.74175348249 -2.94 -2.50 3.0 133ms
18 -36.74068049402 + -2.97 -2.42 2.0 113ms
19 -36.74204248170 -2.87 -2.64 3.0 126ms
20 -36.74245563699 -3.38 -3.24 2.0 112ms
21 -36.74245905557 -5.47 -3.08 2.0 131ms
22 -36.74247141470 -4.91 -3.45 2.0 113ms
23 -36.74247979888 -5.08 -3.78 2.0 103ms
24 -36.74248017863 -6.42 -3.99 2.0 127ms
25 -36.74248060603 -6.37 -4.21 1.0 98.7ms
26 -36.74248049263 + -6.95 -4.34 2.0 109ms
27 -36.74248046980 + -7.64 -4.28 3.0 128ms
28 -36.74247983757 + -6.20 -4.08 2.0 123ms
29 -36.74248052438 -6.16 -4.42 3.0 134ms
30 -36.74248066635 -6.85 -4.93 2.0 112ms
31 -36.74248066547 + -9.05 -5.09 3.0 136ms
32 -36.74248063790 + -7.56 -4.77 3.0 138ms
33 -36.74248061920 + -7.73 -4.61 3.0 138ms
34 -36.74248066760 -7.32 -5.03 3.0 134ms
35 -36.74248067050 -8.54 -5.25 2.0 112ms
36 -36.74248067251 -8.70 -5.85 1.0 93.8ms
37 -36.74248067258 -10.13 -5.91 3.0 147ms
38 -36.74248067261 -10.48 -6.04 2.0 113ms
39 -36.74248067245 + -9.80 -5.80 3.0 134ms
40 -36.74248067264 -9.74 -6.17 2.0 114ms
41 -36.74248067266 -10.66 -6.33 3.0 128ms
42 -36.74248067264 + -10.71 -6.18 3.0 138ms
43 -36.74248067268 -10.37 -6.78 2.0 112ms
44 -36.74248067267 + -11.20 -6.54 3.0 136ms
45 -36.74248067268 -11.36 -6.78 3.0 138ms
46 -36.74248067268 -11.55 -7.42 2.0 101ms
47 -36.74248067268 + -13.45 -7.42 3.0 147ms
48 -36.74248067268 -12.81 -7.76 2.0 113ms
49 -36.74248067268 + -13.07 -7.54 3.0 128ms
50 -36.74248067268 -13.03 -7.97 3.0 122ms
51 -36.74248067268 + -Inf -8.23 3.0 115ms
52 -36.74248067268 + -13.85 -8.07 3.0 139ms
53 -36.74248067268 -13.85 -8.41 2.0 106ms
54 -36.74248067268 + -14.15 -8.58 2.0 108ms
55 -36.74248067268 -13.85 -8.62 2.0 106ms
56 -36.74248067268 + -14.15 -8.62 1.0 93.4ms
57 -36.74248067268 + -Inf -9.03 2.0 107ms
58 -36.74248067268 -14.15 -8.85 3.0 136ms
59 -36.74248067268 + -14.15 -9.62 2.0 113ms
60 -36.74248067268 + -Inf -9.28 2.0 130ms
61 -36.74248067268 + -Inf -9.08 3.0 143ms
62 -36.74248067268 -14.15 -9.05 3.0 148ms
63 -36.74248067268 + -14.15 -9.43 3.0 133ms
64 -36.74248067268 + -Inf -9.39 3.0 138ms
65 -36.74248067268 + -Inf -10.21 2.0 112ms
66 -36.74248067268 + -Inf -10.33 2.0 127ms
67 -36.74248067268 + -14.15 -10.31 3.0 123ms
68 -36.74248067268 -14.15 -10.74 2.0 103ms
69 -36.74248067268 + -Inf -10.88 1.0 97.6ms
70 -36.74248067268 + -Inf -11.10 3.0 132ms
71 -36.74248067268 + -Inf -11.46 2.0 127ms
72 -36.74248067268 + -Inf -11.40 2.0 115ms
73 -36.74248067268 -13.85 -11.51 2.0 105ms
74 -36.74248067268 + -13.85 -11.68 1.0 97.2ms
75 -36.74248067268 + -Inf -11.86 3.0 115ms
76 -36.74248067268 + -Inf -11.41 3.0 147ms
77 -36.74248067268 + -Inf -11.80 3.0 137ms
78 -36.74248067268 + -Inf -11.91 2.0 109ms
79 -36.74248067268 + -Inf -12.43 1.0 195ms
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.73330200708 -0.88 11.0 1.81s
2 -36.73950564521 -2.21 -1.36 1.0 626ms
3 -36.73672668468 + -2.56 -1.53 3.0 129ms
4 -36.74214900873 -2.27 -2.25 1.0 90.9ms
5 -36.74196171128 + -3.73 -2.47 4.0 123ms
6 -36.74238719616 -3.37 -2.43 4.0 107ms
7 -36.74241255856 -4.60 -3.02 1.0 91.6ms
8 -36.74247922674 -4.18 -3.34 5.0 124ms
9 -36.74247895313 + -6.56 -3.67 2.0 116ms
10 -36.74248009114 -5.94 -3.71 2.0 131ms
11 -36.74248049148 -6.40 -4.12 1.0 104ms
12 -36.74248063512 -6.84 -4.57 3.0 142ms
13 -36.74248066225 -7.57 -4.80 2.0 156ms
14 -36.74248067062 -8.08 -5.07 1.0 118ms
15 -36.74248067250 -8.73 -5.62 2.0 159ms
16 -36.74248067267 -9.76 -5.97 3.0 134ms
17 -36.74248067268 -11.02 -6.54 3.0 109ms
18 -36.74248067268 -11.79 -6.74 4.0 138ms
19 -36.74248067268 -12.50 -7.24 4.0 111ms
20 -36.74248067268 -13.67 -7.47 3.0 140ms
21 -36.74248067268 -13.85 -7.91 1.0 94.9ms
22 -36.74248067268 + -14.15 -8.27 4.0 143ms
23 -36.74248067268 -14.15 -8.52 5.0 117ms
24 -36.74248067268 + -14.15 -8.84 2.0 105ms
25 -36.74248067268 -13.67 -8.98 2.0 129ms
26 -36.74248067268 + -13.67 -9.37 1.0 95.1ms
27 -36.74248067268 -13.67 -9.47 4.0 122ms
28 -36.74248067268 + -Inf -9.72 2.0 107ms
29 -36.74248067268 + -13.67 -10.30 1.0 94.9ms
30 -36.74248067268 + -Inf -10.47 3.0 143ms
31 -36.74248067268 -13.85 -10.78 2.0 132ms
32 -36.74248067268 + -14.15 -11.00 2.0 110ms
33 -36.74248067268 + -Inf -11.43 1.0 99.3ms
34 -36.74248067268 + -Inf -11.83 2.0 128ms
35 -36.74248067268 + -Inf -12.14 2.0 109ms
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.02448898041253The 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.244211113878706This 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.723562796213271Since 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).