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.73123404947 -0.88 13.0 381ms
2 -36.63881599238 + -1.03 -1.38 1.0 87.3ms
3 +27.01994189603 + 1.80 -0.16 7.0 241ms
4 -36.21514380117 1.80 -0.92 6.0 234ms
5 -36.33414164448 -0.92 -1.13 4.0 156ms
6 -34.79459001492 + 0.19 -0.92 4.0 158ms
7 -36.69064552086 0.28 -1.59 3.0 148ms
8 -36.74035792297 -1.30 -1.96 2.0 109ms
9 -36.74142368607 -2.97 -2.15 2.0 111ms
10 -36.73849990343 + -2.53 -2.00 2.0 132ms
11 -36.74138389679 -2.54 -2.18 2.0 125ms
12 -36.74225790584 -3.06 -2.42 2.0 109ms
13 -36.74230531646 -4.32 -2.67 1.0 98.1ms
14 -36.74180199788 + -3.30 -2.55 2.0 122ms
15 -36.72838110690 + -1.87 -1.99 3.0 151ms
16 -36.74239137767 -1.85 -2.62 4.0 159ms
17 -36.74197099295 + -3.38 -2.63 2.0 128ms
18 -36.74247573366 -3.30 -3.39 2.0 108ms
19 -36.74162931479 + -3.07 -2.58 4.0 169ms
20 -36.74246398039 -3.08 -3.42 4.0 171ms
21 -36.74247568146 -4.93 -3.39 2.0 111ms
22 -36.74247969124 -5.40 -3.91 1.0 96.0ms
23 -36.74247891710 + -6.11 -3.71 3.0 145ms
24 -36.74248022830 -5.88 -3.99 2.0 126ms
25 -36.74248045299 -6.65 -4.27 2.0 112ms
26 -36.74248066541 -6.67 -4.73 2.0 108ms
27 -36.74248008519 + -6.24 -4.17 3.0 152ms
28 -36.74248062398 -6.27 -4.49 4.0 154ms
29 -36.74248065851 -7.46 -4.91 2.0 112ms
30 -36.74248066873 -7.99 -5.17 2.0 117ms
31 -36.74248066132 + -8.13 -5.02 3.0 149ms
32 -36.74248065802 + -8.48 -4.97 3.0 143ms
33 -36.74248067218 -7.85 -5.58 2.0 114ms
34 -36.74248067267 -9.31 -6.07 3.0 140ms
35 -36.74248067244 + -9.65 -5.82 3.0 146ms
36 -36.74248067267 -9.65 -6.27 2.0 126ms
37 -36.74248067268 -11.04 -6.38 2.0 124ms
38 -36.74248067268 -11.52 -6.64 1.0 95.6ms
39 -36.74248067268 -12.34 -6.69 2.0 136ms
40 -36.74248067268 -11.78 -6.98 1.0 97.3ms
41 -36.74248067267 + -11.11 -6.58 3.0 132ms
42 -36.74248067268 -11.10 -7.55 3.0 265ms
43 -36.74248067268 + -12.17 -7.10 3.0 797ms
44 -36.74248067268 -13.00 -7.15 4.0 164ms
45 -36.74248067268 -12.23 -7.94 3.0 146ms
46 -36.74248067268 + -13.45 -7.63 3.0 135ms
47 -36.74248067268 -13.45 -8.21 3.0 140ms
48 -36.74248067268 + -Inf -8.30 2.0 149ms
49 -36.74248067268 -14.15 -8.39 2.0 127ms
50 -36.74248067268 -14.15 -8.65 2.0 113ms
51 -36.74248067268 + -Inf -9.02 1.0 95.6ms
52 -36.74248067268 + -14.15 -9.07 3.0 145ms
53 -36.74248067268 -14.15 -8.68 3.0 145ms
54 -36.74248067268 + -Inf -9.40 3.0 145ms
55 -36.74248067268 + -Inf -8.94 3.0 145ms
56 -36.74248067268 + -Inf -9.72 3.0 135ms
57 -36.74248067268 -14.15 -9.46 3.0 142ms
58 -36.74248067268 + -14.15 -9.83 2.0 105ms
59 -36.74248067268 + -Inf -9.65 3.0 126ms
60 -36.74248067268 + -Inf -10.04 2.0 117ms
61 -36.74248067268 + -Inf -10.50 2.0 102ms
62 -36.74248067268 + -14.15 -11.04 2.0 124ms
63 -36.74248067268 + -14.15 -10.80 3.0 137ms
64 -36.74248067268 -13.85 -11.22 2.0 108ms
65 -36.74248067268 + -Inf -10.95 3.0 139ms
66 -36.74248067268 + -13.85 -11.30 3.0 145ms
67 -36.74248067268 -14.15 -11.25 3.0 125ms
68 -36.74248067268 -14.15 -11.69 2.0 108ms
69 -36.74248067268 + -Inf -11.41 3.0 135ms
70 -36.74248067268 + -14.15 -11.76 2.0 120ms
71 -36.74248067268 -14.15 -11.58 3.0 128ms
72 -36.74248067268 + -Inf -12.39 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.73221337770 -0.88 11.0 328ms
2 -36.73973484893 -2.12 -1.36 1.0 89.3ms
3 -36.74059296637 -3.07 -1.80 4.0 118ms
4 -36.74207172559 -2.83 -2.12 3.0 109ms
5 -36.74242857809 -3.45 -2.73 2.0 104ms
6 -36.74241853605 + -5.00 -2.51 6.0 144ms
7 -36.74247144903 -4.28 -3.01 1.0 91.8ms
8 -36.74247908017 -5.12 -3.51 1.0 92.0ms
9 -36.74247928995 -6.68 -3.42 4.0 136ms
10 -36.74248053580 -5.90 -4.02 1.0 93.0ms
11 -36.74248063937 -6.98 -4.26 2.0 125ms
12 -36.74248060771 + -7.50 -4.39 2.0 103ms
13 -36.74248066345 -7.25 -4.76 1.0 94.1ms
14 -36.74248067201 -8.07 -5.10 2.0 127ms
15 -36.74248067259 -9.23 -5.39 2.0 101ms
16 -36.74248067258 + -11.01 -5.63 3.0 124ms
17 -36.74248067262 -10.48 -5.77 1.0 95.6ms
18 -36.74248067267 -10.27 -6.16 2.0 103ms
19 -36.74248067268 -10.93 -6.76 3.0 133ms
20 -36.74248067268 -11.96 -7.06 5.0 112ms
21 -36.74248067268 + -13.55 -7.12 2.0 128ms
22 -36.74248067268 -13.37 -7.25 2.0 110ms
23 -36.74248067268 -13.15 -7.59 1.0 97.6ms
24 -36.74248067268 + -14.15 -7.74 2.0 110ms
25 -36.74248067268 -13.55 -8.29 2.0 105ms
26 -36.74248067268 + -14.15 -8.46 2.0 130ms
27 -36.74248067268 -14.15 -8.63 1.0 96.5ms
28 -36.74248067268 -14.15 -9.01 2.0 97.6ms
29 -36.74248067268 + -13.85 -9.28 2.0 121ms
30 -36.74248067268 -13.67 -9.52 1.0 96.3ms
31 -36.74248067268 + -14.15 -10.19 3.0 124ms
32 -36.74248067268 + -Inf -10.18 3.0 147ms
33 -36.74248067268 + -Inf -10.37 1.0 96.3ms
34 -36.74248067268 + -Inf -10.73 2.0 118ms
35 -36.74248067268 + -14.15 -11.14 1.0 92.8ms
36 -36.74248067268 -14.15 -11.17 3.0 137ms
37 -36.74248067268 + -14.15 -11.42 1.0 95.3ms
38 -36.74248067268 -14.15 -11.84 2.0 108ms
39 -36.74248067268 + -Inf -12.15 2.0 128ms
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.024488980393166
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.24421111385837
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.72358138249738
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).