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.73221841801 -0.88 11.0 432ms
2 -36.58873535212 + -0.84 -1.37 1.0 91.9ms
┌ Warning: Eigensolver not converged
│ n_iter =
│ 1-element Vector{Int64}:
│ 22
└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/scf/self_consistent_field.jl:76
3 +50.69650308874 + 1.94 -0.09 22.0 354ms
4 -36.56037747978 1.94 -1.10 9.0 349ms
5 -33.65826785918 + 0.46 -0.81 4.0 229ms
6 -36.42291139455 0.44 -1.28 4.0 180ms
7 -36.71655214650 -0.53 -1.65 3.0 148ms
8 -36.73019809285 -1.86 -1.90 2.0 239ms
9 -36.72547908948 + -2.33 -1.77 2.0 159ms
10 -36.74114415150 -1.81 -2.13 2.0 734ms
11 -36.74156660454 -3.37 -2.31 2.0 125ms
12 -36.74232326261 -3.12 -2.51 1.0 102ms
13 -36.74239550598 -4.14 -2.65 1.0 104ms
14 -36.74244446412 -4.31 -2.72 1.0 107ms
15 -36.73815771764 + -2.37 -2.24 3.0 192ms
16 -36.74224292843 -2.39 -2.77 3.0 153ms
17 -36.72241603436 + -1.70 -1.89 4.0 193ms
18 -36.74248212892 -1.70 -3.24 4.0 188ms
19 -36.74213757453 + -3.46 -2.72 3.0 168ms
20 -36.74230451996 -3.78 -2.80 2.0 134ms
21 -36.74251101340 -3.69 -3.57 3.0 203ms
22 -36.74251403239 -5.52 -4.01 2.0 126ms
23 -36.74251465944 -6.20 -4.13 2.0 141ms
24 -36.74251473337 -7.13 -4.40 2.0 129ms
25 -36.74251462657 + -6.97 -4.22 2.0 147ms
26 -36.74251477091 -6.84 -5.01 2.0 128ms
27 -36.74251476238 + -8.07 -5.01 3.0 157ms
28 -36.74251475122 + -7.95 -4.87 3.0 149ms
29 -36.74251477199 -7.68 -5.45 2.0 161ms
30 -36.74251477092 + -8.97 -5.37 3.0 153ms
31 -36.74251476605 + -8.31 -5.10 3.0 157ms
32 -36.74251477177 -8.24 -5.46 3.0 157ms
33 -36.74251477303 -8.90 -6.48 2.0 127ms
34 -36.74251477303 -11.60 -6.53 3.0 188ms
35 -36.74251477303 + -11.43 -6.39 2.0 121ms
36 -36.74251477304 -11.28 -6.78 2.0 152ms
37 -36.74251477304 -13.55 -6.81 2.0 127ms
38 -36.74251477304 -11.77 -7.15 2.0 110ms
39 -36.74251477304 -12.39 -7.73 2.0 126ms
40 -36.74251477304 + -12.64 -7.28 3.0 166ms
41 -36.74251477304 -12.72 -7.57 3.0 161ms
42 -36.74251477304 + -12.42 -7.23 3.0 153ms
43 -36.74251477304 -12.51 -7.51 3.0 189ms
44 -36.74251477304 -13.03 -7.77 3.0 157ms
45 -36.74251477304 -13.67 -8.49 1.0 106ms
46 -36.74251477304 + -14.15 -8.48 2.0 141ms
47 -36.74251477304 + -Inf -8.16 3.0 156ms
48 -36.74251477304 + -13.85 -8.94 2.0 133ms
49 -36.74251477304 + -Inf -8.57 2.0 145ms
50 -36.74251477304 + -Inf -8.97 3.0 195ms
51 -36.74251477304 -13.85 -8.91 2.0 128ms
52 -36.74251477304 + -Inf -9.53 1.0 105ms
53 -36.74251477304 + -14.15 -9.81 3.0 198ms
54 -36.74251477304 -14.15 -10.26 1.0 101ms
55 -36.74251477304 + -14.15 -9.47 4.0 189ms
56 -36.74251477304 + -14.15 -10.29 3.0 168ms
57 -36.74251477304 -14.15 -10.46 2.0 156ms
58 -36.74251477304 -14.15 -10.10 2.0 135ms
59 -36.74251477304 -14.15 -10.43 2.0 160ms
60 -36.74251477304 + -13.85 -10.34 3.0 158ms
61 -36.74251477304 -14.15 -11.04 2.0 121ms
62 -36.74251477304 + -Inf -10.68 3.0 156ms
63 -36.74251477304 + -Inf -10.90 2.0 140ms
64 -36.74251477304 + -Inf -11.02 3.0 206ms
65 -36.74251477304 + -Inf -11.52 1.0 101ms
66 -36.74251477304 + -Inf -11.61 2.0 146ms
67 -36.74251477304 + -13.85 -11.14 3.0 158ms
68 -36.74251477304 -14.15 -11.58 3.0 169ms
69 -36.74251477304 -14.15 -11.62 2.0 133ms
70 -36.74251477304 + -Inf -11.85 1.0 105ms
71 -36.74251477304 + -Inf -12.07 2.0 127ms
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.73462738510 -0.88 11.0 398ms
2 -36.74034138949 -2.24 -1.36 1.0 98.8ms
3 -36.74145732551 -2.95 -2.10 3.0 115ms
4 -36.74217291917 -3.15 -2.07 4.0 150ms
5 -36.74242204792 -3.60 -2.58 1.0 99.2ms
6 -36.74247175529 -4.30 -2.61 5.0 162ms
7 -36.74248760447 -4.80 -2.94 1.0 100ms
8 -36.74251224660 -4.61 -3.42 2.0 107ms
9 -36.74251388414 -5.79 -3.52 4.0 149ms
10 -36.74251466267 -6.11 -3.77 1.0 102ms
11 -36.74251463363 + -7.54 -3.85 2.0 123ms
12 -36.74251476392 -6.89 -4.57 1.0 104ms
13 -36.74251477232 -8.08 -4.81 3.0 162ms
14 -36.74251476979 + -8.60 -4.74 1.0 103ms
15 -36.74251477271 -8.53 -5.54 1.0 142ms
16 -36.74251477300 -9.54 -5.69 4.0 155ms
17 -36.74251477303 -10.57 -5.92 1.0 100ms
18 -36.74251477304 -11.07 -6.37 2.0 118ms
19 -36.74251477304 -11.70 -6.80 6.0 140ms
20 -36.74251477304 -12.62 -7.09 2.0 141ms
21 -36.74251477304 + -13.85 -7.27 2.0 106ms
22 -36.74251477304 + -13.85 -7.37 2.0 110ms
23 -36.74251477304 -13.19 -7.74 2.0 159ms
24 -36.74251477304 + -14.15 -8.13 2.0 115ms
25 -36.74251477304 + -13.67 -8.55 3.0 131ms
26 -36.74251477304 -13.85 -8.74 3.0 149ms
27 -36.74251477304 + -14.15 -9.17 2.0 123ms
28 -36.74251477304 + -14.15 -9.47 3.0 148ms
29 -36.74251477304 -14.15 -9.96 2.0 106ms
30 -36.74251477304 -14.15 -10.18 3.0 167ms
31 -36.74251477304 + -13.85 -10.67 3.0 167ms
32 -36.74251477304 -13.85 -10.89 3.0 168ms
33 -36.74251477304 + -14.15 -11.59 2.0 106ms
34 -36.74251477304 + -Inf -11.52 4.0 161ms
35 -36.74251477304 + -14.15 -11.84 1.0 103ms
36 -36.74251477304 -14.15 -12.08 3.0 123ms
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.02448824375489
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.244210340078666
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.723581069536744
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).