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.73236537861                   -0.88   12.0    1.27s
  2   -36.72913387965   +   -2.49       -1.59    1.0    242ms
  3   -25.78846525899   +    1.04       -0.52    5.0    209ms
  4   -32.35621007314        0.82       -0.71    5.0    183ms
  5   -36.45878627397        0.61       -1.20    3.0    149ms
  6   -35.33651127387   +    0.05       -0.98    3.0    139ms
  7   -36.72046830986        0.14       -1.78    3.0    130ms
  8   -36.73635221147       -1.80       -1.80    2.0    116ms
  9   -36.74154077164       -2.28       -2.12    2.0    104ms
 10   -36.74161805910       -4.11       -2.35    1.0   95.1ms
 11   -36.74239205602       -3.11       -2.65    2.0    125ms
 12   -36.74221484890   +   -3.75       -2.48    2.0    108ms
 13   -36.74239080961       -3.75       -2.99    1.0   92.7ms
 14   -36.74247345588       -4.08       -3.22    2.0    110ms
 15   -36.74246220735   +   -4.95       -3.35    2.0    129ms
 16   -36.74247799634       -4.80       -3.62    2.0    101ms
 17   -36.74242265117   +   -4.26       -3.18    3.0    137ms
 18   -36.74234901979   +   -4.13       -3.01    3.0    146ms
 19   -36.74247715760       -3.89       -3.78    3.0    132ms
 20   -36.74248056082       -5.47       -4.39    2.0    205ms
 21   -36.74248048149   +   -7.10       -4.23    3.0    141ms
 22   -36.74248052906       -7.32       -4.38    2.0    1.00s
 23   -36.74248065303       -6.91       -4.78    1.0   93.0ms
 24   -36.74248066066       -8.12       -4.74    3.0    136ms
 25   -36.74248066652       -8.23       -4.98    2.0    118ms
 26   -36.74248067009       -8.45       -5.08    2.0    100ms
 27   -36.74248066698   +   -8.51       -4.82    2.0    114ms
 28   -36.74248067216       -8.29       -5.49    2.0    104ms
 29   -36.74248067133   +   -9.08       -5.44    3.0    132ms
 30   -36.74248066905   +   -8.64       -5.28    3.0    130ms
 31   -36.74248067076       -8.77       -5.40    3.0    151ms
 32   -36.74248067260       -8.74       -6.01    2.0    141ms
 33   -36.74248067268      -10.12       -6.42    2.0    140ms
 34   -36.74248067265   +  -10.51       -6.11    3.0    160ms
 35   -36.74248067268      -10.46       -6.82    2.0    119ms
 36   -36.74248067268      -12.19       -7.13    2.0    126ms
 37   -36.74248067268      -12.97       -7.38    2.0    110ms
 38   -36.74248067268      -13.07       -7.46    3.0    109ms
 39   -36.74248067268   +    -Inf       -7.59    1.0   92.4ms
 40   -36.74248067268      -13.55       -7.73    2.0    131ms
 41   -36.74248067268   +  -13.15       -7.53    3.0    133ms
 42   -36.74248067268      -13.19       -7.77    3.0    125ms
 43   -36.74248067268      -14.15       -7.85    3.0    120ms
 44   -36.74248067268      -13.85       -8.38    2.0    102ms
 45   -36.74248067268      -14.15       -8.43    3.0    136ms
 46   -36.74248067268   +  -14.15       -8.97    2.0    115ms
 47   -36.74248067268      -14.15       -8.91    3.0    152ms
 48   -36.74248067268      -14.15       -9.28    2.0    108ms
 49   -36.74248067268   +    -Inf       -9.07    3.0    131ms
 50   -36.74248067268   +    -Inf       -9.68    2.0    108ms
 51   -36.74248067268   +  -13.85       -9.47    3.0    145ms
 52   -36.74248067268      -14.15       -9.77    2.0    110ms
 53   -36.74248067268   +  -14.15      -10.05    2.0    129ms
 54   -36.74248067268   +  -14.15      -10.30    2.0    101ms
 55   -36.74248067268      -14.15      -10.27    3.0    138ms
 56   -36.74248067268   +    -Inf      -10.52    2.0    121ms
 57   -36.74248067268   +    -Inf      -10.44    3.0    128ms
 58   -36.74248067268   +    -Inf      -11.05    2.0    106ms
 59   -36.74248067268      -13.85      -10.73    3.0    145ms
 60   -36.74248067268   +  -13.85      -11.44    3.0    133ms
 61   -36.74248067268   +    -Inf      -10.74    4.0    164ms
 62   -36.74248067268      -14.15      -11.93    4.0    163ms
 63   -36.74248067268   +  -14.15      -11.67    3.0    132ms
 64   -36.74248067268   +    -Inf      -11.83    3.0    137ms
 65   -36.74248067268      -13.85      -12.13    1.0   96.3ms

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.73245619493                   -0.88   12.0    920ms
  2   -36.73974055435       -2.14       -1.36    1.0    611ms
  3   -36.74185523301       -2.67       -2.10    3.0    148ms
  4   -36.74212254955       -3.57       -2.07    6.0    128ms
  5   -36.74237761804       -3.59       -2.55    1.0   91.0ms
  6   -36.74245167953       -4.13       -2.82    2.0    107ms
  7   -36.74247026182       -4.73       -2.97    3.0    117ms
  8   -36.74247972378       -5.02       -3.59    1.0   92.8ms
  9   -36.74247982705       -6.99       -3.47    4.0    154ms
 10   -36.74248048793       -6.18       -3.82    2.0   98.4ms
 11   -36.74248054989       -7.21       -4.22    2.0    104ms
 12   -36.74248062785       -7.11       -4.60    5.0    126ms
 13   -36.74248067151       -7.36       -4.96    3.0    133ms
 14   -36.74248067211       -9.22       -5.07    2.0    104ms
 15   -36.74248067254       -9.36       -5.41    3.0    113ms
 16   -36.74248067266       -9.91       -5.80    2.0    107ms
 17   -36.74248067268      -10.88       -6.00    3.0    138ms
 18   -36.74248067268      -11.57       -6.26    2.0    104ms
 19   -36.74248067268      -11.69       -6.60    3.0    113ms
 20   -36.74248067268      -12.18       -6.94    2.0    104ms
 21   -36.74248067268      -13.03       -7.19    4.0    127ms
 22   -36.74248067268      -13.85       -7.57    5.0    115ms
 23   -36.74248067268      -13.85       -7.57    2.0    224ms
 24   -36.74248067268      -13.85       -8.12    1.0   94.9ms
 25   -36.74248067268   +    -Inf       -8.27    2.0    1.02s
 26   -36.74248067268   +  -14.15       -8.62    2.0    102ms
 27   -36.74248067268      -13.85       -9.09    1.0   95.0ms
 28   -36.74248067268   +  -13.85       -9.29    4.0    141ms
 29   -36.74248067268   +    -Inf       -9.47    2.0    100ms
 30   -36.74248067268      -14.15       -9.87    2.0    101ms
 31   -36.74248067268   +  -14.15      -10.08    3.0    134ms
 32   -36.74248067268   +    -Inf      -10.44    2.0    102ms
 33   -36.74248067268   +    -Inf      -10.70    3.0    154ms
 34   -36.74248067268      -14.15      -10.91    1.0    117ms
 35   -36.74248067268   +  -13.67      -11.16    2.0    123ms
 36   -36.74248067268      -13.85      -11.46    3.0    138ms
 37   -36.74248067268   +    -Inf      -11.81    3.0    160ms
 38   -36.74248067268   +    -Inf      -12.15    2.0    115ms

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

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

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

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

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