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.73437371046 -0.88 11.0 383ms
2 -36.65385140711 + -1.09 -1.48 1.0 88.5ms
3 +30.48889329474 + 1.83 -0.15 8.0 247ms
4 -36.61466156812 1.83 -1.24 6.0 227ms
5 -36.44326387658 + -0.77 -1.30 4.0 154ms
6 -36.07583982608 + -0.43 -1.15 3.0 144ms
7 -36.72754489330 -0.19 -1.74 3.0 138ms
8 -36.73884445109 -1.95 -2.10 2.0 106ms
9 -36.73755670746 + -2.89 -2.12 2.0 131ms
10 -36.74121937176 -2.44 -2.24 2.0 113ms
11 -36.74123558586 -4.79 -2.26 2.0 124ms
12 -36.74244706552 -2.92 -2.73 2.0 109ms
13 -36.74245848680 -4.94 -2.83 2.0 113ms
14 -36.74246137403 -5.54 -2.95 1.0 96.4ms
15 -36.74044907478 + -2.70 -2.37 3.0 142ms
16 -36.74217071674 -2.76 -2.77 3.0 143ms
17 -36.74238170371 -3.68 -2.97 2.0 111ms
18 -36.74205900050 + -3.49 -2.75 3.0 136ms
19 -36.74244832982 -3.41 -3.16 3.0 142ms
20 -36.74247474306 -4.58 -3.49 2.0 114ms
21 -36.74247644368 -5.77 -3.70 2.0 132ms
22 -36.74247981859 -5.47 -3.95 2.0 112ms
23 -36.74248038469 -6.25 -4.19 1.0 97.0ms
24 -36.74248056656 -6.74 -4.30 2.0 132ms
25 -36.74248065820 -7.04 -4.64 2.0 114ms
26 -36.74248067068 -7.90 -4.98 1.0 94.1ms
27 -36.74248066332 + -8.13 -4.97 2.0 133ms
28 -36.74248067028 -8.16 -5.27 2.0 107ms
29 -36.74248066443 + -8.23 -5.07 3.0 143ms
30 -36.74248066343 + -9.00 -5.05 3.0 144ms
31 -36.74248067243 -8.05 -5.79 2.0 122ms
32 -36.74248067227 + -9.78 -5.64 3.0 141ms
33 -36.74248067263 -9.44 -6.11 2.0 114ms
34 -36.74248067239 + -9.62 -5.78 3.0 143ms
35 -36.74248067266 -9.57 -6.07 2.0 125ms
36 -36.74248067249 + -9.77 -5.85 3.0 143ms
37 -36.74248067268 -9.71 -6.79 3.0 118ms
38 -36.74248067268 -12.77 -6.85 3.0 161ms
39 -36.74248067268 -12.13 -7.07 2.0 115ms
40 -36.74248067268 + -12.97 -7.10 1.0 95.9ms
41 -36.74248067268 + -13.85 -7.10 2.0 102ms
42 -36.74248067267 + -10.81 -6.44 3.0 144ms
43 -36.74248067268 -10.80 -7.41 3.0 152ms
44 -36.74248067268 -12.97 -7.67 2.0 133ms
45 -36.74248067268 + -Inf -7.97 1.0 96.4ms
46 -36.74248067268 -13.85 -8.21 2.0 105ms
47 -36.74248067268 + -13.85 -7.95 3.0 145ms
48 -36.74248067268 -13.67 -8.27 3.0 142ms
49 -36.74248067268 + -Inf -8.21 2.0 125ms
50 -36.74248067268 -14.15 -8.68 2.0 109ms
51 -36.74248067268 + -14.15 -8.71 3.0 124ms
52 -36.74248067268 + -Inf -8.65 2.0 106ms
53 -36.74248067268 + -Inf -9.04 2.0 114ms
54 -36.74248067268 + -Inf -9.02 3.0 142ms
55 -36.74248067268 + -Inf -9.58 2.0 112ms
56 -36.74248067268 + -14.15 -9.59 2.0 110ms
57 -36.74248067268 + -Inf -9.27 3.0 136ms
58 -36.74248067268 + -Inf -9.63 3.0 133ms
59 -36.74248067268 + -Inf -9.72 2.0 122ms
60 -36.74248067268 -14.15 -9.78 1.0 96.9ms
61 -36.74248067268 + -Inf -10.17 2.0 108ms
62 -36.74248067268 + -Inf -10.55 3.0 142ms
63 -36.74248067268 + -14.15 -10.45 2.0 112ms
64 -36.74248067268 -14.15 -10.69 2.0 115ms
65 -36.74248067268 + -Inf -10.79 1.0 96.7ms
66 -36.74248067268 + -Inf -10.94 1.0 93.5ms
67 -36.74248067268 -14.15 -10.81 3.0 137ms
68 -36.74248067268 + -14.15 -11.52 2.0 114ms
69 -36.74248067268 + -14.15 -11.34 3.0 152ms
70 -36.74248067268 -13.85 -11.36 2.0 113ms
71 -36.74248067268 + -13.85 -11.28 3.0 132ms
72 -36.74248067268 -14.15 -12.02 2.0 115ms
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.73395235646 -0.88 11.0 384ms
2 -36.74001965493 -2.22 -1.36 1.0 90.2ms
3 -36.73793298459 + -2.68 -1.54 3.0 136ms
4 -36.74235572605 -2.35 -2.31 1.0 94.0ms
5 -36.74240510370 -4.31 -2.45 8.0 152ms
6 -36.74240619859 -5.96 -2.41 1.0 94.5ms
7 -36.74247642148 -4.15 -3.13 1.0 96.0ms
8 -36.74247825895 -5.74 -3.15 3.0 141ms
9 -36.74248035770 -5.68 -3.68 1.0 93.9ms
10 -36.74248051978 -6.79 -3.92 2.0 125ms
11 -36.74248063718 -6.93 -4.47 1.0 97.3ms
12 -36.74248067021 -7.48 -4.84 3.0 142ms
13 -36.74248067085 -9.19 -4.92 3.0 142ms
14 -36.74248067224 -8.86 -5.42 1.0 94.9ms
15 -36.74248067265 -9.39 -5.87 3.0 142ms
16 -36.74248067268 -10.53 -6.34 4.0 116ms
17 -36.74248067268 -12.06 -6.79 3.0 141ms
18 -36.74248067268 -12.46 -7.04 6.0 133ms
19 -36.74248067268 -13.11 -7.35 2.0 135ms
20 -36.74248067268 -13.45 -7.84 1.0 95.5ms
21 -36.74248067268 + -13.85 -8.21 3.0 129ms
22 -36.74248067268 + -Inf -8.34 3.0 142ms
23 -36.74248067268 + -Inf -8.72 1.0 98.3ms
24 -36.74248067268 -13.85 -8.87 2.0 114ms
25 -36.74248067268 + -Inf -9.44 2.0 102ms
26 -36.74248067268 + -Inf -9.71 4.0 148ms
27 -36.74248067268 -13.85 -9.88 2.0 105ms
28 -36.74248067268 + -13.67 -10.30 2.0 118ms
29 -36.74248067268 + -14.15 -10.44 3.0 142ms
30 -36.74248067268 -13.85 -11.05 2.0 102ms
31 -36.74248067268 + -Inf -10.95 3.0 153ms
32 -36.74248067268 + -Inf -11.31 1.0 98.5ms
33 -36.74248067268 + -Inf -11.93 3.0 127ms
34 -36.74248067268 + -13.85 -11.83 3.0 163ms
35 -36.74248067268 -14.15 -12.15 1.0 98.2ms
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.02449278661272
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.244215111988154
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.723638343557099
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).