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.73236690101 -0.88 12.0 352ms
2 -36.36973448476 + -0.44 -1.24 1.0 88.7ms
3 +129.3512968378 + 2.22 0.06 17.0 365ms
4 -32.68375703113 2.21 -0.69 10.0 277ms
5 -33.83685453077 0.06 -0.80 4.0 156ms
6 -32.56203414125 + 0.11 -0.73 4.0 158ms
7 -36.70960264946 0.62 -1.53 3.0 142ms
8 -36.72089355278 -1.95 -1.74 2.0 116ms
9 -36.73563052463 -1.83 -1.85 2.0 112ms
10 -36.74159893973 -2.22 -2.10 1.0 91.3ms
11 -36.73882074020 + -2.56 -2.11 3.0 119ms
12 -36.74070326216 -2.73 -2.26 1.0 93.1ms
13 -36.74045234619 + -3.60 -2.30 1.0 92.6ms
14 -36.74160814624 -2.94 -2.47 1.0 98.3ms
15 -36.73649975112 + -2.29 -2.16 3.0 116ms
16 -36.69239700707 + -1.36 -1.70 4.0 164ms
17 -36.74228493955 -1.30 -2.72 3.0 143ms
18 -36.74208859945 + -3.71 -2.57 2.0 131ms
19 -36.74230186718 -3.67 -2.87 2.0 108ms
20 -36.74242398075 -3.91 -3.16 2.0 108ms
21 -36.74244226957 -4.74 -3.19 3.0 133ms
22 -36.74247653005 -4.47 -3.53 2.0 97.8ms
23 -36.74247663863 -6.96 -3.70 2.0 113ms
24 -36.74247584982 + -6.10 -3.71 3.0 129ms
25 -36.74248023038 -5.36 -4.12 2.0 107ms
26 -36.74248056823 -6.47 -4.47 3.0 134ms
27 -36.74248041985 + -6.83 -4.19 3.0 135ms
28 -36.74248064009 -6.66 -4.74 2.0 113ms
29 -36.74248044956 + -6.72 -4.38 3.0 133ms
30 -36.74248066825 -6.66 -5.09 3.0 138ms
31 -36.74248067094 -8.57 -5.12 2.0 108ms
32 -36.74248067167 -9.14 -5.37 2.0 107ms
33 -36.74248066417 + -8.12 -4.99 3.0 136ms
34 -36.74248067249 -8.08 -5.86 2.0 125ms
35 -36.74248067264 -9.83 -5.87 3.0 135ms
36 -36.74248067173 + -9.04 -5.57 3.0 132ms
37 -36.74248067259 -9.07 -5.99 3.0 138ms
38 -36.74248067266 -10.13 -6.33 1.0 93.0ms
39 -36.74248067267 -10.92 -6.47 2.0 107ms
40 -36.74248067268 -11.27 -6.66 3.0 126ms
41 -36.74248067242 + -9.58 -5.84 3.0 157ms
42 -36.74248067268 -9.58 -6.93 3.0 144ms
43 -36.74248067268 -12.03 -7.32 1.0 99.5ms
44 -36.74248067268 -12.79 -7.27 3.0 136ms
45 -36.74248067268 -13.85 -7.43 2.0 109ms
46 -36.74248067268 -13.45 -7.73 2.0 110ms
47 -36.74248067268 + -Inf -7.98 2.0 102ms
48 -36.74248067268 + -Inf -7.86 2.0 131ms
49 -36.74248067268 + -Inf -7.97 2.0 110ms
50 -36.74248067268 + -14.15 -8.20 2.0 108ms
51 -36.74248067268 -14.15 -8.29 3.0 134ms
52 -36.74248067268 + -Inf -8.78 2.0 102ms
53 -36.74248067268 -13.67 -8.03 4.0 172ms
54 -36.74248067268 + -Inf -9.08 4.0 159ms
55 -36.74248067268 + -13.85 -9.44 1.0 214ms
56 -36.74248067268 -13.85 -9.22 3.0 142ms
57 -36.74248067268 + -13.85 -9.74 2.0 1.21s
58 -36.74248067268 -13.85 -9.61 3.0 134ms
59 -36.74248067268 + -14.15 -9.97 2.0 109ms
60 -36.74248067268 + -14.15 -9.95 2.0 108ms
61 -36.74248067268 + -Inf -10.44 2.0 103ms
62 -36.74248067268 + -Inf -10.46 2.0 127ms
63 -36.74248067268 + -Inf -10.25 2.0 112ms
64 -36.74248067268 + -Inf -10.62 2.0 133ms
65 -36.74248067268 -13.85 -10.19 3.0 164ms
66 -36.74248067268 + -Inf -11.08 3.0 152ms
67 -36.74248067268 + -14.15 -10.81 3.0 161ms
68 -36.74248067268 + -13.85 -11.62 2.0 133ms
69 -36.74248067268 -13.67 -11.62 3.0 136ms
70 -36.74248067268 + -13.85 -11.51 3.0 127ms
71 -36.74248067268 -13.85 -12.12 2.0 108mswhile 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.73232735757 -0.88 11.0 340ms
2 -36.73969462456 -2.13 -1.36 1.0 98.7ms
3 -36.73873150652 + -3.02 -1.58 3.0 122ms
4 -36.74227607858 -2.45 -2.27 1.0 91.2ms
5 -36.74242452012 -3.83 -2.47 3.0 135ms
6 -36.74245501251 -4.52 -2.60 6.0 118ms
7 -36.74246170976 -5.17 -2.75 1.0 98.2ms
8 -36.74247712705 -4.81 -3.15 1.0 93.7ms
9 -36.74247915094 -5.69 -3.33 4.0 122ms
10 -36.74248049861 -5.87 -3.96 2.0 104ms
11 -36.74248063218 -6.87 -4.23 3.0 139ms
12 -36.74248066269 -7.52 -4.64 2.0 106ms
13 -36.74248065871 + -8.40 -4.86 2.0 129ms
14 -36.74248067257 -7.86 -5.51 2.0 117ms
15 -36.74248067247 + -9.99 -5.65 3.0 135ms
16 -36.74248067268 -9.68 -6.24 2.0 110ms
17 -36.74248067265 + -10.47 -6.17 3.0 141ms
18 -36.74248067268 -10.45 -6.79 2.0 116ms
19 -36.74248067268 -12.18 -6.97 3.0 135ms
20 -36.74248067268 -12.55 -7.26 1.0 101ms
21 -36.74248067268 -12.79 -7.58 3.0 118ms
22 -36.74248067268 + -Inf -7.73 4.0 131ms
23 -36.74248067268 -13.85 -8.24 1.0 95.8ms
24 -36.74248067268 + -Inf -8.36 3.0 141ms
25 -36.74248067268 + -13.85 -8.66 1.0 95.7ms
26 -36.74248067268 -13.85 -9.12 3.0 116ms
27 -36.74248067268 + -Inf -9.30 3.0 150ms
28 -36.74248067268 + -Inf -9.93 2.0 101ms
29 -36.74248067268 + -Inf -9.89 4.0 158ms
30 -36.74248067268 + -Inf -10.44 2.0 111ms
31 -36.74248067268 + -Inf -10.47 3.0 145ms
32 -36.74248067268 -14.15 -11.03 1.0 95.8ms
33 -36.74248067268 + -13.85 -11.10 2.0 134ms
34 -36.74248067268 -14.15 -11.41 1.0 95.8ms
35 -36.74248067268 + -13.85 -11.76 4.0 124ms
36 -36.74248067268 -13.85 -12.06 3.0 136msGiven 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.02448898024415The 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.244211113701844This 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.723586500758022Since 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).