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
using LazyArtifacts
al_supercell = bulk(:Al; cubic=true) * (4, 1, 1)
system_Al = attach_psp(al_supercell;
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_DFT(system_Al; functionals=LDA(), 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.98161598310 -0.86 11.0 484ms
2 -35.89310345687 + -1.05 -1.46 1.0 102ms
3 +34.19023978773 + 1.85 -0.14 19.0 380ms
4 -35.73036916574 1.84 -1.21 7.0 315ms
5 -34.54607486613 + 0.07 -0.99 4.0 175ms
6 -35.81517831575 0.10 -1.39 4.0 180ms
7 -35.97849675751 -0.79 -1.89 3.0 144ms
8 -35.98763692054 -2.04 -2.11 2.0 146ms
9 -35.98676336305 + -3.06 -1.99 2.0 140ms
10 -35.98908705665 -2.63 -2.34 2.0 138ms
11 -35.98775412555 + -2.88 -2.29 2.0 170ms
12 -35.98935352324 -2.80 -2.65 1.0 109ms
13 -35.98962317921 -3.57 -3.03 2.0 131ms
14 -35.98959449040 + -4.54 -3.00 2.0 153ms
15 -35.98911337319 + -3.32 -2.66 3.0 146ms
16 -35.98960052823 -3.31 -3.13 5.0 155ms
17 -35.98924049716 + -3.44 -2.75 3.0 166ms
18 -35.98962575574 -3.41 -3.35 3.0 192ms
19 -35.98951805842 + -3.97 -3.01 3.0 165ms
20 -35.98962789922 -3.96 -3.70 3.0 149ms
21 -35.98963125118 -5.47 -4.09 2.0 148ms
22 -35.98963159752 -6.46 -4.23 3.0 159ms
23 -35.98963193907 -6.47 -4.43 2.0 150ms
24 -35.98963181270 + -6.90 -4.32 3.0 138ms
25 -35.98963197982 -6.78 -5.00 1.0 141ms
26 -35.98963198310 -8.48 -5.21 3.0 171ms
27 -35.98963198503 -8.71 -5.17 2.0 129ms
28 -35.98963198600 -9.01 -5.88 1.0 107ms
29 -35.98963198549 + -9.29 -5.64 3.0 166ms
30 -35.98963198356 + -8.71 -5.32 4.0 197ms
31 -35.98963198604 -8.61 -5.99 3.0 166ms
32 -35.98963198557 + -9.33 -5.65 3.0 189ms
33 -35.98963198612 -9.26 -6.56 4.0 158ms
34 -35.98963198613 -11.10 -7.14 3.0 153ms
35 -35.98963198613 -13.00 -7.35 4.0 170ms
36 -35.98963198613 + -12.83 -7.32 3.0 139ms
37 -35.98963198613 -14.15 -7.38 3.0 151ms
38 -35.98963198613 -12.83 -7.73 3.0 173ms
39 -35.98963198613 -13.85 -7.79 2.0 116ms
40 -35.98963198613 -13.85 -8.05 2.0 118ms
41 -35.98963198613 + -13.85 -8.18 2.0 148ms
42 -35.98963198613 + -13.85 -7.83 3.0 149ms
43 -35.98963198613 -13.67 -8.95 5.0 147ms
44 -35.98963198613 -14.15 -8.67 4.0 197ms
45 -35.98963198613 + -13.85 -8.91 3.0 198ms
46 -35.98963198613 + -Inf -9.15 3.0 155ms
47 -35.98963198613 -13.85 -9.39 1.0 109ms
48 -35.98963198613 + -Inf -9.39 2.0 149ms
49 -35.98963198613 + -14.15 -9.96 2.0 111ms
50 -35.98963198613 + -13.85 -9.74 3.0 165ms
51 -35.98963198613 -14.15 -10.09 3.0 147ms
52 -35.98963198613 -14.15 -10.14 2.0 128ms
53 -35.98963198613 -14.15 -10.43 2.0 137ms
54 -35.98963198613 + -13.85 -10.65 2.0 149ms
55 -35.98963198613 -13.85 -10.79 1.0 106ms
56 -35.98963198613 + -13.85 -10.75 2.0 142ms
57 -35.98963198613 -13.85 -10.84 2.0 120ms
58 -35.98963198613 + -13.85 -10.58 3.0 163ms
59 -35.98963198613 + -Inf -11.02 4.0 154ms
60 -35.98963198613 + -Inf -11.04 2.0 153ms
61 -35.98963198613 + -Inf -11.77 2.0 125ms
62 -35.98963198613 + -Inf -11.61 4.0 183ms
63 -35.98963198613 -14.15 -11.61 1.0 108ms
64 -35.98963198613 + -14.15 -12.29 2.0 131ms
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.98246302732 -0.86 13.0 479ms
2 -35.98785214251 -2.27 -1.34 1.0 103ms
3 -35.98810276555 -3.60 -1.76 3.0 145ms
4 -35.98944647371 -2.87 -2.22 1.0 104ms
5 -35.98957058679 -3.91 -2.55 10.0 201ms
6 -35.98959119792 -4.69 -2.53 2.0 117ms
7 -35.98962207643 -4.51 -3.05 2.0 119ms
8 -35.98962867079 -5.18 -3.19 5.0 208ms
9 -35.98963057367 -5.72 -3.48 1.0 103ms
10 -35.98963184536 -5.90 -4.27 1.0 109ms
11 -35.98963185496 -8.02 -4.32 4.0 250ms
12 -35.98963198075 -6.90 -4.77 2.0 710ms
13 -35.98963196586 + -7.83 -4.81 4.0 151ms
14 -35.98963198603 -7.70 -5.64 3.0 125ms
15 -35.98963198601 + -10.69 -5.84 8.0 180ms
16 -35.98963198609 -10.11 -6.15 2.0 108ms
17 -35.98963198612 -10.45 -6.51 3.0 151ms
18 -35.98963198613 -11.45 -6.79 2.0 138ms
19 -35.98963198613 -11.66 -7.15 6.0 187ms
20 -35.98963198613 -12.79 -7.40 3.0 149ms
21 -35.98963198613 -12.94 -7.92 2.0 107ms
22 -35.98963198613 + -14.15 -8.11 3.0 144ms
23 -35.98963198613 + -Inf -8.59 3.0 130ms
24 -35.98963198613 + -Inf -8.94 5.0 157ms
25 -35.98963198613 + -14.15 -9.12 6.0 134ms
26 -35.98963198613 + -Inf -9.46 2.0 124ms
27 -35.98963198613 + -Inf -9.90 3.0 153ms
28 -35.98963198613 -14.15 -10.13 5.0 137ms
29 -35.98963198613 + -Inf -10.45 2.0 111ms
30 -35.98963198613 + -Inf -10.79 4.0 151ms
31 -35.98963198613 + -13.85 -11.33 5.0 142ms
32 -35.98963198613 -14.15 -11.45 8.0 184ms
33 -35.98963198613 + -Inf -11.82 1.0 104ms
34 -35.98963198613 -14.15 -12.03 3.0 154ms
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.387461708332786
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.62548498774453
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.685462194152469
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).