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.73381335966                   -0.88   11.0    352ms
  2   -36.64462503839   +   -1.05       -1.43    1.0   96.0ms
  3   +31.37221362818   +    1.83       -0.15    7.0    222ms
  4   -36.58479363786        1.83       -1.10    7.0    237ms
  5   -33.75456415744   +    0.45       -0.83    4.0    160ms
  6   -36.60902259688        0.46       -1.43    4.0    157ms
  7   -36.73397290960       -0.90       -1.78    2.0    114ms
  8   -36.73752624236       -2.45       -2.07    2.0    114ms
  9   -36.74144229148       -2.41       -2.04    2.0    128ms
 10   -36.74240588272       -3.02       -2.39    1.0    100ms
 11   -36.74234773612   +   -4.24       -2.40    1.0   94.1ms
 12   -36.74234502399   +   -5.57       -2.77    1.0   95.4ms
 13   -36.74246660703       -3.92       -2.94    1.0    101ms
 14   -36.74242715645   +   -4.40       -2.98    1.0   95.5ms
 15   -36.74216349741   +   -3.58       -2.73    2.0    129ms
 16   -36.74237636404       -3.67       -2.93    2.0    124ms
 17   -36.74201723054   +   -3.44       -2.73    3.0    138ms
 18   -36.74247657827       -3.34       -3.48    2.0    123ms
 19   -36.74247330560   +   -5.49       -3.50    2.0    135ms
 20   -36.74247595297       -5.58       -3.47    2.0    124ms
 21   -36.74244531280   +   -4.51       -3.28    3.0    128ms
 22   -36.74247967914       -4.46       -3.94    2.0    121ms
 23   -36.74248065122       -6.01       -4.35    2.0    111ms
 24   -36.74248058258   +   -7.16       -4.48    3.0    138ms
 25   -36.74248063795       -7.26       -4.74    2.0    127ms
 26   -36.74248067007       -7.49       -5.20    2.0    131ms
 27   -36.74248067101       -9.03       -5.36    2.0    107ms
 28   -36.74248067067   +   -9.47       -5.26    2.0    130ms
 29   -36.74248060776   +   -7.20       -4.64    3.0    152ms
 30   -36.74248067259       -7.19       -5.94    3.0    147ms
 31   -36.74248067262      -10.52       -6.01    3.0    151ms
 32   -36.74248067263      -10.76       -5.99    2.0    122ms
 33   -36.74248067267      -10.41       -6.42    1.0   95.2ms
 34   -36.74248067264   +  -10.50       -6.21    3.0    143ms
 35   -36.74248067268      -10.40       -6.85    2.0    121ms
 36   -36.74248067268      -11.77       -7.17    2.0    135ms
 37   -36.74248067268   +  -12.94       -7.09    3.0    118ms
 38   -36.74248067268      -12.63       -7.41    2.0    126ms
 39   -36.74248067268      -12.85       -7.91    2.0    130ms
 40   -36.74248067268   +  -13.67       -7.71    3.0    134ms
 41   -36.74248067268      -14.15       -7.96    2.0    123ms
 42   -36.74248067268      -14.15       -8.22    2.0    111ms
 43   -36.74248067268      -14.15       -8.48    2.0    130ms
 44   -36.74248067268      -14.15       -8.61    3.0    116ms
 45   -36.74248067268   +  -14.15       -8.74    3.0    134ms
 46   -36.74248067268   +    -Inf       -9.15    1.0    100ms
 47   -36.74248067268   +    -Inf       -9.36    3.0    146ms
 48   -36.74248067268   +    -Inf       -9.42    2.0    116ms
 49   -36.74248067268   +    -Inf      -10.11    2.0    106ms
 50   -36.74248067268   +    -Inf       -9.86    3.0    149ms
 51   -36.74248067268   +  -14.15      -10.18    2.0    135ms
 52   -36.74248067268      -14.15      -10.22    2.0    111ms
 53   -36.74248067268   +    -Inf      -10.54    2.0    110ms
 54   -36.74248067268      -13.85      -10.42    2.0    132ms
 55   -36.74248067268   +  -13.85      -10.09    3.0    142ms
 56   -36.74248067268   +    -Inf      -10.99    3.0    138ms
 57   -36.74248067268      -14.15      -10.99    2.0    115ms
 58   -36.74248067268   +  -14.15      -11.46    2.0    112ms
 59   -36.74248067268   +    -Inf      -11.37    3.0    146ms
 60   -36.74248067268   +  -14.15      -11.53    3.0    138ms
 61   -36.74248067268      -14.15      -12.09    2.0    110ms

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.73296667142                   -0.88   12.0    351ms
  2   -36.73987986524       -2.16       -1.36    1.0   90.3ms
  3   -36.74036519327       -3.31       -1.71    4.0    129ms
  4   -36.74199710039       -2.79       -2.04    1.0   90.9ms
  5   -36.74240565455       -3.39       -2.58    4.0    107ms
  6   -36.74244219454       -4.44       -2.51    3.0    139ms
  7   -36.74247752057       -4.45       -3.06    1.0   93.4ms
  8   -36.74247962380       -5.68       -3.39    1.0   97.8ms
  9   -36.74247954376   +   -7.10       -3.36    5.0    126ms
 10   -36.74248061003       -5.97       -4.04    1.0   98.2ms
 11   -36.74248064742       -7.43       -4.18    4.0    138ms
 12   -36.74248066750       -7.70       -4.62    1.0    100ms
 13   -36.74248067221       -8.33       -5.06    3.0    126ms
 14   -36.74248067252       -9.50       -5.26    3.0    131ms
 15   -36.74248067260      -10.10       -5.46    2.0    120ms
 16   -36.74248067267      -10.20       -5.83    2.0    106ms
 17   -36.74248067268      -10.97       -6.07    3.0    135ms
 18   -36.74248067268      -11.44       -6.46    2.0    106ms
 19   -36.74248067268      -11.96       -7.15    3.0    127ms
 20   -36.74248067268      -13.45       -7.35    7.0    161ms
 21   -36.74248067268      -13.85       -7.58    2.0    107ms
 22   -36.74248067268   +    -Inf       -7.86    1.0    100ms
 23   -36.74248067268   +  -13.85       -8.42    2.0    111ms
 24   -36.74248067268      -13.67       -8.27    4.0    144ms
 25   -36.74248067268   +  -14.15       -8.53    1.0   97.2ms
 26   -36.74248067268      -13.85       -8.87    2.0    121ms
 27   -36.74248067268   +  -13.85       -9.29    1.0   95.3ms
 28   -36.74248067268      -14.15       -9.55    4.0    136ms
 29   -36.74248067268   +  -13.67       -9.84    1.0   95.8ms
 30   -36.74248067268      -14.15      -10.07    4.0    137ms
 31   -36.74248067268      -13.85      -10.68    1.0   95.3ms
 32   -36.74248067268   +  -13.67      -10.88    4.0    149ms
 33   -36.74248067268      -13.85      -11.21    2.0    101ms
 34   -36.74248067268   +  -14.15      -11.54    6.0    135ms
 35   -36.74248067268      -13.67      -11.65    2.0    128ms
 36   -36.74248067268   +  -14.15      -12.11    1.0   99.2ms

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.024488980277866

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.24421111373726

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))
Example block output

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))
Example block output

Clearly the charge-sloshing mode is no longer dominating.

The largest eigenvalue is now

maximum(real.(λ_Kerker))
4.723535012295406

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).