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.oncvpsp3.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.73428021615                   -0.88   11.0    398ms
  2   -36.71524842247   +   -1.72       -1.57    1.0   91.5ms
  3   -7.304051169376   +    1.47       -0.33    6.0    211ms
  4   -36.46068681068        1.46       -1.09    8.0    224ms
  5   -36.53554558489       -1.13       -1.32    4.0    153ms
  6   -36.19382764659   +   -0.47       -1.20    4.0    166ms
  7   -36.72215678971       -0.28       -1.83    3.0    163ms
  8   -36.73982008731       -1.75       -2.02    2.0    116ms
  9   -36.73963951969   +   -3.74       -2.16    2.0    135ms
 10   -36.74089717462       -2.90       -2.18    2.0    130ms
 11   -36.74183900087       -3.03       -2.37    1.0   96.3ms
 12   -36.74203236118       -3.71       -3.03    1.0    101ms
 13   -36.74207082424       -4.41       -3.21    3.0    146ms
 14   -36.74208208335       -4.95       -3.47    2.0    116ms
 15   -36.74203697229   +   -4.35       -3.22    3.0    158ms
 16   -36.74205913442       -4.65       -3.19    4.0    153ms
 17   -36.74191071209   +   -3.83       -2.92    3.0    149ms
 18   -36.74205674756       -3.84       -3.35    4.0    153ms
 19   -36.74207306843       -4.79       -3.51    3.0    134ms
 20   -36.74208361005       -4.98       -4.35    2.0    112ms
 21   -36.74208363290       -7.64       -4.33    4.0    175ms
 22   -36.74208365929       -7.58       -4.21    2.0    139ms
 23   -36.74208375297       -7.03       -4.71    2.0    130ms
 24   -36.74208375688       -8.41       -5.02    2.0    137ms
 25   -36.74208376444       -8.12       -5.28    1.0    101ms
 26   -36.74208376446      -10.57       -5.41    2.0    136ms
 27   -36.74208376434   +   -9.92       -5.65    1.0   96.2ms
 28   -36.74208376469       -9.45       -5.99    1.0    103ms
 29   -36.74208376173   +   -8.53       -5.31    4.0    172ms
 30   -36.74208376458       -8.55       -5.88    4.0    189ms
 31   -36.74208376478       -9.70       -6.43    3.0    122ms
 32   -36.74208376478   +  -11.56       -6.47    3.0    149ms
 33   -36.74208376478      -11.89       -6.52    5.0    139ms
 34   -36.74208376479      -11.08       -6.90    2.0    111ms
 35   -36.74208376479      -11.94       -7.15    2.0    103ms
 36   -36.74208376479   +  -14.15       -7.33    3.0    143ms
 37   -36.74208376479      -12.97       -7.42    1.0    103ms
 38   -36.74208376479      -13.85       -7.79    1.0   96.6ms
 39   -36.74208376479   +    -Inf       -8.06    2.0    135ms
 40   -36.74208376479   +  -13.85       -7.78    2.0    138ms
 41   -36.74208376479   +  -13.85       -7.73    4.0    149ms
 42   -36.74208376479   +  -13.85       -7.61    4.0    156ms
 43   -36.74208376479      -13.37       -8.46    3.0    147ms
 44   -36.74208376479      -13.85       -9.13    2.0    107ms
 45   -36.74208376479   +  -13.85       -8.81    4.0    174ms
 46   -36.74208376479   +    -Inf       -9.09    3.0    159ms
 47   -36.74208376479   +    -Inf       -9.78    2.0    125ms
 48   -36.74208376479      -13.85       -9.61    3.0    143ms
 49   -36.74208376479   +    -Inf       -9.72    2.0    121ms
 50   -36.74208376479   +  -13.85       -9.94    1.0   96.4ms
 51   -36.74208376479   +    -Inf      -10.39    1.0    101ms
 52   -36.74208376479      -14.15      -10.36    2.0    139ms
 53   -36.74208376479   +    -Inf      -10.88    1.0   96.2ms
 54   -36.74208376479   +  -14.15      -10.44    3.0    152ms
 55   -36.74208376479   +    -Inf      -10.32    4.0    221ms
 56   -36.74208376479   +    -Inf      -11.03    4.0    170ms
 57   -36.74208376479      -13.85      -11.21    2.0    116ms
 58   -36.74208376479   +  -13.85      -11.31    2.0    130ms
 59   -36.74208376479      -14.15      -11.24    3.0    147ms
 60   -36.74208376479   +    -Inf      -11.59    5.0    137ms
 61   -36.74208376479   +  -14.15      -12.25    2.0    135ms

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.73415368139                   -0.88   10.0    366ms
  2   -36.74003466815       -2.23       -1.36    1.0   97.0ms
  3   -36.73944685822   +   -3.23       -1.71    4.0    132ms
  4   -36.74182517334       -2.62       -2.23    1.0   98.1ms
  5   -36.74184538592       -4.69       -2.49    4.0    150ms
  6   -36.74203763190       -3.72       -2.54    3.0    108ms
  7   -36.74205089467       -4.88       -2.90    1.0    126ms
  8   -36.74207764017       -4.57       -3.09    2.0    109ms
  9   -36.74208169789       -5.39       -3.46    2.0    105ms
 10   -36.74208345541       -5.76       -4.09    6.0    137ms
 11   -36.74208352826       -7.14       -4.25    3.0    154ms
 12   -36.74208374879       -6.66       -4.55    2.0    121ms
 13   -36.74208373402   +   -7.83       -4.78    2.0    215ms
 14   -36.74208376446       -7.52       -5.61    2.0    117ms
 15   -36.74208376460       -9.85       -5.85    6.0    837ms
 16   -36.74208376473       -9.89       -6.11    3.0    123ms
 17   -36.74208376477      -10.42       -6.38    2.0    137ms
 18   -36.74208376479      -10.71       -6.95    2.0    118ms
 19   -36.74208376479      -13.00       -7.19    4.0    150ms
 20   -36.74208376479      -12.94       -7.50    2.0    108ms
 21   -36.74208376479      -13.45       -7.90    6.0    144ms
 22   -36.74208376479   +  -13.85       -8.05    2.0    151ms
 23   -36.74208376479   +    -Inf       -8.47    1.0    109ms
 24   -36.74208376479   +    -Inf       -8.79    5.0    128ms
 25   -36.74208376479      -13.85       -9.25    2.0    117ms
 26   -36.74208376479   +  -13.85       -9.53    4.0    159ms
 27   -36.74208376479      -13.85       -9.86    1.0    105ms
 28   -36.74208376479   +  -13.67      -10.07    3.0    145ms
 29   -36.74208376479      -13.85      -10.56    2.0    106ms
 30   -36.74208376479   +    -Inf      -10.86    3.0    144ms
 31   -36.74208376479   +    -Inf      -11.15    9.0    157ms
 32   -36.74208376479      -14.15      -11.81    2.0    107ms
 33   -36.74208376479   +    -Inf      -11.92    4.0    156ms
 34   -36.74208376479   +  -14.15      -12.15    2.0    130ms

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

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

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

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