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.73323730013                   -0.88   12.0    400ms
  2   -36.65612793990   +   -1.11       -1.45    1.0   90.0ms
  3   +28.98931407855   +    1.82       -0.15    7.0    295ms
  4   -36.45545231838        1.82       -1.06    6.0    243ms
  5   -36.61025039698       -0.81       -1.35    3.0    146ms
  6   -35.56157323675   +    0.02       -1.03    4.0    160ms
  7   -36.69864963710        0.06       -1.65    3.0    139ms
  8   -36.73670297328       -1.42       -2.01    2.0    118ms
  9   -36.73901169936       -2.64       -2.06    2.0    178ms
 10   -36.73575513762   +   -2.49       -1.97    2.0    121ms
 11   -36.74152569886       -2.24       -2.23    2.0    126ms
 12   -36.74215551913       -3.20       -2.53    1.0    101ms
 13   -36.74242617505       -3.57       -2.75    1.0   96.3ms
 14   -36.74240691278   +   -4.72       -2.77    2.0    111ms
 15   -36.74106349352   +   -2.87       -2.47    3.0    143ms
 16   -36.74174886310       -3.16       -2.54    3.0    143ms
 17   -36.74244445047       -3.16       -3.11    2.0    120ms
 18   -36.74218192581   +   -3.58       -2.79    3.0    189ms
 19   -36.74246834522       -3.54       -3.22    3.0    142ms
 20   -36.74250812968       -4.40       -3.48    2.0    114ms
 21   -36.74251270876       -5.34       -3.81    1.0    101ms
 22   -36.74251404421       -5.87       -3.97    2.0    138ms
 23   -36.74251472013       -6.17       -4.28    2.0    115ms
 24   -36.74251456482   +   -6.81       -4.26    3.0    134ms
 25   -36.74251471140       -6.83       -4.52    1.0    123ms
 26   -36.74251475521       -7.36       -4.80    3.0    155ms
 27   -36.74251471523   +   -7.40       -4.46    3.0    151ms
 28   -36.74251476670       -7.29       -5.11    2.0    130ms
 29   -36.74251477211       -8.27       -5.42    2.0    108ms
 30   -36.74251473601   +   -7.44       -4.76    4.0    170ms
 31   -36.74251467534   +   -7.22       -4.56    4.0    177ms
 32   -36.74251477268       -7.01       -5.59    3.0    158ms
 33   -36.74251476564   +   -8.15       -5.11    3.0    183ms
 34   -36.74251477290       -8.14       -5.86    3.0    148ms
 35   -36.74251477302       -9.90       -6.22    2.0    113ms
 36   -36.74251477303      -10.99       -6.47    2.0    139ms
 37   -36.74251477304      -11.62       -6.58    2.0    118ms
 38   -36.74251477304   +  -12.36       -6.82    2.0    121ms
 39   -36.74251477303   +  -11.00       -6.39    3.0    150ms
 40   -36.74251477304      -10.95       -6.95    2.0    173ms
 41   -36.74251477304   +    -Inf       -7.21    2.0    113ms
 42   -36.74251477304   +  -12.55       -7.10    2.0    140ms
 43   -36.74251477304      -12.29       -7.80    2.0    115ms
 44   -36.74251477304   +  -12.70       -7.36    3.0    167ms
 45   -36.74251477304      -12.70       -7.87    3.0    148ms
 46   -36.74251477304   +  -13.07       -7.50    3.0    147ms
 47   -36.74251477304      -13.07       -8.06    3.0    146ms
 48   -36.74251477304   +    -Inf       -7.87    3.0    175ms
 49   -36.74251477304   +    -Inf       -8.61    2.0    120ms
 50   -36.74251477304      -14.15       -8.61    3.0    157ms
 51   -36.74251477304      -14.15       -8.88    2.0    114ms
 52   -36.74251477304   +    -Inf       -8.80    2.0    126ms
 53   -36.74251477304   +  -14.15       -9.20    1.0    101ms
 54   -36.74251477304      -14.15       -9.37    3.0    144ms
 55   -36.74251477304   +    -Inf       -9.68    1.0    101ms
 56   -36.74251477304   +    -Inf       -9.60    2.0    161ms
 57   -36.74251477304   +  -13.85       -9.43    2.0    114ms
 58   -36.74251477304      -13.85       -9.47    3.0    147ms
 59   -36.74251477304   +    -Inf      -10.05    2.0    119ms
 60   -36.74251477304      -14.15      -10.13    3.0    140ms
 61   -36.74251477304   +  -14.15       -9.71    3.0    148ms
 62   -36.74251477304   +    -Inf      -10.34    3.0    148ms
 63   -36.74251477304   +    -Inf      -10.30    2.0    171ms
 64   -36.74251477304      -13.85      -10.81    2.0    119ms
 65   -36.74251477304   +  -13.85      -10.93    3.0    133ms
 66   -36.74251477304   +  -13.85      -11.16    1.0   96.2ms
 67   -36.74251477304      -13.85      -11.48    3.0    142ms
 68   -36.74251477304      -14.15      -11.75    2.0    139ms
 69   -36.74251477304   +  -13.67      -11.49    3.0    156ms
 70   -36.74251477304      -13.85      -11.40    3.0    143ms
 71   -36.74251477304   +  -14.15      -11.60    3.0    183ms
 72   -36.74251477304      -14.15      -11.86    2.0    129ms
 73   -36.74251477304   +  -13.85      -12.11    2.0    113ms

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.73331700589                   -0.88   12.0    391ms
  2   -36.74011291797       -2.17       -1.36    1.0   91.9ms
  3   -36.74025920102       -3.83       -1.67    3.0    134ms
  4   -36.74231924412       -2.69       -2.27    1.0    132ms
  5   -36.74241528551       -4.02       -2.56    6.0    148ms
  6   -36.74246308553       -4.32       -2.52    2.0    135ms
  7   -36.74249480403       -4.50       -2.86    1.0   98.8ms
  8   -36.74250905465       -4.85       -3.10    1.0   95.4ms
  9   -36.74251173199       -5.57       -3.34    1.0    100ms
 10   -36.74251435631       -5.58       -3.90    1.0    100ms
 11   -36.74251469664       -6.47       -4.31    3.0    140ms
 12   -36.74251476038       -7.20       -4.54    3.0    161ms
 13   -36.74251476666       -8.20       -4.87    5.0    150ms
 14   -36.74251476940       -8.56       -5.04    3.0    148ms
 15   -36.74251477091       -8.82       -5.30    2.0    116ms
 16   -36.74251477284       -8.71       -5.82    1.0    101ms
 17   -36.74251477301       -9.77       -6.21    3.0    153ms
 18   -36.74251477304      -10.62       -6.55    2.0    140ms
 19   -36.74251477304   +  -12.17       -6.68    2.0    108ms
 20   -36.74251477304      -11.99       -6.97    1.0    143ms
 21   -36.74251477304      -12.58       -7.53    2.0    142ms
 22   -36.74251477304      -13.67       -7.71    2.0    109ms
 23   -36.74251477304   +  -14.15       -7.68    3.0    141ms
 24   -36.74251477304      -14.15       -7.98    1.0    101ms
 25   -36.74251477304   +  -14.15       -8.47    2.0    115ms
 26   -36.74251477304   +    -Inf       -8.77    3.0    142ms
 27   -36.74251477304   +    -Inf       -9.10    2.0    135ms
 28   -36.74251477304      -14.15       -9.33    2.0    116ms
 29   -36.74251477304   +    -Inf       -9.81    2.0    147ms
 30   -36.74251477304   +  -14.15      -10.13    3.0    152ms
 31   -36.74251477304      -14.15      -10.48    2.0    104ms
 32   -36.74251477304   +  -14.15      -10.57    3.0    147ms
 33   -36.74251477304      -14.15      -11.08    1.0    102ms
 34   -36.74251477304   +  -14.15      -11.22    3.0    148ms
 35   -36.74251477304   +    -Inf      -11.85    2.0    120ms
 36   -36.74251477304      -14.15      -11.92    3.0    200ms
 37   -36.74251477304   +  -14.15      -12.28    1.0    103ms

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

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

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

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