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 ASEconvert
using DFTK
using LazyArtifacts
ase_Al = ase.build.bulk("Al"; cubic=true) * pytuple((4, 1, 1))
system_Al = attach_psp(pyconvert(AbstractSystem, ase_Al);
Al=artifact"pd_nc_sr_pbe_standard_0.4.1_upf/Al.upf")
FlexibleSystem(Al₁₆, periodic = TTT):
bounding_box : [ 16.2 0 0;
0 4.05 0;
0 0 4.05]u"Å"
.---------------------------------------.
/| |
* | |
|Al Al Al Al |
| | |
| .--Al--------Al--------Al--------Al-----.
|/ Al Al Al Al /
Al--------Al--------Al--------Al--------*
and we discretise:
model_Al = model_LDA(system_Al; temperature=1e-3, symmetries=false)
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 -35.97981658525 -0.86 13.0 366ms
2 -35.93804429713 + -1.38 -1.53 1.0 88.4ms
3 +11.09753977603 + 1.67 -0.23 7.0 212ms
4 -35.92948387565 1.67 -1.24 9.0 246ms
5 -35.94449046271 -1.82 -1.61 2.0 118ms
6 -35.95420237704 -2.01 -1.69 2.0 106ms
7 -35.44722887040 + -0.30 -1.20 4.0 146ms
8 -35.98457164922 -0.27 -2.07 3.0 151ms
9 -35.98856221782 -2.40 -2.12 2.0 131ms
10 -35.98837160542 + -3.72 -2.20 2.0 119ms
11 -35.98947926686 -2.96 -2.67 2.0 108ms
12 -35.98961195725 -3.88 -2.94 1.0 92.1ms
13 -35.98962393316 -4.92 -3.30 2.0 123ms
14 -35.98962058209 + -5.47 -3.14 3.0 136ms
15 -35.98956738266 + -4.27 -3.16 3.0 126ms
16 -35.98961969383 -4.28 -3.49 3.0 138ms
17 -35.98963094148 -4.95 -3.81 2.0 109ms
18 -35.98963185038 -6.04 -4.37 3.0 116ms
19 -35.98963083347 + -5.99 -3.99 3.0 146ms
20 -35.98963181984 -6.01 -4.39 4.0 144ms
21 -35.98963187588 -7.25 -4.39 2.0 108ms
22 -35.98963197661 -7.00 -4.74 2.0 100ms
23 -35.98963198038 -8.42 -5.05 3.0 134ms
24 -35.98963198395 -8.45 -5.25 2.0 113ms
25 -35.98963198433 -9.42 -5.37 1.0 106ms
26 -35.98963198427 + -10.24 -5.18 3.0 129ms
27 -35.98963198607 -8.75 -6.05 2.0 98.4ms
28 -35.98963198608 -10.84 -6.02 3.0 141ms
29 -35.98963198604 + -10.35 -5.91 2.0 112ms
30 -35.98963198605 -10.83 -6.10 3.0 126ms
31 -35.98963198611 -10.21 -6.27 2.0 110ms
32 -35.98963198613 -10.75 -6.89 2.0 100ms
33 -35.98963198612 + -11.11 -6.58 3.0 143ms
34 -35.98963198613 -11.15 -6.94 3.0 136ms
35 -35.98963198613 -12.50 -7.01 2.0 120ms
36 -35.98963198613 -12.19 -7.16 2.0 100ms
37 -35.98963198613 -12.81 -7.50 2.0 100ms
38 -35.98963198613 + -13.85 -7.54 2.0 127ms
39 -35.98963198613 -13.67 -7.70 3.0 116ms
40 -35.98963198613 -13.55 -7.83 1.0 96.8ms
41 -35.98963198613 -14.15 -8.57 1.0 93.2ms
42 -35.98963198613 -14.15 -8.70 3.0 140ms
43 -35.98963198613 + -14.15 -8.66 5.0 144ms
44 -35.98963198613 + -Inf -9.04 2.0 112ms
45 -35.98963198613 -14.15 -9.01 3.0 137ms
46 -35.98963198613 + -Inf -8.93 2.0 120ms
47 -35.98963198613 + -Inf -9.42 2.0 109ms
48 -35.98963198613 + -Inf -9.77 2.0 110ms
49 -35.98963198613 + -Inf -9.72 3.0 134ms
50 -35.98963198613 + -Inf -9.76 2.0 113ms
51 -35.98963198613 + -Inf -10.22 2.0 110ms
52 -35.98963198613 + -14.15 -10.45 2.0 117ms
53 -35.98963198613 + -Inf -10.66 2.0 130ms
54 -35.98963198613 -14.15 -10.83 5.0 129ms
55 -35.98963198613 + -Inf -10.85 1.0 93.9ms
56 -35.98963198613 + -Inf -11.24 2.0 105ms
57 -35.98963198613 + -Inf -11.47 3.0 136ms
58 -35.98963198613 + -14.15 -11.49 3.0 117ms
59 -35.98963198613 -14.15 -11.79 2.0 110ms
60 -35.98963198613 + -Inf -12.08 1.0 92.6ms
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 -35.98238251815 -0.86 11.0 351ms
2 -35.98783086439 -2.26 -1.34 1.0 88.8ms
3 -35.98738844474 + -3.35 -1.70 6.0 130ms
4 -35.98939347500 -2.70 -2.25 1.0 88.6ms
5 -35.98954556104 -3.82 -2.49 6.0 156ms
6 -35.98960904463 -4.20 -2.65 2.0 98.1ms
7 -35.98961279580 -5.43 -2.94 2.0 97.8ms
8 -35.98962415713 -4.94 -3.10 2.0 101ms
9 -35.98962969502 -5.26 -3.50 8.0 136ms
10 -35.98963171931 -5.69 -4.18 6.0 132ms
11 -35.98963176905 -7.30 -4.26 4.0 150ms
12 -35.98963197204 -6.69 -4.55 2.0 110ms
13 -35.98963196033 + -7.93 -4.80 2.0 127ms
14 -35.98963198575 -7.59 -5.51 2.0 101ms
15 -35.98963198602 -9.56 -5.79 4.0 160ms
16 -35.98963198612 -10.02 -6.15 2.0 118ms
17 -35.98963198611 + -11.26 -6.38 3.0 137ms
18 -35.98963198613 -10.83 -6.64 5.0 129ms
19 -35.98963198613 -11.62 -7.17 2.0 101ms
20 -35.98963198613 + -12.66 -7.22 3.0 141ms
21 -35.98963198613 -12.57 -7.85 2.0 100ms
22 -35.98963198613 -13.67 -7.92 3.0 136ms
23 -35.98963198613 + -14.15 -8.45 1.0 93.2ms
24 -35.98963198613 -14.15 -8.74 5.0 126ms
25 -35.98963198613 + -13.85 -8.89 3.0 142ms
26 -35.98963198613 -13.85 -9.25 1.0 93.6ms
27 -35.98963198613 + -14.15 -9.58 3.0 116ms
28 -35.98963198613 + -13.85 -9.89 3.0 135ms
29 -35.98963198613 -13.85 -10.31 2.0 103ms
30 -35.98963198613 + -13.85 -10.70 3.0 154ms
31 -35.98963198613 -13.85 -10.89 10.0 206ms
32 -35.98963198613 + -Inf -11.17 5.0 141ms
33 -35.98963198613 + -Inf -11.50 2.0 181ms
34 -35.98963198613 + -Inf -11.86 2.0 106ms
35 -35.98963198613 + -Inf -12.42 2.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.387452180856535
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.62547497989132
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.68545844693817
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).