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.73207104408 -0.88 11.0 329ms
2 -36.71853256641 + -1.87 -1.59 1.0 88.4ms
3 -9.743597454124 + 1.43 -0.35 7.0 192ms
4 -36.53208899782 1.43 -1.15 6.0 182ms
5 -36.32488383184 + -0.68 -1.24 3.0 125ms
6 -36.40557517805 -1.09 -1.29 5.0 154ms
7 -36.72845708057 -0.49 -1.74 2.0 115ms
8 -36.74146595893 -1.89 -2.19 2.0 104ms
9 -36.73229405477 + -2.04 -1.96 2.0 120ms
10 -36.73917233251 -2.16 -1.97 3.0 129ms
11 -36.74217032245 -2.52 -2.42 2.0 117ms
12 -36.74241101710 -3.62 -2.70 2.0 109ms
13 -36.74245360665 -4.37 -2.99 1.0 90.2ms
14 -36.74247360390 -4.70 -3.22 2.0 103ms
15 -36.74219278778 + -3.55 -2.78 3.0 135ms
16 -36.74223150004 -4.41 -2.79 4.0 146ms
17 -36.74243041781 -3.70 -3.14 3.0 120ms
18 -36.74244520183 -4.83 -3.25 3.0 122ms
19 -36.74246384566 -4.73 -3.43 2.0 107ms
20 -36.74245174710 + -4.92 -3.33 3.0 126ms
21 -36.74248013191 -4.55 -4.10 2.0 105ms
22 -36.74248041388 -6.55 -4.14 3.0 140ms
23 -36.74248034582 + -7.17 -4.28 2.0 108ms
24 -36.74248065486 -6.51 -4.66 2.0 99.0ms
25 -36.74248065646 -8.80 -4.75 2.0 102ms
26 -36.74248066821 -7.93 -4.91 2.0 108ms
27 -36.74248066366 + -8.34 -4.90 2.0 116ms
28 -36.74248067146 -8.11 -5.36 1.0 90.0ms
29 -36.74248065606 + -7.81 -4.93 3.0 151ms
30 -36.74248067210 -7.79 -5.63 3.0 132ms
31 -36.74248067053 + -8.80 -5.30 3.0 131ms
32 -36.74248067266 -8.67 -6.13 3.0 124ms
33 -36.74248067263 + -10.58 -6.09 3.0 133ms
34 -36.74248067253 + -9.99 -5.91 3.0 126ms
35 -36.74248067240 + -9.88 -5.81 3.0 132ms
36 -36.74248067266 -9.57 -6.38 3.0 129ms
37 -36.74248067268 -10.75 -6.86 2.0 116ms
38 -36.74248067268 + -12.03 -6.86 2.0 134ms
39 -36.74248067268 -11.68 -7.39 2.0 118ms
40 -36.74248067268 -13.67 -7.56 3.0 141ms
41 -36.74248067268 + -12.64 -7.31 3.0 134ms
42 -36.74248067268 -12.64 -7.83 3.0 137ms
43 -36.74248067268 -13.85 -7.68 3.0 137ms
44 -36.74248067268 -13.85 -8.10 2.0 113ms
45 -36.74248067268 + -13.85 -8.09 3.0 148ms
46 -36.74248067268 -13.85 -8.51 1.0 97.9ms
47 -36.74248067268 + -Inf -8.30 3.0 145ms
48 -36.74248067268 + -14.15 -8.60 2.0 113ms
49 -36.74248067268 -14.15 -8.89 2.0 116ms
50 -36.74248067268 + -Inf -8.58 3.0 147ms
51 -36.74248067268 + -14.15 -9.32 3.0 152ms
52 -36.74248067268 + -Inf -9.48 2.0 125ms
53 -36.74248067268 -14.15 -9.31 3.0 142ms
54 -36.74248067268 + -Inf -9.20 3.0 139ms
55 -36.74248067268 -14.15 -9.39 3.0 145ms
56 -36.74248067268 + -13.85 -9.85 2.0 105ms
57 -36.74248067268 -14.15 -9.89 3.0 145ms
58 -36.74248067268 + -Inf -9.76 2.0 129ms
59 -36.74248067268 + -Inf -9.95 2.0 117ms
60 -36.74248067268 + -Inf -10.05 2.0 119ms
61 -36.74248067268 + -14.15 -10.57 1.0 98.6ms
62 -36.74248067268 + -Inf -10.57 3.0 143ms
63 -36.74248067268 -14.15 -10.75 2.0 116ms
64 -36.74248067268 + -13.85 -10.95 1.0 98.3ms
65 -36.74248067268 -13.85 -10.86 3.0 140ms
66 -36.74248067268 + -13.85 -11.15 2.0 113ms
67 -36.74248067268 -13.85 -11.02 3.0 136ms
68 -36.74248067268 + -14.15 -11.58 2.0 116ms
69 -36.74248067268 -14.15 -11.05 3.0 151ms
70 -36.74248067268 + -Inf -11.29 3.0 156ms
71 -36.74248067268 + -14.15 -11.51 2.0 127ms
72 -36.74248067268 + -Inf -12.14 2.0 118ms
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.73240620110 -0.88 11.0 377ms
2 -36.73980129374 -2.13 -1.36 1.0 93.5ms
3 -36.74070509400 -3.04 -1.83 6.0 132ms
4 -36.74198850438 -2.89 -2.04 2.0 121ms
5 -36.74241950692 -3.37 -2.61 1.0 97.6ms
6 -36.74246153737 -4.38 -2.71 2.0 129ms
7 -36.74247334098 -4.93 -2.90 1.0 97.4ms
8 -36.74247978412 -5.19 -3.56 1.0 98.8ms
9 -36.74247954799 + -6.63 -3.36 8.0 169ms
10 -36.74248043660 -6.05 -3.76 1.0 94.3ms
11 -36.74248064241 -6.69 -4.28 1.0 100ms
12 -36.74248066326 -7.68 -4.57 3.0 140ms
13 -36.74248067118 -8.10 -4.79 2.0 139ms
14 -36.74248066963 + -8.81 -5.14 1.0 96.2ms
15 -36.74248067261 -8.53 -5.54 2.0 119ms
16 -36.74248067267 -10.19 -5.80 3.0 142ms
17 -36.74248067267 + -12.64 -5.95 2.0 103ms
18 -36.74248067268 -10.95 -6.23 2.0 118ms
19 -36.74248067268 -11.90 -6.64 2.0 107ms
20 -36.74248067268 -12.59 -7.06 4.0 141ms
21 -36.74248067268 -13.30 -7.38 3.0 143ms
22 -36.74248067268 -13.55 -7.76 2.0 103ms
23 -36.74248067268 + -13.85 -7.87 3.0 144ms
24 -36.74248067268 -13.67 -8.45 2.0 107ms
25 -36.74248067268 + -13.85 -8.43 3.0 147ms
26 -36.74248067268 -14.15 -8.82 2.0 103ms
27 -36.74248067268 -14.15 -9.06 2.0 138ms
28 -36.74248067268 + -14.15 -9.58 1.0 99.5ms
29 -36.74248067268 -14.15 -9.66 3.0 146ms
30 -36.74248067268 + -Inf -9.93 1.0 96.3ms
31 -36.74248067268 + -Inf -10.24 2.0 116ms
32 -36.74248067268 + -Inf -10.57 2.0 136ms
33 -36.74248067268 + -13.85 -10.89 2.0 117ms
34 -36.74248067268 -13.85 -11.28 2.0 104ms
35 -36.74248067268 + -13.85 -11.39 3.0 146ms
36 -36.74248067268 -13.85 -11.55 2.0 106ms
37 -36.74248067268 + -Inf -11.89 1.0 99.2ms
38 -36.74248067268 + -Inf -11.87 3.0 133ms
39 -36.74248067268 + -Inf -12.14 1.0 101ms
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.02448898023324
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.244211113690376
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.723584343543123
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).