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.73283658161 -0.88 12.0 388ms
2 -36.65325604863 + -1.10 -1.42 1.0 92.0ms
3 +27.19429573113 + 1.81 -0.16 7.0 235ms
4 -36.17594804254 1.80 -0.99 7.0 249ms
5 -33.58164405801 + 0.41 -0.80 4.0 169ms
6 -36.44641860560 0.46 -1.27 4.0 169ms
7 -36.70558162824 -0.59 -1.63 3.0 131ms
8 -36.73720999520 -1.50 -1.97 2.0 113ms
9 -36.73557089660 + -2.79 -1.96 2.0 120ms
10 -36.73982771993 -2.37 -2.04 2.0 113ms
11 -36.74227153896 -2.61 -2.41 1.0 103ms
12 -36.74132341511 + -3.02 -2.34 3.0 147ms
13 -36.74232816713 -3.00 -2.71 1.0 105ms
14 -36.74231917484 + -5.05 -2.70 3.0 139ms
15 -36.74142958478 + -3.05 -2.52 3.0 127ms
16 -36.74244976764 -2.99 -3.09 2.0 132ms
17 -36.73616277783 + -2.20 -2.15 4.0 188ms
18 -36.74246464714 -2.20 -3.32 4.0 173ms
19 -36.74246428316 + -6.44 -3.36 2.0 116ms
20 -36.74247469617 -4.98 -3.62 3.0 130ms
21 -36.74247920436 -5.35 -3.77 2.0 111ms
22 -36.74248061300 -5.85 -4.33 2.0 139ms
23 -36.74247981408 + -6.10 -4.02 3.0 147ms
24 -36.74248061957 -6.09 -4.52 3.0 140ms
25 -36.74248043540 + -6.73 -4.17 3.0 153ms
26 -36.74248065880 -6.65 -4.78 2.0 138ms
27 -36.74248065291 + -8.23 -4.85 2.0 117ms
28 -36.74248066694 -7.85 -5.04 2.0 115ms
29 -36.74248034598 + -6.49 -4.29 4.0 161ms
30 -36.74248067170 -6.49 -5.36 4.0 168ms
31 -36.74248067109 + -9.22 -5.31 3.0 139ms
32 -36.74248067222 -8.94 -5.69 2.0 114ms
33 -36.74248067244 -9.68 -5.71 3.0 145ms
34 -36.74248067255 -9.93 -5.94 1.0 98.9ms
35 -36.74248067074 + -8.74 -5.41 4.0 167ms
36 -36.74248067267 -8.72 -6.34 3.0 144ms
37 -36.74248067267 -12.08 -6.35 2.0 138ms
38 -36.74248067268 -10.88 -6.81 1.0 100ms
39 -36.74248067268 + -12.01 -6.86 3.0 134ms
40 -36.74248067268 -11.84 -7.06 3.0 125ms
41 -36.74248067268 -12.44 -7.20 2.0 141ms
42 -36.74248067268 -12.47 -7.60 2.0 117ms
43 -36.74248067268 + -12.46 -7.16 3.0 152ms
44 -36.74248067268 -12.47 -7.56 3.0 148ms
45 -36.74248067268 -14.15 -7.47 2.0 136ms
46 -36.74248067268 -13.55 -8.09 2.0 116ms
47 -36.74248067268 + -14.15 -7.83 3.0 155ms
48 -36.74248067268 -13.85 -8.31 3.0 146ms
49 -36.74248067268 -13.85 -8.60 1.0 100ms
50 -36.74248067268 + -14.15 -8.76 2.0 143ms
51 -36.74248067268 + -Inf -8.87 2.0 130ms
52 -36.74248067268 -14.15 -8.74 3.0 132ms
53 -36.74248067268 + -14.15 -9.04 2.0 111ms
54 -36.74248067268 + -Inf -9.08 3.0 150ms
55 -36.74248067268 + -14.15 -9.50 2.0 117ms
56 -36.74248067268 -14.15 -9.74 2.0 141ms
57 -36.74248067268 -14.15 -9.55 3.0 159ms
58 -36.74248067268 + -14.15 -10.24 2.0 118ms
59 -36.74248067268 + -13.85 -9.75 3.0 154ms
60 -36.74248067268 -13.55 -9.97 3.0 149ms
61 -36.74248067268 + -13.85 -10.32 2.0 117ms
62 -36.74248067268 -14.15 -10.72 2.0 131ms
63 -36.74248067268 + -14.15 -11.06 2.0 136ms
64 -36.74248067268 + -14.15 -11.29 2.0 117ms
65 -36.74248067268 -13.67 -11.24 2.0 140ms
66 -36.74248067268 + -13.85 -11.34 2.0 128ms
67 -36.74248067268 -13.85 -11.15 3.0 140ms
68 -36.74248067268 + -13.85 -11.43 3.0 135ms
69 -36.74248067268 + -Inf -11.94 2.0 116ms
70 -36.74248067268 + -14.15 -12.04 3.0 148mswhile 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.73277580366 -0.88 12.0 367ms
2 -36.73997639692 -2.14 -1.36 1.0 98.2ms
3 -36.73973460955 + -3.62 -1.69 5.0 135ms
4 -36.74198578349 -2.65 -2.06 1.0 94.0ms
5 -36.74235417469 -3.43 -2.61 2.0 127ms
6 -36.74242512810 -4.15 -2.67 3.0 137ms
7 -36.74246434571 -4.41 -2.91 1.0 101ms
8 -36.74247689209 -4.90 -3.36 4.0 114ms
9 -36.74247965589 -5.56 -3.49 3.0 143ms
10 -36.74248052539 -6.06 -3.96 1.0 98.8ms
11 -36.74248063281 -6.97 -4.40 4.0 118ms
12 -36.74248060254 + -7.52 -4.54 3.0 141ms
13 -36.74248067167 -7.16 -5.13 2.0 109ms
14 -36.74248067253 -9.07 -5.57 3.0 149ms
15 -36.74248067258 -10.23 -5.64 3.0 123ms
16 -36.74248067265 -10.21 -5.98 2.0 110ms
17 -36.74248067266 -10.81 -6.30 3.0 140ms
18 -36.74248067268 -10.80 -6.49 1.0 104ms
19 -36.74248067268 + -12.10 -6.69 2.0 125ms
20 -36.74248067268 -11.31 -7.50 2.0 139ms
21 -36.74248067268 + -13.00 -7.43 5.0 171ms
22 -36.74248067268 -12.97 -8.08 2.0 115ms
23 -36.74248067268 -14.15 -8.42 2.0 139ms
24 -36.74248067268 + -14.15 -8.59 2.0 110ms
25 -36.74248067268 + -Inf -8.91 1.0 98.3ms
26 -36.74248067268 + -Inf -9.19 2.0 139ms
27 -36.74248067268 + -Inf -9.59 2.0 105ms
28 -36.74248067268 + -Inf -10.02 4.0 136ms
29 -36.74248067268 + -13.85 -10.00 3.0 145ms
30 -36.74248067268 -13.67 -10.58 2.0 120ms
31 -36.74248067268 + -14.15 -10.46 3.0 158ms
32 -36.74248067268 + -Inf -10.81 1.0 105ms
33 -36.74248067268 + -Inf -11.47 2.0 117ms
34 -36.74248067268 + -Inf -11.32 3.0 167ms
35 -36.74248067268 -14.15 -11.78 2.0 105ms
36 -36.74248067268 + -Inf -12.19 3.0 131msGiven 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
endepsilon (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.02448898025313The smallest eigenvalue is a bit more tricky to obtain, so we will just assume
λ_Simple_min = 0.9520.952This makes the condition number around 30:
cond_Simple = λ_Simple_max / λ_Simple_min46.24421111371127This 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.7235936329219035Since 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).