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.73593375254 -0.88 11.0 339ms
2 -36.71952832145 + -1.79 -1.62 1.0 89.4ms
3 -8.177546947962 + 1.46 -0.34 6.0 184ms
4 -36.67605664700 1.45 -1.27 5.0 174ms
5 -36.62876808072 + -1.33 -1.51 2.0 116ms
6 -36.39003050349 + -0.62 -1.29 4.0 142ms
7 -36.73737239748 -0.46 -1.91 3.0 126ms
8 -36.74202583102 -2.33 -2.31 2.0 104ms
9 -36.74114632955 + -3.06 -2.21 2.0 124ms
10 -36.74200281620 -3.07 -2.41 2.0 113ms
11 -36.74239593466 -3.41 -2.71 2.0 111ms
12 -36.74246222576 -4.18 -3.00 2.0 105ms
13 -36.74246933480 -5.15 -3.39 2.0 126ms
14 -36.74247416896 -5.32 -3.33 2.0 104ms
15 -36.74192392591 + -3.26 -2.70 4.0 148ms
16 -36.74247574061 -3.26 -3.67 3.0 133ms
17 -36.74247976784 -5.39 -3.69 2.0 108ms
18 -36.74230625516 + -3.76 -2.95 3.0 138ms
19 -36.74248057532 -3.76 -4.28 4.0 147ms
20 -36.74248063526 -7.22 -4.50 3.0 141ms
21 -36.74248057386 + -7.21 -4.45 2.0 108ms
22 -36.74248062913 -7.26 -4.69 1.0 91.6ms
23 -36.74248066174 -7.49 -4.72 2.0 109ms
24 -36.74248066864 -8.16 -5.21 1.0 94.2ms
25 -36.74248067193 -8.48 -5.41 3.0 132ms
26 -36.74248067226 -9.49 -5.61 1.0 93.6ms
27 -36.74248066468 + -8.12 -5.10 3.0 210ms
28 -36.74248067170 -8.15 -5.54 3.0 722ms
29 -36.74248067267 -9.01 -6.10 2.0 106ms
30 -36.74248067259 + -10.09 -5.95 3.0 131ms
31 -36.74248067267 -10.09 -6.36 2.0 107ms
32 -36.74248067268 -11.03 -6.42 2.0 122ms
33 -36.74248067268 -11.47 -6.99 1.0 89.7ms
34 -36.74248067268 + -12.63 -7.09 3.0 135ms
35 -36.74248067268 -12.42 -7.14 2.0 105ms
36 -36.74248067268 -13.00 -7.62 1.0 90.4ms
37 -36.74248067268 -14.15 -7.61 2.0 122ms
38 -36.74248067268 + -12.49 -7.28 4.0 137ms
39 -36.74248067268 -12.50 -7.79 3.0 129ms
40 -36.74248067268 -13.85 -7.87 2.0 105ms
41 -36.74248067268 + -13.55 -7.67 3.0 127ms
42 -36.74248067268 -13.55 -7.96 3.0 122ms
43 -36.74248067268 + -Inf -8.12 1.0 94.2ms
44 -36.74248067268 -13.85 -8.65 2.0 106ms
45 -36.74248067268 + -Inf -8.92 3.0 138ms
46 -36.74248067268 + -14.15 -8.94 2.0 104ms
47 -36.74248067268 -14.15 -9.13 2.0 110ms
48 -36.74248067268 + -Inf -9.15 2.0 106ms
49 -36.74248067268 -14.15 -8.89 3.0 124ms
50 -36.74248067268 + -14.15 -8.97 2.0 116ms
51 -36.74248067268 + -13.85 -9.32 3.0 126ms
52 -36.74248067268 + -Inf -9.43 3.0 127ms
53 -36.74248067268 -13.85 -9.77 1.0 90.2ms
54 -36.74248067268 + -14.15 -9.81 3.0 137ms
55 -36.74248067268 -14.15 -10.07 2.0 110ms
56 -36.74248067268 + -Inf -10.39 1.0 93.9ms
57 -36.74248067268 + -Inf -10.63 2.0 124ms
58 -36.74248067268 + -Inf -10.90 1.0 94.8ms
59 -36.74248067268 + -13.85 -10.70 3.0 134ms
60 -36.74248067268 + -Inf -10.70 3.0 134ms
61 -36.74248067268 -14.15 -11.27 2.0 117ms
62 -36.74248067268 -14.15 -11.22 2.0 123ms
63 -36.74248067268 + -Inf -10.71 4.0 148ms
64 -36.74248067268 + -Inf -11.22 3.0 134ms
65 -36.74248067268 + -Inf -11.41 2.0 108ms
66 -36.74248067268 + -Inf -10.92 3.0 130ms
67 -36.74248067268 + -13.85 -11.67 3.0 131ms
68 -36.74248067268 -13.85 -11.37 3.0 135ms
69 -36.74248067268 + -Inf -12.22 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.72884006558 -0.88 12.0 336ms
2 -36.73795394623 -2.04 -1.36 1.0 89.0ms
3 -36.73916241928 -2.92 -1.69 4.0 137ms
4 -36.74197011389 -2.55 -2.02 2.0 97.5ms
5 -36.74222833672 -3.59 -2.55 2.0 94.5ms
6 -36.74237644352 -3.83 -2.41 6.0 146ms
7 -36.74244050399 -4.19 -2.92 1.0 93.9ms
8 -36.74247861754 -4.42 -3.32 2.0 102ms
9 -36.74247869772 -7.10 -3.37 3.0 116ms
10 -36.74248034539 -5.78 -3.72 4.0 106ms
11 -36.74248059049 -6.61 -3.90 1.0 96.2ms
12 -36.74248064315 -7.28 -4.46 2.0 106ms
13 -36.74248066880 -7.59 -4.57 3.0 141ms
14 -36.74248067076 -8.71 -4.69 1.0 97.2ms
15 -36.74248067207 -8.88 -5.12 2.0 109ms
16 -36.74248067249 -9.38 -5.33 2.0 102ms
17 -36.74248067261 -9.94 -5.43 1.0 97.0ms
18 -36.74248067267 -10.22 -6.00 4.0 113ms
19 -36.74248067268 -11.01 -6.18 3.0 137ms
20 -36.74248067268 -12.01 -6.45 1.0 95.1ms
21 -36.74248067268 -12.50 -6.48 3.0 127ms
22 -36.74248067268 -12.47 -6.66 1.0 96.2ms
23 -36.74248067268 -12.57 -7.36 2.0 111ms
24 -36.74248067268 + -13.85 -7.21 4.0 151ms
25 -36.74248067268 -13.67 -7.78 2.0 120ms
26 -36.74248067268 + -Inf -8.04 2.0 118ms
27 -36.74248067268 -13.85 -8.52 1.0 96.4ms
28 -36.74248067268 + -13.67 -8.93 3.0 146ms
29 -36.74248067268 -13.85 -9.41 3.0 117ms
30 -36.74248067268 + -Inf -9.28 4.0 140ms
31 -36.74248067268 + -Inf -9.54 1.0 92.7ms
32 -36.74248067268 + -13.85 -10.18 2.0 112ms
33 -36.74248067268 -13.85 -10.40 4.0 151ms
34 -36.74248067268 + -13.85 -10.67 1.0 97.2ms
35 -36.74248067268 -14.15 -10.84 3.0 127ms
36 -36.74248067268 + -Inf -11.00 1.0 96.4ms
37 -36.74248067268 -14.15 -11.21 1.0 94.0ms
38 -36.74248067268 + -Inf -11.46 5.0 116ms
39 -36.74248067268 + -13.85 -12.11 2.0 122ms
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.02448898022948
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.244211113686426
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.723669726590018
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).