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.73102093623 -0.88 11.0 1.50s
2 -36.27736853257 + -0.34 -1.20 1.0 262ms
3 +138.2230318173 + 2.24 0.07 41.0 513ms
4 -32.83661323521 2.23 -0.71 11.0 310ms
5 -36.00930200948 0.50 -1.03 4.0 171ms
6 -29.00073631365 + 0.85 -0.61 4.0 153ms
7 -36.66780890160 0.88 -1.40 4.0 157ms
8 -36.73797669855 -1.15 -1.93 2.0 103ms
9 -36.73349362184 + -2.35 -1.85 2.0 119ms
10 -36.74126629260 -2.11 -2.05 2.0 112ms
11 -36.73907701013 + -2.66 -2.10 1.0 89.3ms
12 -36.73874627424 + -3.48 -2.19 1.0 94.6ms
13 -36.73959603380 -3.07 -2.24 1.0 90.6ms
14 -36.74088336604 -2.89 -2.37 1.0 96.2ms
15 -36.73958619105 + -2.89 -2.31 1.0 90.3ms
16 -36.28721991287 + -0.34 -1.24 4.0 159ms
17 -36.74202249741 -0.34 -2.59 4.0 154ms
18 -36.74202075787 + -5.76 -2.73 2.0 120ms
19 -36.73652590285 + -2.26 -2.17 3.0 139ms
20 -36.74247341279 -2.23 -3.39 3.0 135ms
21 -36.74245088471 + -4.65 -3.20 3.0 133ms
22 -36.74247552857 -4.61 -3.59 2.0 111ms
23 -36.74248022399 -5.33 -3.93 3.0 106ms
24 -36.74248030491 -7.09 -4.07 2.0 128ms
25 -36.74248047022 -6.78 -3.96 2.0 107ms
26 -36.74248065772 -6.73 -4.62 1.0 94.9ms
27 -36.74248049461 + -6.79 -4.35 3.0 133ms
28 -36.74248031041 + -6.73 -4.27 2.0 120ms
29 -36.74248059799 -6.54 -4.59 3.0 130ms
30 -36.74248066746 -7.16 -4.91 2.0 110ms
31 -36.74248064115 + -7.58 -4.78 3.0 123ms
32 -36.74248066432 -7.64 -5.07 2.0 229ms
33 -36.74248066912 -8.32 -5.10 2.0 106ms
34 -36.74248067230 -8.50 -5.69 2.0 1.27s
35 -36.74248067228 + -10.70 -5.61 3.0 137ms
36 -36.74248067252 -9.62 -5.88 1.0 89.2ms
37 -36.74248067266 -9.84 -6.16 2.0 104ms
38 -36.74248067266 + -11.30 -6.27 3.0 120ms
39 -36.74248067267 -10.86 -6.39 2.0 103ms
40 -36.74248067268 -11.03 -6.65 2.0 121ms
41 -36.74248067268 + -11.29 -6.62 2.0 105ms
42 -36.74248067268 -11.19 -7.03 2.0 105ms
43 -36.74248067268 + -11.68 -6.82 3.0 129ms
44 -36.74248067268 + -11.49 -6.62 3.0 137ms
45 -36.74248067268 -11.26 -7.58 3.0 151ms
46 -36.74248067268 + -12.75 -7.35 3.0 159ms
47 -36.74248067268 -12.73 -7.76 2.0 133ms
48 -36.74248067268 -13.67 -8.00 2.0 122ms
49 -36.74248067268 + -14.15 -8.38 2.0 145ms
50 -36.74248067268 -14.15 -8.86 2.0 104ms
51 -36.74248067268 + -14.15 -8.51 3.0 145ms
52 -36.74248067268 -13.85 -9.16 2.0 120ms
53 -36.74248067268 + -14.15 -8.80 3.0 139ms
54 -36.74248067268 + -14.15 -9.44 3.0 136ms
55 -36.74248067268 -14.15 -9.63 2.0 132ms
56 -36.74248067268 + -Inf -9.28 3.0 138ms
57 -36.74248067268 + -14.15 -9.44 3.0 128ms
58 -36.74248067268 + -Inf -10.18 2.0 108ms
59 -36.74248067268 -13.85 -10.11 3.0 136ms
60 -36.74248067268 + -14.15 -10.69 2.0 104ms
61 -36.74248067268 + -14.15 -10.76 3.0 130ms
62 -36.74248067268 + -Inf -10.95 2.0 103ms
63 -36.74248067268 -14.15 -11.20 1.0 89.7ms
64 -36.74248067268 -14.15 -11.04 2.0 119ms
65 -36.74248067268 + -14.15 -10.82 3.0 129ms
66 -36.74248067268 -14.15 -11.52 3.0 121ms
67 -36.74248067268 + -13.85 -11.38 2.0 126ms
68 -36.74248067268 + -13.85 -11.88 2.0 116ms
69 -36.74248067268 -13.85 -11.36 3.0 128ms
70 -36.74248067268 -14.15 -11.99 3.0 134ms
71 -36.74248067268 + -14.15 -12.06 2.0 122ms
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.73171339542 -0.88 10.0 944ms
2 -36.73967821558 -2.10 -1.37 1.0 704ms
3 -36.74118240188 -2.82 -1.85 3.0 118ms
4 -36.74214145814 -3.02 -2.10 1.0 87.0ms
5 -36.74233555660 -3.71 -2.54 6.0 136ms
6 -36.74242657389 -4.04 -2.46 3.0 118ms
7 -36.74246421900 -4.42 -3.18 1.0 88.7ms
8 -36.74247900611 -4.83 -3.47 3.0 134ms
9 -36.74247977983 -6.11 -3.61 1.0 94.5ms
10 -36.74248037194 -6.23 -3.92 1.0 91.0ms
11 -36.74248063439 -6.58 -4.49 2.0 129ms
12 -36.74248067211 -7.42 -5.14 7.0 126ms
13 -36.74248067183 + -9.55 -5.24 3.0 142ms
14 -36.74248067263 -9.10 -5.84 1.0 91.3ms
15 -36.74248067267 -10.34 -5.99 3.0 139ms
16 -36.74248067267 + -11.32 -6.30 2.0 96.9ms
17 -36.74248067268 -10.88 -6.70 3.0 111ms
18 -36.74248067268 -12.04 -7.08 3.0 133ms
19 -36.74248067268 -14.15 -7.40 3.0 104ms
20 -36.74248067268 -13.45 -7.52 3.0 134ms
21 -36.74248067268 -13.55 -8.00 2.0 100ms
22 -36.74248067268 + -13.85 -8.31 3.0 117ms
23 -36.74248067268 + -14.15 -8.75 2.0 124ms
24 -36.74248067268 -13.85 -9.03 3.0 117ms
25 -36.74248067268 + -Inf -9.26 3.0 129ms
26 -36.74248067268 + -Inf -9.61 3.0 116ms
27 -36.74248067268 + -Inf -10.08 2.0 100ms
28 -36.74248067268 -13.85 -10.45 3.0 137ms
29 -36.74248067268 + -14.15 -10.46 5.0 118ms
30 -36.74248067268 + -14.15 -11.11 1.0 91.3ms
31 -36.74248067268 + -Inf -11.27 4.0 149ms
32 -36.74248067268 + -Inf -11.59 1.0 91.0ms
33 -36.74248067268 + -Inf -11.91 2.0 128ms
34 -36.74248067268 + -13.85 -12.04 2.0 100ms
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.024489088058615The 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.24421122695233This 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.723569123549863Since 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).