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.73404414421                   -0.88   11.0    354ms
  2   -36.67597664890   +   -1.24       -1.47    1.0   92.0ms
  3   +20.99499526877   +    1.76       -0.18    7.0    216ms
  4   -36.00924729617        1.76       -0.93    6.0    221ms
  5   -35.45218695340   +   -0.25       -0.99    4.0    175ms
  6   -35.63768230182       -0.73       -1.03    4.0    165ms
  7   -36.73564802716        0.04       -1.80    3.0    136ms
  8   -36.73837828570       -2.56       -1.91    2.0    116ms
  9   -36.73708142925   +   -2.89       -2.03    2.0    114ms
 10   -36.74207902392       -2.30       -2.21    1.0   95.4ms
 11   -36.74213255168       -4.27       -2.29    1.0   91.5ms
 12   -36.74219379585       -4.21       -2.39    2.0    114ms
 13   -36.74245347953       -3.59       -2.74    1.0   97.6ms
 14   -36.74233757005   +   -3.94       -2.78    2.0    129ms
 15   -36.73856270185   +   -2.42       -2.26    3.0    154ms
 16   -36.74233661347       -2.42       -2.81    3.0    139ms
 17   -36.74237260939       -4.44       -2.97    2.0   99.5ms
 18   -36.73467576266   +   -2.11       -2.12    4.0    165ms
 19   -36.74251353933       -2.11       -3.63    3.0    149ms
 20   -36.74251394794       -6.39       -3.61    2.0    131ms
 21   -36.74251452137       -6.24       -4.10    1.0   92.8ms
 22   -36.74251457664       -7.26       -4.18    2.0    113ms
 23   -36.74251468280       -6.97       -4.51    1.0   99.0ms
 24   -36.74251476124       -7.11       -4.58    2.0    127ms
 25   -36.74251476724       -8.22       -4.94    1.0   97.1ms
 26   -36.74251477004       -8.55       -4.94    2.0    117ms
 27   -36.74251477196       -8.72       -5.29    1.0   92.7ms
 28   -36.74251476653   +   -8.27       -5.12    3.0    140ms
 29   -36.74251477243       -8.23       -5.50    3.0    131ms
 30   -36.74251477138   +   -8.98       -5.43    2.0    118ms
 31   -36.74251477223       -9.07       -5.59    3.0    132ms
 32   -36.74251477283       -9.23       -5.64    2.0    138ms
 33   -36.74251477296       -9.87       -6.07    2.0    109ms
 34   -36.74251477300      -10.51       -6.22    1.0   97.4ms
 35   -36.74251477303      -10.44       -6.60    2.0    132ms
 36   -36.74251477304      -11.30       -6.88    2.0    109ms
 37   -36.74251477304      -12.22       -7.20    2.0    132ms
 38   -36.74251477304   +  -13.67       -7.16    2.0    119ms
 39   -36.74251477304      -12.94       -7.35    2.0    109ms
 40   -36.74251477304      -13.19       -7.54    1.0   98.4ms
 41   -36.74251477304   +  -14.15       -7.71    1.0    101ms
 42   -36.74251477304   +  -14.15       -7.80    2.0    109ms
 43   -36.74251477304   +  -14.15       -7.72    2.0    114ms
 44   -36.74251477304      -13.85       -8.04    2.0    123ms
 45   -36.74251477304      -13.85       -8.55    2.0    128ms
 46   -36.74251477304      -13.85       -8.91    2.0    132ms
 47   -36.74251477304   +  -13.67       -8.04    4.0    167ms
 48   -36.74251477304      -13.85       -8.91    4.0    157ms
 49   -36.74251477304   +    -Inf       -9.14    1.0   94.2ms
 50   -36.74251477304   +    -Inf       -9.52    3.0    128ms
 51   -36.74251477304   +    -Inf       -9.69    2.0    133ms
 52   -36.74251477304   +    -Inf       -9.89    1.0   93.5ms
 53   -36.74251477304   +  -14.15       -9.97    2.0    106ms
 54   -36.74251477304      -14.15      -10.31    1.0   92.8ms
 55   -36.74251477304   +  -14.15       -9.73    4.0    154ms
 56   -36.74251477304      -14.15      -10.39    3.0    140ms
 57   -36.74251477304   +  -13.85      -10.64    2.0    124ms
 58   -36.74251477304   +    -Inf      -11.22    1.0    103ms
 59   -36.74251477304      -13.85      -10.55    4.0    178ms
 60   -36.74251477304   +    -Inf      -11.23    3.0    150ms
 61   -36.74251477304   +    -Inf      -11.69    2.0    110ms
 62   -36.74251477304   +    -Inf      -11.83    2.0    133ms
 63   -36.74251477304   +    -Inf      -11.66    2.0    116ms
 64   -36.74251477304   +  -13.85      -12.02    2.0    102ms

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.73174177943                   -0.88   11.0    359ms
  2   -36.73979149242       -2.09       -1.36    1.0   88.2ms
  3   -36.73894763856   +   -3.07       -1.61    3.0    119ms
  4   -36.74224053186       -2.48       -2.25    1.0   93.1ms
  5   -36.74241027244       -3.77       -2.49    3.0    135ms
  6   -36.74243710889       -4.57       -2.41    5.0    115ms
  7   -36.74249397482       -4.25       -2.92    1.0   94.4ms
  8   -36.74251252781       -4.73       -3.33    1.0   90.6ms
  9   -36.74251390386       -5.86       -3.47    4.0    129ms
 10   -36.74251465584       -6.12       -3.96    1.0    116ms
 11   -36.74251472360       -7.17       -4.40    3.0    118ms
 12   -36.74251472688       -8.48       -4.59    3.0    138ms
 13   -36.74251476044       -7.47       -4.94    3.0    118ms
 14   -36.74251477283       -7.91       -5.64    1.0   97.5ms
 15   -36.74251477262   +   -9.70       -5.59    4.0    148ms
 16   -36.74251477302       -9.41       -5.96    2.0    238ms
 17   -36.74251477303      -10.92       -6.25    3.0    783ms
 18   -36.74251477304      -11.06       -6.62    1.0   96.5ms
 19   -36.74251477304      -12.70       -6.96    2.0    130ms
 20   -36.74251477304      -12.75       -7.28    3.0    114ms
 21   -36.74251477304      -14.15       -7.52    4.0    140ms
 22   -36.74251477304   +    -Inf       -7.78    2.0    122ms
 23   -36.74251477304   +    -Inf       -7.98    2.0    119ms
 24   -36.74251477304   +    -Inf       -8.22    2.0    120ms
 25   -36.74251477304   +  -14.15       -8.57    2.0    136ms
 26   -36.74251477304      -14.15       -9.00    3.0    158ms
 27   -36.74251477304   +    -Inf       -9.13    3.0    117ms
 28   -36.74251477304   +  -13.85       -9.61    1.0   94.0ms
 29   -36.74251477304   +    -Inf       -9.97    4.0    140ms
 30   -36.74251477304   +    -Inf      -10.16    5.0    119ms
 31   -36.74251477304      -13.85      -10.31    2.0    134ms
 32   -36.74251477304   +    -Inf      -11.04    2.0    108ms
 33   -36.74251477304   +    -Inf      -11.10    3.0    137ms
 34   -36.74251477304   +  -14.15      -11.49    1.0   93.3ms
 35   -36.74251477304      -14.15      -11.62    3.0    139ms
 36   -36.74251477304   +    -Inf      -11.91    1.0   93.6ms
 37   -36.74251477304   +  -13.85      -12.14    3.0    120ms

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

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

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

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