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.73341725471 -0.88 12.0 369ms
2 -36.50812345319 + -0.65 -1.30 1.0 97.7ms
3 +102.3027491901 + 2.14 0.02 21.0 349ms
4 -35.93993413728 2.14 -0.96 9.0 300ms
5 -26.36731724714 + 0.98 -0.55 5.0 176ms
6 -36.51598429397 1.01 -1.27 4.0 170ms
7 -36.66270656606 -0.83 -1.51 2.0 120ms
8 -36.67216770206 -2.02 -1.59 2.0 122ms
9 -36.73291132745 -1.22 -1.85 2.0 103ms
10 -36.74130275946 -2.08 -2.09 1.0 223ms
11 -36.73638328958 + -2.31 -2.02 2.0 121ms
12 -36.73835971603 -2.70 -2.17 1.0 1.40s
13 -36.74152373725 -2.50 -2.46 1.0 98.0ms
14 -36.74229034472 -3.12 -2.64 3.0 123ms
15 -36.73670975824 + -2.25 -2.16 3.0 140ms
16 -36.73843884860 -2.76 -2.20 3.0 148ms
17 -36.73871180014 -3.56 -2.27 3.0 138ms
18 -36.74246907801 -2.43 -3.30 3.0 146ms
19 -36.74245205428 + -4.77 -3.15 3.0 158ms
20 -36.74239845133 + -4.27 -3.05 3.0 159ms
21 -36.74247901345 -4.09 -3.76 2.0 120ms
22 -36.74247667080 + -5.63 -3.61 2.0 128ms
23 -36.74248038552 -5.43 -4.22 2.0 102ms
24 -36.74248037168 + -7.86 -4.20 3.0 139ms
25 -36.74248065825 -6.54 -4.64 2.0 112ms
26 -36.74248065631 + -8.71 -4.65 2.0 120ms
27 -36.74248065815 -8.73 -4.92 2.0 103ms
28 -36.74248048862 + -6.77 -4.39 3.0 144ms
29 -36.74248061983 -6.88 -4.69 4.0 159ms
30 -36.74248065043 -7.51 -4.85 2.0 121ms
31 -36.74248067170 -7.67 -5.41 2.0 110ms
32 -36.74248067218 -9.32 -5.66 3.0 135ms
33 -36.74248067165 + -9.27 -5.48 2.0 128ms
34 -36.74248067255 -9.05 -5.99 3.0 136ms
35 -36.74248067261 -10.22 -5.99 3.0 136ms
36 -36.74248067268 -10.17 -6.37 1.0 94.7ms
37 -36.74248067264 + -10.48 -6.18 3.0 145ms
38 -36.74248067268 -10.43 -6.77 2.0 110ms
39 -36.74248067268 + -12.26 -6.81 2.0 106ms
40 -36.74248067268 -11.76 -7.16 2.0 110ms
41 -36.74248067268 + -12.51 -7.14 3.0 143ms
42 -36.74248067268 + -11.42 -6.74 3.0 136ms
43 -36.74248067268 -11.40 -7.24 3.0 142ms
44 -36.74248067268 -12.94 -7.60 1.0 94.8ms
45 -36.74248067268 -13.67 -7.82 3.0 123ms
46 -36.74248067268 -13.85 -7.96 3.0 130ms
47 -36.74248067268 -13.67 -8.14 2.0 134ms
48 -36.74248067268 + -13.67 -8.25 2.0 113ms
49 -36.74248067268 -13.85 -8.37 2.0 125ms
50 -36.74248067268 -13.67 -8.68 2.0 110ms
51 -36.74248067268 + -Inf -8.85 1.0 102ms
52 -36.74248067268 + -13.85 -8.92 2.0 103ms
53 -36.74248067268 + -Inf -9.32 1.0 100ms
54 -36.74248067268 + -Inf -9.42 3.0 139ms
55 -36.74248067268 + -Inf -9.69 2.0 116ms
56 -36.74248067268 + -13.85 -9.87 2.0 104ms
57 -36.74248067268 -14.15 -10.02 2.0 111ms
58 -36.74248067268 -14.15 -10.28 2.0 128ms
59 -36.74248067268 + -Inf -10.43 2.0 116ms
60 -36.74248067268 + -13.85 -10.67 2.0 107ms
61 -36.74248067268 -14.15 -10.79 3.0 139ms
62 -36.74248067268 -14.15 -10.71 2.0 130ms
63 -36.74248067268 + -Inf -11.29 2.0 116ms
64 -36.74248067268 + -Inf -11.40 2.0 120ms
65 -36.74248067268 + -14.15 -11.62 3.0 124ms
66 -36.74248067268 -14.15 -11.65 2.0 128ms
67 -36.74248067268 + -Inf -11.50 3.0 141ms
68 -36.74248067268 + -Inf -11.92 2.0 110ms
69 -36.74248067268 + -13.85 -12.22 3.0 122ms
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.73494120206 -0.88 11.0 361ms
2 -36.74049226626 -2.26 -1.36 1.0 91.8ms
3 -36.74164897467 -2.94 -2.30 2.0 102ms
4 -36.74199958264 -3.46 -2.12 4.0 162ms
5 -36.74238552944 -3.41 -2.54 1.0 92.5ms
6 -36.74244401081 -4.23 -2.84 2.0 104ms
7 -36.74247126010 -4.56 -2.82 2.0 98.4ms
8 -36.74247462426 -5.47 -3.41 1.0 100ms
9 -36.74248030830 -5.25 -3.62 3.0 134ms
10 -36.74248046239 -6.81 -3.89 2.0 109ms
11 -36.74248007649 + -6.41 -3.96 2.0 111ms
12 -36.74248066500 -6.23 -4.50 2.0 112ms
13 -36.74248067104 -8.22 -5.00 5.0 128ms
14 -36.74248067138 -9.47 -5.14 3.0 145ms
15 -36.74248067252 -8.94 -5.34 2.0 102ms
16 -36.74248067259 -10.15 -5.95 1.0 102ms
17 -36.74248067268 -10.07 -6.11 3.0 136ms
18 -36.74248067268 -11.82 -6.48 1.0 102ms
19 -36.74248067268 -11.77 -6.61 5.0 122ms
20 -36.74248067268 -11.89 -6.92 2.0 105ms
21 -36.74248067268 -12.94 -7.17 2.0 135ms
22 -36.74248067268 -14.15 -7.37 2.0 102ms
23 -36.74248067268 -13.85 -7.75 1.0 101ms
24 -36.74248067268 -13.85 -8.00 2.0 129ms
25 -36.74248067268 -14.15 -8.38 4.0 119ms
26 -36.74248067268 + -14.15 -8.66 2.0 130ms
27 -36.74248067268 + -Inf -9.07 2.0 111ms
28 -36.74248067268 + -Inf -9.42 2.0 130ms
29 -36.74248067268 + -13.85 -9.75 2.0 115ms
30 -36.74248067268 -13.85 -10.03 3.0 111ms
31 -36.74248067268 + -Inf -10.41 3.0 135ms
32 -36.74248067268 + -13.85 -10.47 2.0 130ms
33 -36.74248067268 -13.85 -11.04 2.0 117ms
34 -36.74248067268 + -Inf -11.22 2.0 135ms
35 -36.74248067268 + -Inf -11.52 1.0 95.8ms
36 -36.74248067268 + -Inf -11.83 2.0 111ms
37 -36.74248067268 + -14.15 -11.88 2.0 129ms
38 -36.74248067268 -13.67 -12.23 1.0 101ms
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.024488980448886
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.2442111139169
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.723582812210099
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).