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.73310452359 -0.88 11.0 1.37s
2 -36.72208333559 + -1.96 -1.54 1.0 276ms
3 -17.36531823474 + 1.29 -0.42 6.0 236ms
4 -35.57727373762 1.26 -0.91 5.0 177ms
5 -36.40112054894 -0.08 -1.20 4.0 146ms
6 -35.56056015107 + -0.08 -1.03 4.0 167ms
7 -36.69748921447 0.06 -1.66 3.0 127ms
8 -36.73861194594 -1.39 -2.02 2.0 120ms
9 -36.74085723463 -2.65 -2.20 2.0 136ms
10 -36.74029286020 + -3.25 -2.02 2.0 123ms
11 -36.74199826165 -2.77 -2.41 2.0 117ms
12 -36.74208163819 -4.08 -2.44 2.0 116ms
13 -36.74246024797 -3.42 -3.01 1.0 91.9ms
14 -36.74246201842 -5.75 -3.24 2.0 124ms
15 -36.74247969368 -4.75 -3.60 1.0 97.2ms
16 -36.74247594358 + -5.43 -3.62 2.0 124ms
17 -36.74246120469 + -4.83 -3.41 3.0 130ms
18 -36.74244929604 + -4.92 -3.30 4.0 145ms
19 -36.74247905454 -4.53 -3.88 2.0 122ms
20 -36.74247965344 -6.22 -4.01 2.0 115ms
21 -36.74248037366 -6.14 -4.26 2.0 121ms
22 -36.74248057688 -6.69 -4.45 3.0 114ms
23 -36.74248064636 -7.16 -4.71 2.0 258ms
24 -36.74248066821 -7.66 -4.79 2.0 100ms
25 -36.74248067145 -8.49 -5.18 1.0 1.31s
26 -36.74248067194 -9.31 -5.36 2.0 130ms
27 -36.74248067224 -9.53 -5.46 1.0 90.8ms
28 -36.74248066969 + -8.59 -5.18 2.0 116ms
29 -36.74248067229 -8.58 -5.71 2.0 105ms
30 -36.74248066764 + -8.33 -5.20 3.0 140ms
31 -36.74248067246 -8.32 -5.87 3.0 130ms
32 -36.74248067263 -9.77 -6.08 2.0 109ms
33 -36.74248067255 + -10.10 -5.95 3.0 139ms
34 -36.74248067268 -9.88 -6.41 2.0 127ms
35 -36.74248067268 -11.63 -6.99 2.0 128ms
36 -36.74248067268 + -12.51 -6.94 2.0 146ms
37 -36.74248067268 + -11.78 -6.73 2.0 146ms
38 -36.74248067268 -11.83 -7.01 2.0 138ms
39 -36.74248067268 -12.19 -7.35 2.0 120ms
40 -36.74248067268 -13.11 -7.67 3.0 133ms
41 -36.74248067268 -13.37 -7.69 3.0 132ms
42 -36.74248067268 + -13.55 -7.71 1.0 90.5ms
43 -36.74248067268 + -13.37 -7.54 3.0 130ms
44 -36.74248067268 -13.11 -8.05 2.0 106ms
45 -36.74248067268 + -12.87 -7.48 3.0 141ms
46 -36.74248067268 -12.92 -7.94 4.0 153ms
47 -36.74248067268 -14.15 -8.48 2.0 105ms
48 -36.74248067268 + -Inf -8.58 2.0 124ms
49 -36.74248067268 + -Inf -8.87 1.0 90.9ms
50 -36.74248067268 -14.15 -9.14 2.0 100ms
51 -36.74248067268 + -Inf -9.40 2.0 130ms
52 -36.74248067268 + -Inf -9.37 2.0 106ms
53 -36.74248067268 + -Inf -9.59 2.0 101ms
54 -36.74248067268 -13.85 -9.76 2.0 126ms
55 -36.74248067268 + -13.85 -9.59 2.0 121ms
56 -36.74248067268 + -Inf -10.25 2.0 100ms
57 -36.74248067268 + -Inf -10.52 2.0 124ms
58 -36.74248067268 -13.85 -10.23 3.0 146ms
59 -36.74248067268 + -Inf -10.45 3.0 130ms
60 -36.74248067268 + -13.85 -10.68 2.0 106ms
61 -36.74248067268 + -Inf -10.67 2.0 111ms
62 -36.74248067268 + -Inf -11.06 1.0 90.6ms
63 -36.74248067268 -13.85 -11.42 3.0 130ms
64 -36.74248067268 + -14.15 -11.55 2.0 125ms
65 -36.74248067268 + -14.15 -11.22 3.0 127ms
66 -36.74248067268 + -Inf -11.78 3.0 125ms
67 -36.74248067268 + -Inf -11.71 3.0 136ms
68 -36.74248067268 + -Inf -11.89 2.0 111ms
69 -36.74248067268 + -Inf -11.83 3.0 129ms
70 -36.74248067268 + -Inf -12.24 1.0 91.4ms
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.73326513035 -0.88 11.0 936ms
2 -36.74015160140 -2.16 -1.37 1.0 643ms
3 -36.74050165161 -3.46 -1.71 3.0 145ms
4 -36.74223787070 -2.76 -2.20 1.0 87.6ms
5 -36.74235206449 -3.94 -2.49 4.0 110ms
6 -36.74242315059 -4.15 -2.50 5.0 127ms
7 -36.74247105505 -4.32 -3.19 1.0 89.2ms
8 -36.74247732764 -5.20 -3.15 4.0 145ms
9 -36.74247957116 -5.65 -3.47 1.0 90.3ms
10 -36.74248038335 -6.09 -3.97 1.0 96.3ms
11 -36.74248063512 -6.60 -4.36 3.0 130ms
12 -36.74248066880 -7.47 -4.57 5.0 117ms
13 -36.74248066926 -9.34 -4.97 1.0 92.3ms
14 -36.74248067216 -8.54 -5.27 2.0 132ms
15 -36.74248067260 -9.35 -5.45 1.0 92.3ms
16 -36.74248067267 -10.18 -5.85 2.0 117ms
17 -36.74248067268 -11.01 -6.30 4.0 107ms
18 -36.74248067268 -11.72 -6.47 2.0 132ms
19 -36.74248067268 -11.99 -6.89 2.0 101ms
20 -36.74248067268 -12.87 -7.16 3.0 128ms
21 -36.74248067268 -12.79 -7.69 3.0 107ms
22 -36.74248067268 + -13.85 -7.71 4.0 151ms
23 -36.74248067268 -13.67 -8.28 2.0 103ms
24 -36.74248067268 + -Inf -8.45 3.0 138ms
25 -36.74248067268 + -13.85 -8.81 2.0 113ms
26 -36.74248067268 -13.67 -9.05 2.0 125ms
27 -36.74248067268 + -13.67 -9.33 1.0 97.7ms
28 -36.74248067268 -14.15 -9.92 2.0 102ms
29 -36.74248067268 + -Inf -10.24 3.0 143ms
30 -36.74248067268 + -Inf -10.35 2.0 97.7ms
31 -36.74248067268 -14.15 -10.29 2.0 120ms
32 -36.74248067268 + -14.15 -10.67 1.0 92.5ms
33 -36.74248067268 -14.15 -11.30 2.0 118ms
34 -36.74248067268 + -Inf -11.62 3.0 147ms
35 -36.74248067268 + -14.15 -11.67 2.0 132ms
36 -36.74248067268 + -Inf -12.10 1.0 92.5ms
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.0244891268031The 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.244211267650314This 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.723615615286123Since 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).