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.73260571261 -0.88 12.0 340ms
2 -36.73426117503 -2.78 -1.62 1.0 92.6ms
3 -30.13148225263 + 0.82 -0.65 4.0 161ms
4 -35.95141590417 0.76 -1.01 4.0 160ms
5 -36.61536145233 -0.18 -1.46 3.0 134ms
6 -35.93275514924 + -0.17 -1.11 4.0 144ms
7 -36.71642685032 -0.11 -1.66 3.0 141ms
8 -36.74047107031 -1.62 -2.11 2.0 102ms
9 -36.73925089459 + -2.91 -2.08 2.0 126ms
10 -36.74222963528 -2.53 -2.30 1.0 93.2ms
11 -36.74215279477 + -4.11 -2.61 2.0 101ms
12 -36.74226658969 -3.94 -2.76 2.0 107ms
13 -36.74245446406 -3.73 -3.17 1.0 96.7ms
14 -36.74246526339 -4.97 -3.31 2.0 131ms
15 -36.74245779559 + -5.13 -3.19 2.0 121ms
16 -36.74232935071 + -3.89 -2.95 3.0 138ms
17 -36.74247757010 -3.83 -3.61 3.0 131ms
18 -36.74245735104 + -4.69 -3.37 3.0 129ms
19 -36.74246813206 -4.97 -3.49 3.0 136ms
20 -36.74248051488 -4.91 -4.29 2.0 109ms
21 -36.74247982029 + -6.16 -3.96 3.0 152ms
22 -36.74247961720 + -6.69 -3.91 3.0 134ms
23 -36.74248065300 -5.98 -4.70 3.0 122ms
24 -36.74248065410 -8.96 -4.80 2.0 129ms
25 -36.74248066881 -7.83 -5.10 2.0 104ms
26 -36.74248066985 -8.98 -5.04 2.0 108ms
27 -36.74248067209 -8.65 -5.55 1.0 93.1ms
28 -36.74248064709 + -7.60 -4.83 4.0 156ms
29 -36.74248067211 -7.60 -5.65 4.0 158ms
30 -36.74248066824 + -8.41 -5.21 3.0 135ms
31 -36.74248067252 -8.37 -5.82 3.0 132ms
32 -36.74248067264 -9.90 -5.96 2.0 116ms
33 -36.74248067000 + -8.58 -5.33 3.0 142ms
34 -36.74248067251 -8.60 -5.81 3.0 137ms
35 -36.74248067268 -9.79 -6.26 2.0 112ms
36 -36.74248067267 + -11.00 -6.38 3.0 120ms
37 -36.74248067268 -10.89 -6.67 1.0 94.0ms
38 -36.74248067268 -11.60 -6.91 2.0 128ms
39 -36.74248067268 + -13.85 -7.23 1.0 94.9ms
40 -36.74248067268 -12.75 -7.56 2.0 105ms
41 -36.74248067268 + -11.85 -6.95 4.0 148ms
42 -36.74248067268 -11.87 -7.44 3.0 135ms
43 -36.74248067268 -12.85 -7.94 3.0 129ms
44 -36.74248067268 -14.15 -8.23 2.0 102ms
45 -36.74248067268 + -13.85 -8.01 3.0 133ms
46 -36.74248067268 -14.15 -8.30 3.0 127ms
47 -36.74248067268 -14.15 -8.57 2.0 104ms
48 -36.74248067268 + -14.15 -8.93 2.0 111ms
49 -36.74248067268 -13.85 -9.01 2.0 125ms
50 -36.74248067268 + -Inf -9.02 3.0 113ms
51 -36.74248067268 + -Inf -9.35 1.0 98.9ms
52 -36.74248067268 + -Inf -9.71 3.0 118ms
53 -36.74248067268 + -13.85 -9.52 3.0 132ms
54 -36.74248067268 -13.85 -9.72 2.0 121ms
55 -36.74248067268 + -Inf -9.38 3.0 150ms
56 -36.74248067268 + -Inf -9.70 3.0 140ms
57 -36.74248067268 + -13.85 -10.08 2.0 109ms
58 -36.74248067268 -13.85 -10.20 1.0 97.0ms
59 -36.74248067268 + -Inf -10.42 3.0 132ms
60 -36.74248067268 + -14.15 -10.57 2.0 138ms
61 -36.74248067268 -14.15 -10.75 2.0 118ms
62 -36.74248067268 + -Inf -10.84 3.0 137ms
63 -36.74248067268 + -14.15 -11.02 1.0 128ms
64 -36.74248067268 + -Inf -11.17 2.0 115ms
65 -36.74248067268 + -Inf -11.20 2.0 109ms
66 -36.74248067268 + -14.15 -11.57 1.0 100ms
67 -36.74248067268 -13.55 -11.12 3.0 156ms
68 -36.74248067268 + -13.85 -11.53 3.0 156ms
69 -36.74248067268 + -Inf -12.14 2.0 108ms
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.73386796139 -0.88 12.0 400ms
2 -36.74014855843 -2.20 -1.36 1.0 95.8ms
3 -36.73953580366 + -3.21 -1.62 4.0 145ms
4 -36.74230116214 -2.56 -2.26 1.0 95.5ms
5 -36.74240873499 -3.97 -2.56 8.0 161ms
6 -36.74244636568 -4.42 -2.57 2.0 115ms
7 -36.74247109434 -4.61 -2.92 1.0 100ms
8 -36.74247644584 -5.27 -3.13 1.0 100ms
9 -36.74247886418 -5.62 -3.30 2.0 119ms
10 -36.74248041227 -5.81 -4.15 1.0 98.1ms
11 -36.74248062083 -6.68 -4.26 5.0 179ms
12 -36.74248066666 -7.34 -4.65 2.0 107ms
13 -36.74248066874 -8.68 -4.98 2.0 138ms
14 -36.74248067237 -8.44 -5.35 2.0 109ms
15 -36.74248067261 -9.62 -5.79 3.0 123ms
16 -36.74248067267 -10.24 -6.05 3.0 150ms
17 -36.74248067267 -11.35 -6.21 2.0 110ms
18 -36.74248067268 -10.99 -6.46 1.0 104ms
19 -36.74248067268 -11.96 -6.82 2.0 105ms
20 -36.74248067268 -11.98 -7.09 7.0 163ms
21 -36.74248067268 -13.37 -7.24 2.0 122ms
22 -36.74248067268 -13.25 -7.54 2.0 132ms
23 -36.74248067268 -13.67 -7.85 2.0 115ms
24 -36.74248067268 + -Inf -8.20 2.0 137ms
25 -36.74248067268 + -13.85 -8.53 2.0 117ms
26 -36.74248067268 -13.85 -8.63 2.0 142ms
27 -36.74248067268 + -Inf -9.43 1.0 101ms
28 -36.74248067268 + -Inf -9.30 4.0 188ms
29 -36.74248067268 + -Inf -9.74 2.0 119ms
30 -36.74248067268 + -13.85 -10.04 2.0 129ms
31 -36.74248067268 -14.15 -10.46 1.0 104ms
32 -36.74248067268 + -Inf -10.70 2.0 158ms
33 -36.74248067268 -14.15 -11.06 1.0 126ms
34 -36.74248067268 + -Inf -11.19 3.0 149ms
35 -36.74248067268 + -Inf -11.61 2.0 139ms
36 -36.74248067268 + -13.85 -11.82 1.0 106ms
37 -36.74248067268 -13.85 -12.36 2.0 133ms
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.02449311217027
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.24421545396037
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.723576776606016
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).