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.73316158809                   -0.88   11.0    1.34s
  2   -36.63011169611   +   -0.99       -1.43    1.0    252ms
  3   +37.16030281076   +    1.87       -0.13    7.0    267ms
  4   -36.48829654671        1.87       -1.03    9.0    371ms
  5   -36.66018198585       -0.76       -1.41    3.0    132ms
  6   -36.70869232930       -1.31       -1.65    2.0    1.12s
  7   -36.05820460338   +   -0.19       -1.15    3.0    130ms
  8   -36.73635159248       -0.17       -2.04    3.0    140ms
  9   -36.74160055200       -2.28       -2.03    2.0    124ms
 10   -36.74164584657       -4.34       -2.19    1.0   91.1ms
 11   -36.74201411823       -3.43       -2.42    1.0   91.7ms
 12   -36.74196151377   +   -4.28       -2.56    1.0   92.4ms
 13   -36.74243358440       -3.33       -2.91    1.0   92.8ms
 14   -36.74244424302       -4.97       -3.05    2.0    118ms
 15   -36.74183499188   +   -3.22       -2.56    3.0    137ms
 16   -36.74247603733       -3.19       -3.52    3.0    134ms
 17   -36.74247710104       -5.97       -3.52    3.0    130ms
 18   -36.74247759248       -6.31       -3.72    2.0    108ms
 19   -36.74247400484   +   -5.45       -3.64    3.0    127ms
 20   -36.74247209549   +   -5.72       -3.59    3.0    133ms
 21   -36.74248008469       -5.10       -4.14    2.0    117ms
 22   -36.74248016750       -7.08       -4.13    3.0    147ms
 23   -36.74248064856       -6.32       -4.74    2.0    109ms
 24   -36.74248066538       -7.77       -5.00    2.0    126ms
 25   -36.74248065893   +   -8.19       -4.88    2.0    112ms
 26   -36.74248067024       -7.95       -5.17    1.0   93.2ms
 27   -36.74248067225       -8.70       -5.59    3.0    126ms
 28   -36.74248067162   +   -9.20       -5.39    2.0    126ms
 29   -36.74248067249       -9.06       -5.88    2.0    107ms
 30   -36.74248067266       -9.78       -6.07    3.0    128ms
 31   -36.74248067185   +   -9.10       -5.57    4.0    154ms
 32   -36.74248067245       -9.22       -5.83    3.0    135ms
 33   -36.74248067260       -9.83       -6.07    2.0    118ms
 34   -36.74248067268      -10.10       -6.58    2.0    123ms
 35   -36.74248067268      -11.69       -6.79    2.0    127ms
 36   -36.74248067268      -12.10       -7.11    2.0    109ms
 37   -36.74248067268   +  -13.55       -6.97    2.0    126ms
 38   -36.74248067268      -12.36       -7.40    2.0    114ms
 39   -36.74248067268   +  -13.45       -7.48    3.0    121ms
 40   -36.74248067268      -13.19       -7.60    2.0    109ms
 41   -36.74248067268      -14.15       -7.95    1.0   97.4ms
 42   -36.74248067268   +  -13.67       -7.70    3.0    140ms
 43   -36.74248067268      -13.45       -8.29    3.0    134ms
 44   -36.74248067268   +    -Inf       -8.33    2.0    107ms
 45   -36.74248067268      -14.15       -8.74    2.0    108ms
 46   -36.74248067268   +    -Inf       -8.28    3.0    147ms
 47   -36.74248067268   +    -Inf       -8.99    3.0    135ms
 48   -36.74248067268   +  -13.85       -8.83    3.0    126ms
 49   -36.74248067268      -14.15       -9.44    3.0    119ms
 50   -36.74248067268   +    -Inf       -9.29    3.0    147ms
 51   -36.74248067268      -14.15       -9.66    2.0    127ms
 52   -36.74248067268      -14.15       -9.39    3.0    138ms
 53   -36.74248067268   +    -Inf       -9.77    3.0    137ms
 54   -36.74248067268   +    -Inf      -10.13    2.0    110ms
 55   -36.74248067268   +    -Inf      -10.17    2.0    131ms
 56   -36.74248067268   +  -13.85      -10.43    2.0    112ms
 57   -36.74248067268   +    -Inf      -10.85    2.0    104ms
 58   -36.74248067268      -14.15      -10.36    3.0    166ms
 59   -36.74248067268      -14.15      -11.06    4.0    154ms
 60   -36.74248067268   +    -Inf      -10.92    2.0    121ms
 61   -36.74248067268   +  -13.85      -11.60    2.0    113ms
 62   -36.74248067268   +    -Inf      -11.34    3.0    137ms
 63   -36.74248067268      -14.15      -11.49    3.0    137ms
 64   -36.74248067268   +    -Inf      -11.20    3.0    129ms
 65   -36.74248067268   +  -13.67      -11.97    3.0    133ms
 66   -36.74248067268      -13.85      -12.34    2.0    130ms

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.73401554671                   -0.88   11.0    969ms
  2   -36.74020815422       -2.21       -1.36    1.0    662ms
  3   -36.74116672210       -3.02       -1.88    3.0    120ms
  4   -36.74209216711       -3.03       -2.09    2.0    109ms
  5   -36.74244311421       -3.45       -2.86    2.0    104ms
  6   -36.74245845612       -4.81       -2.78    4.0    155ms
  7   -36.74247532500       -4.77       -3.01    1.0   96.4ms
  8   -36.74247995255       -5.33       -3.54    1.0    197ms
  9   -36.74248038250       -6.37       -3.79    2.0    107ms
 10   -36.74248064216       -6.59       -4.11    2.0    1.07s
 11   -36.74248065228       -7.99       -4.33    3.0    116ms
 12   -36.74248066818       -7.80       -4.58    5.0    117ms
 13   -36.74248066948       -8.89       -4.67    2.0    128ms
 14   -36.74248067170       -8.65       -5.05    1.0   94.9ms
 15   -36.74248067211       -9.39       -5.31    2.0    120ms
 16   -36.74248067261       -9.30       -5.52    1.0   95.0ms
 17   -36.74248067265      -10.36       -5.85    2.0    104ms
 18   -36.74248067268      -10.57       -6.33    3.0    126ms
 19   -36.74248067268      -11.59       -6.71    2.0    129ms
 20   -36.74248067268      -12.79       -6.97    5.0    137ms
 21   -36.74248067268   +  -13.45       -6.89    2.0    143ms
 22   -36.74248067268      -12.87       -7.27    1.0    119ms
 23   -36.74248067268      -13.25       -7.55    1.0    120ms
 24   -36.74248067268   +    -Inf       -8.03    3.0    123ms
 25   -36.74248067268   +    -Inf       -8.31    3.0    125ms
 26   -36.74248067268   +    -Inf       -8.38    2.0    129ms
 27   -36.74248067268      -13.85       -8.63    1.0   94.8ms
 28   -36.74248067268   +  -13.85       -9.21    2.0    110ms
 29   -36.74248067268   +  -13.85       -9.41    3.0    132ms
 30   -36.74248067268      -13.67       -9.79    2.0    100ms
 31   -36.74248067268   +  -14.15      -10.14    3.0    134ms
 32   -36.74248067268   +  -14.15      -10.51    2.0    100ms
 33   -36.74248067268      -14.15      -10.66    6.0    135ms
 34   -36.74248067268      -14.15      -10.88    2.0    115ms
 35   -36.74248067268      -14.15      -11.29    1.0   94.7ms
 36   -36.74248067268   +  -13.85      -11.73    3.0    135ms
 37   -36.74248067268      -13.85      -12.05    2.0    109ms

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

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

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

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