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.73443713260 -0.88 11.0 1.31s
2 -36.65993148799 + -1.13 -1.48 1.0 242ms
3 +28.23101028305 + 1.81 -0.16 7.0 251ms
4 -36.56075651407 1.81 -1.14 7.0 233ms
5 -35.17451685101 + 0.14 -0.97 4.0 172ms
6 -36.53588798180 0.13 -1.36 4.0 163ms
7 -36.72874409709 -0.71 -1.77 2.0 118ms
8 -36.73845578621 -2.01 -2.09 2.0 103ms
9 -36.73831155996 + -3.84 -1.97 2.0 127ms
10 -36.74113579043 -2.55 -2.17 2.0 124ms
11 -36.74103487639 + -4.00 -2.32 2.0 100ms
12 -36.74243208026 -2.85 -2.77 1.0 98.6ms
13 -36.74242564036 + -5.19 -2.67 2.0 134ms
14 -36.74243098167 -5.27 -2.94 1.0 96.2ms
15 -36.74181057522 + -3.21 -2.59 2.0 125ms
16 -36.74243381450 -3.21 -3.09 3.0 137ms
17 -36.73902819858 + -2.47 -2.28 3.0 154ms
18 -36.74200431075 -2.53 -2.69 4.0 168ms
19 -36.74242958425 -3.37 -3.15 3.0 138ms
20 -36.74247483275 -4.34 -3.61 2.0 111ms
21 -36.74247920372 -5.36 -3.57 2.0 133ms
22 -36.74247811584 + -5.96 -3.80 2.0 110ms
23 -36.74248030511 -5.66 -4.20 2.0 116ms
24 -36.74248065431 -6.46 -4.33 2.0 139ms
25 -36.74248056377 + -7.04 -4.46 1.0 97.0ms
26 -36.74248059611 -7.49 -4.53 1.0 193ms
27 -36.74248065413 -7.24 -4.88 1.0 95.4ms
28 -36.74248066747 -7.87 -5.14 2.0 1.09s
29 -36.74248064337 + -7.62 -4.80 3.0 139ms
30 -36.74248066604 -7.64 -5.14 3.0 137ms
31 -36.74248065958 + -8.19 -4.90 3.0 138ms
32 -36.74248067259 -7.89 -5.66 2.0 111ms
33 -36.74248067227 + -9.49 -5.68 3.0 139ms
34 -36.74248067264 -9.42 -6.05 1.0 94.5ms
35 -36.74248067263 + -10.93 -6.16 2.0 111ms
36 -36.74248067267 -10.45 -6.16 3.0 134ms
37 -36.74248067268 -11.03 -6.46 1.0 115ms
38 -36.74248067268 -12.09 -6.51 2.0 153ms
39 -36.74248067268 -12.50 -6.45 1.0 115ms
40 -36.74248067268 -11.48 -6.98 1.0 115ms
41 -36.74248067268 + -11.21 -6.62 3.0 154ms
42 -36.74248067268 -11.15 -7.27 3.0 138ms
43 -36.74248067268 -12.94 -7.57 2.0 130ms
44 -36.74248067268 + -13.03 -7.19 3.0 136ms
45 -36.74248067268 + -13.15 -7.25 3.0 141ms
46 -36.74248067268 -12.92 -7.60 3.0 143ms
47 -36.74248067268 -13.45 -7.80 2.0 112ms
48 -36.74248067268 -14.15 -7.87 2.0 111ms
49 -36.74248067268 -14.15 -8.45 1.0 94.9ms
50 -36.74248067268 + -Inf -8.29 3.0 145ms
51 -36.74248067268 -14.15 -8.79 2.0 110ms
52 -36.74248067268 + -13.85 -8.74 2.0 129ms
53 -36.74248067268 + -14.15 -8.61 3.0 122ms
54 -36.74248067268 -13.85 -8.90 3.0 131ms
55 -36.74248067268 -13.67 -8.87 3.0 138ms
56 -36.74248067268 + -13.85 -9.47 2.0 116ms
57 -36.74248067268 + -Inf -9.65 3.0 138ms
58 -36.74248067268 -14.15 -9.78 2.0 134ms
59 -36.74248067268 -14.15 -10.35 2.0 106ms
60 -36.74248067268 + -13.85 -10.16 3.0 143ms
61 -36.74248067268 + -13.85 -10.16 2.0 110ms
62 -36.74248067268 -13.85 -10.40 2.0 148ms
63 -36.74248067268 + -14.15 -10.52 2.0 110ms
64 -36.74248067268 -14.15 -10.92 2.0 106ms
65 -36.74248067268 + -Inf -11.04 3.0 146ms
66 -36.74248067268 + -Inf -11.29 2.0 116ms
67 -36.74248067268 + -Inf -11.47 1.0 95.5ms
68 -36.74248067268 + -Inf -11.80 2.0 116ms
69 -36.74248067268 + -Inf -11.98 2.0 133ms
70 -36.74248067268 + -Inf -12.11 1.0 95.5ms
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.73372492229 -0.88 11.0 948ms
2 -36.74003658552 -2.20 -1.36 1.0 644ms
3 -36.74023456483 -3.70 -1.71 3.0 123ms
4 -36.74203615227 -2.74 -2.05 2.0 98.0ms
5 -36.74239494869 -3.45 -2.61 4.0 115ms
6 -36.74244275352 -4.32 -2.51 3.0 144ms
7 -36.74247488450 -4.49 -3.13 1.0 97.3ms
8 -36.74247938656 -5.35 -3.35 2.0 103ms
9 -36.74247989270 -6.30 -3.46 2.0 100ms
10 -36.74248060895 -6.14 -4.11 2.0 104ms
11 -36.74248064430 -7.45 -4.25 5.0 144ms
12 -36.74248065532 -7.96 -4.77 2.0 107ms
13 -36.74248067048 -7.82 -4.77 3.0 139ms
14 -36.74248067199 -8.82 -5.00 1.0 95.2ms
15 -36.74248067243 -9.36 -5.15 3.0 132ms
16 -36.74248067259 -9.77 -5.43 1.0 99.2ms
17 -36.74248067267 -10.15 -5.89 2.0 111ms
18 -36.74248067268 -10.97 -6.01 4.0 158ms
19 -36.74248067268 -11.38 -6.37 1.0 100ms
20 -36.74248067268 -11.78 -7.02 3.0 127ms
21 -36.74248067268 + -13.45 -6.83 5.0 150ms
22 -36.74248067268 -13.15 -7.12 1.0 99.0ms
23 -36.74248067268 -13.55 -7.80 2.0 120ms
24 -36.74248067268 + -13.85 -7.90 3.0 131ms
25 -36.74248067268 -14.15 -8.06 1.0 99.2ms
26 -36.74248067268 -14.15 -8.54 2.0 110ms
27 -36.74248067268 + -Inf -9.02 3.0 138ms
28 -36.74248067268 + -Inf -9.14 3.0 130ms
29 -36.74248067268 + -Inf -9.45 1.0 96.8ms
30 -36.74248067268 + -14.15 -9.77 3.0 223ms
31 -36.74248067268 -14.15 -10.12 2.0 1.06s
32 -36.74248067268 + -Inf -10.72 3.0 106ms
33 -36.74248067268 + -Inf -10.88 4.0 152ms
34 -36.74248067268 + -13.85 -11.25 4.0 111ms
35 -36.74248067268 -13.85 -11.77 2.0 121ms
36 -36.74248067268 -13.85 -11.97 3.0 136ms
37 -36.74248067268 + -Inf -12.15 1.0 94.4ms
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.024488980478985The 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.244211113948516This 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.7236066980025315Since 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).