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.73399000446 -0.88 11.0 359ms
2 -36.68500536825 + -1.31 -1.51 1.0 97.2ms
3 +16.16038012389 + 1.72 -0.20 7.0 213ms
4 -36.57117113419 1.72 -1.10 6.0 225ms
5 -36.70905814846 -0.86 -1.48 3.0 135ms
6 -36.73706930606 -1.55 -1.84 2.0 119ms
7 -36.51935833986 + -0.66 -1.39 3.0 134ms
8 -36.71739091010 -0.70 -1.82 3.0 141ms
9 -36.74053702972 -1.64 -2.21 2.0 105ms
10 -36.73924081992 + -2.89 -2.12 2.0 130ms
11 -36.74067409098 -2.84 -2.23 2.0 127ms
12 -36.74209389628 -2.85 -2.55 1.0 95.9ms
13 -36.74243349331 -3.47 -2.92 1.0 102ms
14 -36.74243954873 -5.22 -2.91 3.0 141ms
15 -36.74038400136 + -2.69 -2.38 3.0 144ms
16 -36.74231908637 -2.71 -2.87 4.0 154ms
17 -36.74247567815 -3.81 -3.51 2.0 110ms
18 -36.74246555397 + -4.99 -3.35 3.0 148ms
19 -36.74247817197 -4.90 -3.67 2.0 117ms
20 -36.74248007757 -5.72 -3.99 2.0 130ms
21 -36.74247790370 + -5.66 -3.84 3.0 142ms
22 -36.74248026866 -5.63 -4.22 2.0 122ms
23 -36.74247980867 + -6.34 -4.04 3.0 138ms
24 -36.74248065072 -6.07 -4.61 2.0 117ms
25 -36.74248065124 -9.29 -4.77 3.0 138ms
26 -36.74248067164 -7.69 -5.29 2.0 110ms
27 -36.74248058345 + -7.05 -4.55 4.0 165ms
28 -36.74248067129 -7.06 -5.30 3.0 153ms
29 -36.74248066662 + -8.33 -5.12 2.0 122ms
30 -36.74248067172 -8.29 -5.30 2.0 115ms
31 -36.74248067249 -9.11 -5.82 1.0 96.4ms
32 -36.74248067268 -9.72 -6.16 3.0 141ms
33 -36.74248067253 + -9.83 -5.91 3.0 142ms
34 -36.74248067183 + -9.16 -5.60 4.0 160ms
35 -36.74248067264 -9.09 -6.14 3.0 148ms
36 -36.74248067268 -10.43 -6.50 2.0 110ms
37 -36.74248067268 -11.63 -6.58 3.0 119ms
38 -36.74248067268 -11.46 -6.86 1.0 95.6ms
39 -36.74248067268 -12.66 -7.09 2.0 135ms
40 -36.74248067268 + -12.17 -7.04 2.0 112ms
41 -36.74248067268 + -12.94 -7.03 3.0 153ms
42 -36.74248067268 -12.01 -7.68 2.0 112ms
43 -36.74248067268 -13.67 -7.72 3.0 148ms
44 -36.74248067268 -13.67 -8.14 1.0 95.8ms
45 -36.74248067268 + -13.85 -8.52 2.0 112ms
46 -36.74248067268 + -14.15 -8.36 3.0 142ms
47 -36.74248067268 + -Inf -7.94 4.0 159ms
48 -36.74248067268 -14.15 -8.21 4.0 155ms
49 -36.74248067268 -13.67 -8.69 3.0 144ms
50 -36.74248067268 + -13.85 -8.76 2.0 104ms
51 -36.74248067268 + -Inf -9.18 2.0 113ms
52 -36.74248067268 -14.15 -9.49 2.0 136ms
53 -36.74248067268 + -14.15 -9.52 2.0 118ms
54 -36.74248067268 + -Inf -10.03 1.0 101ms
55 -36.74248067268 + -Inf -9.68 3.0 141ms
56 -36.74248067268 + -Inf -10.10 2.0 137ms
57 -36.74248067268 + -Inf -10.23 2.0 111ms
58 -36.74248067268 + -14.15 -10.42 1.0 95.6ms
59 -36.74248067268 -13.85 -10.55 3.0 133ms
60 -36.74248067268 + -14.15 -10.34 3.0 139ms
61 -36.74248067268 + -14.15 -10.87 3.0 134ms
62 -36.74248067268 -13.85 -11.13 2.0 130ms
63 -36.74248067268 + -14.15 -11.52 2.0 117ms
64 -36.74248067268 -13.85 -11.39 2.0 130ms
65 -36.74248067268 + -13.85 -11.95 1.0 100ms
66 -36.74248067268 -13.85 -11.81 3.0 139ms
67 -36.74248067268 + -13.55 -11.68 2.0 129ms
68 -36.74248067268 -13.85 -12.17 2.0 119mswhile 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.73290859376 -0.88 11.0 348ms
2 -36.73938279455 -2.19 -1.36 1.0 90.1ms
3 -36.73986521045 -3.32 -1.70 3.0 133ms
4 -36.74208106770 -2.65 -2.09 1.0 92.1ms
5 -36.74230891475 -3.64 -2.58 5.0 123ms
6 -36.74243722271 -3.89 -2.50 3.0 137ms
7 -36.74246077486 -4.63 -2.90 1.0 93.9ms
8 -36.74247568151 -4.83 -3.14 1.0 99.4ms
9 -36.74247975694 -5.39 -3.58 5.0 116ms
10 -36.74248053748 -6.11 -3.83 2.0 132ms
11 -36.74248064673 -6.96 -4.09 1.0 95.3ms
12 -36.74248055700 + -7.05 -4.35 2.0 116ms
13 -36.74248066463 -6.97 -4.40 2.0 131ms
14 -36.74248066381 + -9.09 -4.66 1.0 100ms
15 -36.74248067159 -8.11 -4.98 2.0 122ms
16 -36.74248067251 -9.04 -5.27 1.0 100ms
17 -36.74248067265 -9.84 -5.74 3.0 117ms
18 -36.74248067264 + -10.86 -5.81 3.0 134ms
19 -36.74248067268 -10.41 -6.28 2.0 101ms
20 -36.74248067268 -11.24 -6.62 3.0 141ms
21 -36.74248067268 + -12.94 -6.50 3.0 135ms
22 -36.74248067268 -12.35 -7.15 2.0 116ms
23 -36.74248067268 + -13.67 -7.14 3.0 148ms
24 -36.74248067268 -13.00 -7.43 1.0 100ms
25 -36.74248067268 + -14.15 -7.78 1.0 96.0ms
26 -36.74248067268 -13.85 -8.01 3.0 149ms
27 -36.74248067268 + -14.15 -8.37 1.0 95.2ms
28 -36.74248067268 + -13.85 -8.64 2.0 127ms
29 -36.74248067268 -13.67 -8.81 2.0 130ms
30 -36.74248067268 + -14.15 -9.12 1.0 97.2ms
31 -36.74248067268 + -Inf -9.92 2.0 126ms
32 -36.74248067268 + -Inf -10.01 4.0 152ms
33 -36.74248067268 + -Inf -10.24 2.0 107ms
34 -36.74248067268 + -Inf -10.49 2.0 102ms
35 -36.74248067268 + -Inf -10.91 3.0 139ms
36 -36.74248067268 + -Inf -11.14 3.0 107ms
37 -36.74248067268 + -Inf -11.49 1.0 100ms
38 -36.74248067268 -14.15 -11.67 3.0 139ms
39 -36.74248067268 + -14.15 -11.93 1.0 100ms
40 -36.74248067268 + -Inf -12.20 3.0 131msGiven 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.024488980366854The 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.24421111383073This 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.723585185687596Since 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).