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.97952583063 -0.86 10.0 385ms
2 -35.93502095675 + -1.35 -1.51 1.0 91.0ms
3 +11.71496685359 + 1.68 -0.23 7.0 263ms
4 -35.83503144533 1.68 -1.12 8.0 289ms
5 -35.93028394113 -1.02 -1.50 3.0 146ms
6 -35.97905866242 -1.31 -1.82 2.0 112ms
7 -35.50217257411 + -0.32 -1.22 4.0 172ms
8 -35.97201711349 -0.33 -1.85 4.0 158ms
9 -35.98869552979 -1.78 -2.27 2.0 103ms
10 -35.98769083153 + -3.00 -2.08 2.0 143ms
11 -35.98944018643 -2.76 -2.55 2.0 171ms
12 -35.98971192890 -3.57 -2.76 1.0 97.5ms
13 -35.98975807497 -4.34 -3.12 1.0 99.1ms
14 -35.98976390919 -5.23 -3.24 2.0 139ms
15 -35.98927959496 + -3.31 -2.72 4.0 160ms
16 -35.98969950173 -3.38 -3.05 4.0 160ms
17 -35.98975263529 -4.27 -3.38 2.0 126ms
18 -35.98976053949 -5.10 -3.55 3.0 149ms
19 -35.98976384835 -5.48 -3.81 2.0 175ms
20 -35.98973677450 + -4.57 -3.32 4.0 183ms
21 -35.98976664443 -4.52 -4.28 3.0 152ms
22 -35.98976662238 + -7.66 -4.32 3.0 148ms
23 -35.98976676822 -6.84 -4.80 2.0 107ms
24 -35.98976677429 -8.22 -4.93 3.0 148ms
25 -35.98976676205 + -7.91 -4.74 3.0 133ms
26 -35.98976678317 -7.68 -5.28 2.0 150ms
27 -35.98976678301 + -9.80 -5.39 2.0 138ms
28 -35.98976678368 -9.17 -5.56 5.0 149ms
29 -35.98976678443 -9.13 -5.93 1.0 97.6ms
30 -35.98976678390 + -9.28 -5.54 3.0 153ms
31 -35.98976678427 -9.44 -5.82 4.0 159ms
32 -35.98976678448 -9.68 -6.15 5.0 152ms
33 -35.98976678418 + -9.54 -5.77 4.0 202ms
34 -35.98976678452 -9.47 -6.49 3.0 151ms
35 -35.98976678452 -11.92 -6.79 2.0 117ms
36 -35.98976678452 -11.74 -7.20 2.0 137ms
37 -35.98976678452 + -11.39 -6.73 4.0 160ms
38 -35.98976678452 -11.38 -7.42 4.0 160ms
39 -35.98976678452 -13.55 -7.71 2.0 144ms
40 -35.98976678452 + -12.81 -7.35 3.0 180ms
41 -35.98976678452 -12.81 -7.94 7.0 187ms
42 -35.98976678452 -13.85 -8.49 2.0 107ms
43 -35.98976678452 + -Inf -8.66 4.0 156ms
44 -35.98976678452 + -Inf -8.81 2.0 119ms
45 -35.98976678452 + -13.85 -8.82 2.0 107ms
46 -35.98976678452 -14.15 -9.12 2.0 109ms
47 -35.98976678452 + -14.15 -9.38 2.0 143ms
48 -35.98976678452 + -Inf -9.60 2.0 166ms
49 -35.98976678452 -14.15 -9.47 3.0 141ms
50 -35.98976678452 -14.15 -10.00 2.0 119ms
51 -35.98976678452 + -13.85 -10.11 3.0 150ms
52 -35.98976678452 -13.85 -10.38 1.0 98.1ms
53 -35.98976678452 + -Inf -10.23 4.0 162ms
54 -35.98976678452 + -13.85 -10.51 3.0 164ms
55 -35.98976678452 + -Inf -11.12 2.0 136ms
56 -35.98976678452 -13.85 -11.29 4.0 171ms
57 -35.98976678452 + -Inf -11.84 2.0 116ms
58 -35.98976678452 + -Inf -11.68 3.0 160ms
59 -35.98976678452 + -Inf -11.86 3.0 141ms
60 -35.98976678452 + -13.85 -11.74 2.0 129ms
61 -35.98976678452 -13.67 -12.48 2.0 147ms
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.98287709775 -0.86 10.0 429ms
2 -35.98808969606 -2.28 -1.34 1.0 92.7ms
3 -35.98776893400 + -3.49 -1.75 4.0 146ms
4 -35.98952299799 -2.76 -2.19 1.0 93.7ms
5 -35.98968615438 -3.79 -2.93 4.0 128ms
6 -35.98969750357 -4.95 -2.81 4.0 232ms
7 -35.98975939903 -4.21 -3.24 3.0 120ms
8 -35.98976030549 -6.04 -3.38 3.0 147ms
9 -35.98975946319 + -6.07 -3.52 2.0 118ms
10 -35.98976631600 -5.16 -3.93 1.0 97.0ms
11 -35.98976658601 -6.57 -4.29 2.0 137ms
12 -35.98976676381 -6.75 -4.55 5.0 140ms
13 -35.98976677909 -7.82 -4.88 2.0 117ms
14 -35.98976678436 -8.28 -5.75 8.0 222ms
15 -35.98976678452 -9.80 -5.81 3.0 172ms
16 -35.98976678445 + -10.17 -5.88 2.0 139ms
17 -35.98976678450 -10.28 -6.25 1.0 99.1ms
18 -35.98976678452 -10.65 -6.71 3.0 130ms
19 -35.98976678452 -13.45 -7.05 3.0 144ms
20 -35.98976678452 -12.60 -7.42 3.0 155ms
21 -35.98976678452 -13.37 -7.92 1.0 145ms
22 -35.98976678452 + -Inf -8.18 6.0 183ms
23 -35.98976678452 + -14.15 -8.54 9.0 194ms
24 -35.98976678452 -13.85 -8.97 3.0 151ms
25 -35.98976678452 + -13.85 -9.20 5.0 142ms
26 -35.98976678452 -13.85 -9.51 2.0 143ms
27 -35.98976678452 + -13.85 -10.12 1.0 100ms
28 -35.98976678452 -13.85 -10.55 4.0 214ms
29 -35.98976678452 + -13.85 -10.83 3.0 152ms
30 -35.98976678452 -13.85 -11.28 2.0 111ms
31 -35.98976678452 + -13.85 -11.51 4.0 164ms
32 -35.98976678452 + -Inf -12.01 7.0 163ms
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.38745571346606
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.62547869061561
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.6854599422886105
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).