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.73229578312 -0.88 11.0 1.29s
2 -36.71475643211 + -1.76 -1.49 1.0 375ms
3 -10.00589823584 + 1.43 -0.35 7.0 196ms
4 -35.43129498066 1.41 -0.88 6.0 1.17s
5 -35.84342058730 -0.38 -1.03 3.0 146ms
6 -35.31245081058 + -0.27 -0.97 6.0 182ms
7 -36.73006989703 0.15 -1.76 3.0 132ms
8 -36.73339793553 -2.48 -1.82 2.0 113ms
9 -36.73974248754 -2.20 -2.16 2.0 110ms
10 -36.74015294330 -3.39 -2.12 2.0 125ms
11 -36.74227123383 -2.67 -2.49 1.0 92.0ms
12 -36.74226841558 + -5.55 -2.48 2.0 113ms
13 -36.74245487408 -3.73 -3.08 1.0 96.4ms
14 -36.74242498571 + -4.52 -3.05 3.0 164ms
15 -36.74246791371 -4.37 -3.32 1.0 93.6ms
16 -36.74221558998 + -3.60 -2.85 3.0 144ms
17 -36.74247611825 -3.58 -3.53 3.0 134ms
18 -36.74194250957 + -3.27 -2.70 4.0 161ms
19 -36.74247769882 -3.27 -3.64 3.0 149ms
20 -36.74247975089 -5.69 -3.84 2.0 109ms
21 -36.74248046684 -6.15 -4.34 2.0 127ms
22 -36.74248064890 -6.74 -4.58 3.0 118ms
23 -36.74248066136 -7.90 -4.79 2.0 133ms
24 -36.74248066522 -8.41 -4.87 2.0 101ms
25 -36.74248067201 -8.17 -5.36 2.0 104ms
26 -36.74248067225 -9.63 -5.38 3.0 141ms
27 -36.74248067219 + -10.26 -5.30 2.0 119ms
28 -36.74248066959 + -8.58 -5.28 2.0 111ms
29 -36.74248067156 -8.71 -5.50 3.0 134ms
30 -36.74248067146 + -10.02 -5.52 3.0 133ms
31 -36.74248067181 -9.46 -5.59 2.0 109ms
32 -36.74248067258 -9.11 -5.85 2.0 113ms
33 -36.74248067268 -10.02 -6.42 1.0 93.6ms
34 -36.74248067267 + -11.01 -6.36 3.0 142ms
35 -36.74248067268 -10.98 -6.72 2.0 101ms
36 -36.74248067268 + -11.66 -6.67 2.0 132ms
37 -36.74248067268 -11.34 -7.29 2.0 108ms
38 -36.74248067268 -12.87 -7.36 3.0 142ms
39 -36.74248067268 + -13.45 -7.48 2.0 113ms
40 -36.74248067268 -14.15 -7.53 1.0 94.0ms
41 -36.74248067268 -13.25 -7.86 1.0 97.7ms
42 -36.74248067268 + -13.67 -7.71 3.0 134ms
43 -36.74248067268 + -13.25 -7.55 3.0 139ms
44 -36.74248067268 -13.15 -7.64 3.0 136ms
45 -36.74248067268 -13.85 -8.07 2.0 113ms
46 -36.74248067268 -14.15 -8.48 2.0 109ms
47 -36.74248067268 + -13.85 -8.34 3.0 141ms
48 -36.74248067268 -14.15 -8.52 2.0 125ms
49 -36.74248067268 + -Inf -8.93 2.0 104ms
50 -36.74248067268 + -Inf -8.92 2.0 132ms
51 -36.74248067268 -14.15 -8.87 2.0 119ms
52 -36.74248067268 + -Inf -9.18 2.0 107ms
53 -36.74248067268 + -Inf -9.62 3.0 124ms
54 -36.74248067268 + -Inf -9.33 3.0 142ms
55 -36.74248067268 + -Inf -9.88 3.0 148ms
56 -36.74248067268 + -Inf -9.91 2.0 109ms
57 -36.74248067268 + -13.85 -10.26 1.0 97.8ms
58 -36.74248067268 -14.15 -10.44 3.0 136ms
59 -36.74248067268 -14.15 -10.41 2.0 115ms
60 -36.74248067268 -14.15 -11.05 1.0 93.6ms
61 -36.74248067268 + -14.15 -11.06 3.0 152ms
62 -36.74248067268 + -Inf -11.22 1.0 93.1ms
63 -36.74248067268 -14.15 -11.39 3.0 135ms
64 -36.74248067268 + -13.85 -11.72 1.0 93.4ms
65 -36.74248067268 -14.15 -11.93 3.0 139ms
66 -36.74248067268 + -Inf -11.92 4.0 138ms
67 -36.74248067268 + -Inf -12.18 2.0 109ms
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.73357150311 -0.88 13.0 973ms
2 -36.74029448850 -2.17 -1.36 1.0 647ms
3 -36.74075863621 -3.33 -1.77 4.0 135ms
4 -36.74222673872 -2.83 -2.17 1.0 90.7ms
5 -36.74241698291 -3.72 -2.73 5.0 121ms
6 -36.74241740071 -6.38 -2.48 3.0 132ms
7 -36.74247629855 -4.23 -3.19 1.0 96.9ms
8 -36.74247906030 -5.56 -3.31 2.0 204ms
9 -36.74248024424 -5.93 -3.48 2.0 109ms
10 -36.74248054156 -6.53 -3.90 1.0 1.07s
11 -36.74248061106 -7.16 -4.01 6.0 127ms
12 -36.74248050884 + -6.99 -4.21 2.0 109ms
13 -36.74248066763 -6.80 -4.67 1.0 95.2ms
14 -36.74248067050 -8.54 -4.83 3.0 132ms
15 -36.74248067198 -8.83 -4.90 1.0 94.8ms
16 -36.74248066687 + -8.29 -4.85 2.0 120ms
17 -36.74248067227 -8.27 -5.43 1.0 94.8ms
18 -36.74248067237 -10.02 -5.61 3.0 138ms
19 -36.74248067262 -9.60 -5.91 1.0 103ms
20 -36.74248067268 -10.21 -6.17 2.0 152ms
21 -36.74248067268 -11.84 -6.53 1.0 114ms
22 -36.74248067268 -11.94 -7.00 2.0 143ms
23 -36.74248067268 -12.73 -7.52 2.0 152ms
24 -36.74248067268 -14.15 -7.59 4.0 147ms
25 -36.74248067268 -13.67 -7.78 1.0 95.5ms
26 -36.74248067268 + -13.67 -8.00 2.0 101ms
27 -36.74248067268 -13.67 -8.30 3.0 116ms
28 -36.74248067268 + -14.15 -8.74 2.0 129ms
29 -36.74248067268 + -Inf -9.32 2.0 114ms
30 -36.74248067268 -14.15 -9.24 4.0 156ms
31 -36.74248067268 + -13.85 -9.75 2.0 101ms
32 -36.74248067268 -13.85 -9.99 3.0 137ms
33 -36.74248067268 + -14.15 -10.28 2.0 100ms
34 -36.74248067268 + -Inf -10.44 2.0 125ms
35 -36.74248067268 -13.85 -10.78 1.0 95.3ms
36 -36.74248067268 + -13.85 -11.15 2.0 130ms
37 -36.74248067268 + -Inf -11.50 1.0 96.1ms
38 -36.74248067268 + -Inf -11.79 3.0 140ms
39 -36.74248067268 + -Inf -12.08 2.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
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.0244896539103The 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.24421182133435This 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.723580668406999Since 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).