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.oncvpsp3.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.73531893402 -0.88 11.0 402ms
2 -36.72334171715 + -1.92 -1.63 1.0 96.3ms
3 -12.80539502420 + 1.38 -0.38 8.0 222ms
4 -36.69651556233 1.38 -1.34 8.0 209ms
5 -36.70827751734 -1.93 -1.71 2.0 128ms
6 -36.33102751100 + -0.42 -1.26 4.0 144ms
7 -36.72397256281 -0.41 -1.81 3.0 151ms
8 -36.74134456692 -1.76 -2.28 2.0 131ms
9 -36.74160915865 -3.58 -2.25 3.0 143ms
10 -36.74133007527 + -3.55 -2.28 2.0 114ms
11 -36.74200311688 -3.17 -2.72 1.0 100ms
12 -36.74205126599 -4.32 -2.90 1.0 102ms
13 -36.74207317789 -4.66 -3.09 2.0 138ms
14 -36.74208106177 -5.10 -3.69 1.0 97.5ms
15 -36.74169040375 + -3.41 -2.77 5.0 181ms
16 -36.74208289959 -3.41 -3.93 4.0 307ms
17 -36.74207838538 + -5.35 -3.41 3.0 759ms
18 -36.74200703653 + -4.15 -3.13 4.0 200ms
19 -36.74208308969 -4.12 -4.09 4.0 163ms
20 -36.74208371523 -6.20 -4.35 2.0 121ms
21 -36.74208376195 -7.33 -4.81 2.0 150ms
22 -36.74208374414 + -7.75 -4.65 2.0 116ms
23 -36.74208376211 -7.75 -4.90 1.0 98.0ms
24 -36.74208376449 -8.62 -5.44 2.0 109ms
25 -36.74208376406 + -9.36 -5.37 3.0 142ms
26 -36.74208376463 -9.24 -5.73 1.0 98.9ms
27 -36.74208363055 + -6.87 -4.50 6.0 202ms
28 -36.74208376470 -6.87 -6.00 5.0 207ms
29 -36.74208376478 -10.10 -6.20 2.0 137ms
30 -36.74208376464 + -9.88 -5.96 4.0 154ms
31 -36.74208376371 + -9.03 -5.54 4.0 164ms
32 -36.74208376479 -8.97 -6.41 4.0 169ms
33 -36.74208376479 -11.58 -6.77 2.0 103ms
34 -36.74208376479 -12.07 -6.93 2.0 137ms
35 -36.74208376479 -12.62 -7.48 1.0 101ms
36 -36.74208376479 -13.25 -7.77 2.0 138ms
37 -36.74208376479 + -13.85 -7.58 3.0 128ms
38 -36.74208376479 -13.85 -8.12 1.0 101ms
39 -36.74208376479 + -14.15 -8.10 3.0 145ms
40 -36.74208376479 + -13.67 -7.77 3.0 175ms
41 -36.74208376479 -13.45 -8.66 4.0 140ms
42 -36.74208376479 + -14.15 -8.63 3.0 158ms
43 -36.74208376479 + -13.67 -7.90 4.0 164ms
44 -36.74208376479 -13.55 -8.85 4.0 164ms
45 -36.74208376479 + -13.85 -9.33 2.0 110ms
46 -36.74208376479 -14.15 -9.16 3.0 146ms
47 -36.74208376479 + -14.15 -9.91 2.0 117ms
48 -36.74208376479 -14.15 -10.03 3.0 176ms
49 -36.74208376479 -14.15 -10.15 2.0 120ms
50 -36.74208376479 + -13.85 -10.58 1.0 97.1ms
51 -36.74208376479 -13.85 -10.81 2.0 138ms
52 -36.74208376479 + -14.15 -11.00 2.0 109ms
53 -36.74208376479 + -14.15 -11.14 2.0 119ms
54 -36.74208376479 -13.85 -10.68 4.0 165ms
55 -36.74208376479 + -13.85 -11.27 3.0 147ms
56 -36.74208376479 + -Inf -11.85 3.0 128ms
57 -36.74208376479 -14.15 -11.98 3.0 146ms
58 -36.74208376479 + -14.15 -12.20 2.0 104ms
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.73406806115 -0.88 11.0 371ms
2 -36.73992122803 -2.23 -1.36 1.0 96.8ms
3 -36.73807959166 + -2.73 -1.62 3.0 127ms
4 -36.74186253640 -2.42 -2.35 2.0 119ms
5 -36.74189102830 -4.55 -2.53 4.0 147ms
6 -36.74202521916 -3.87 -2.75 7.0 140ms
7 -36.74206861779 -4.36 -2.77 3.0 114ms
8 -36.74206915422 -6.27 -3.25 2.0 109ms
9 -36.74208326893 -4.85 -3.53 3.0 145ms
10 -36.74208328272 -7.86 -4.00 2.0 110ms
11 -36.74208345679 -6.76 -4.28 5.0 137ms
12 -36.74208374089 -6.55 -4.56 3.0 169ms
13 -36.74208376431 -7.63 -5.09 2.0 105ms
14 -36.74208376466 -9.44 -5.40 8.0 184ms
15 -36.74208376475 -10.08 -5.73 2.0 107ms
16 -36.74208376477 -10.64 -6.17 4.0 127ms
17 -36.74208376479 -10.76 -6.65 3.0 147ms
18 -36.74208376479 -13.15 -6.83 3.0 147ms
19 -36.74208376479 -12.67 -6.95 5.0 142ms
20 -36.74208376479 -11.96 -7.47 2.0 111ms
21 -36.74208376479 -13.30 -7.88 4.0 155ms
22 -36.74208376479 -13.85 -8.18 7.0 144ms
23 -36.74208376479 + -Inf -8.29 3.0 145ms
24 -36.74208376479 + -Inf -8.63 2.0 109ms
25 -36.74208376479 + -Inf -8.85 3.0 127ms
26 -36.74208376479 + -Inf -9.47 1.0 102ms
27 -36.74208376479 + -13.85 -9.67 8.0 199ms
28 -36.74208376479 + -Inf -10.04 2.0 107ms
29 -36.74208376479 -14.15 -10.18 3.0 148ms
30 -36.74208376479 -14.15 -10.62 1.0 103ms
31 -36.74208376479 + -14.15 -10.83 5.0 135ms
32 -36.74208376479 + -Inf -11.23 3.0 146ms
33 -36.74208376479 + -Inf -11.55 5.0 130ms
34 -36.74208376479 -14.15 -11.88 3.0 138ms
35 -36.74208376479 + -13.85 -12.15 2.0 142ms
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.024488546359336
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.24421065794048
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.723592146727547
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).