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.73303395601                   -0.88   11.0    408ms
  2   -36.63513273746   +   -1.01       -1.45    1.0   97.7ms
  3   +34.37190970911   +    1.85       -0.14   22.0    352ms
  4   -36.44543980441        1.85       -1.11    8.0    268ms
  5   -36.21769280663   +   -0.64       -1.18    4.0    164ms
  6   -36.23500221181       -1.76       -1.21    4.0    164ms
  7   -36.73330560709       -0.30       -1.86    3.0    149ms
  8   -36.73942727430       -2.21       -2.15    2.0    138ms
  9   -36.73821824025   +   -2.92       -1.95    3.0    151ms
 10   -36.74142642011       -2.49       -2.25    2.0    136ms
 11   -36.74133932667   +   -4.06       -2.51    2.0    110ms
 12   -36.74197114922       -3.20       -2.65    1.0   98.0ms
 13   -36.74207505550       -3.98       -3.31    1.0    101ms
 14   -36.74206301749   +   -4.92       -3.26    3.0    158ms
 15   -36.74192218058   +   -3.85       -2.93    2.0    124ms
 16   -36.74193894979       -4.78       -2.94    4.0    172ms
 17   -36.74202380816       -4.07       -3.12    3.0    139ms
 18   -36.74192768893   +   -4.02       -2.97    3.0    175ms
 19   -36.74207959410       -3.82       -3.73    3.0    157ms
 20   -36.74208341701       -5.42       -3.97    3.0    150ms
 21   -36.74208361547       -6.70       -4.15    1.0   98.4ms
 22   -36.74208362381       -8.08       -4.17    2.0    134ms
 23   -36.74208374839       -6.90       -4.71    1.0    102ms
 24   -36.74208375850       -8.00       -4.92    3.0    160ms
 25   -36.74208376130       -8.55       -4.83    2.0    131ms
 26   -36.74208376338       -8.68       -5.19    1.0    116ms
 27   -36.74208376460       -8.92       -5.69    1.0   98.1ms
 28   -36.74208376380   +   -9.10       -5.50    2.0    141ms
 29   -36.74208376476       -9.02       -6.11    5.0    147ms
 30   -36.74208376375   +   -9.00       -5.56    4.0    180ms
 31   -36.74208376468       -9.03       -6.01    4.0    169ms
 32   -36.74208376474      -10.19       -6.20    3.0    141ms
 33   -36.74208376479      -10.37       -6.73    2.0    127ms
 34   -36.74208376478   +  -11.50       -6.64    3.0    151ms
 35   -36.74208376479      -11.38       -7.00    2.0    122ms
 36   -36.74208376479      -12.20       -7.49    2.0    106ms
 37   -36.74208376479   +  -13.67       -7.41    4.0    158ms
 38   -36.74208376479      -13.67       -7.58    2.0    135ms
 39   -36.74208376479   +    -Inf       -7.74    2.0    121ms
 40   -36.74208376479      -13.55       -8.35    1.0    105ms
 41   -36.74208376479   +  -14.15       -8.29    3.0    193ms
 42   -36.74208376479   +  -14.15       -8.43    2.0    119ms
 43   -36.74208376479      -14.15       -8.81    2.0    110ms
 44   -36.74208376479   +  -14.15       -8.32    4.0    167ms
 45   -36.74208376479   +    -Inf       -9.08    4.0    160ms
 46   -36.74208376479   +    -Inf       -9.01    3.0    141ms
 47   -36.74208376479      -13.85       -9.19    3.0    131ms
 48   -36.74208376479   +  -13.85       -9.75    1.0    103ms
 49   -36.74208376479      -14.15      -10.01    3.0    170ms
 50   -36.74208376479      -14.15       -9.60    4.0    161ms
 51   -36.74208376479   +  -13.85      -10.29    3.0    154ms
 52   -36.74208376479   +  -14.15      -10.44    3.0    152ms
 53   -36.74208376479      -14.15      -10.43    3.0    129ms
 54   -36.74208376479      -13.85      -11.21    1.0    104ms
 55   -36.74208376479   +  -13.85      -11.08    3.0    149ms
 56   -36.74208376479      -14.15      -11.10    3.0    218ms
 57   -36.74208376479   +    -Inf      -11.40    3.0    145ms
 58   -36.74208376479   +  -14.15      -11.73    3.0    131ms
 59   -36.74208376479      -13.85      -11.87    2.0    142ms
 60   -36.74208376479   +  -13.85      -11.94    2.0    122ms
 61   -36.74208376479   +    -Inf      -12.07    2.0    123ms

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.73400539981                   -0.88   12.0    402ms
  2   -36.73998293407       -2.22       -1.36    1.0   96.6ms
  3   -36.73762451251   +   -2.63       -1.56    4.0    140ms
  4   -36.74190354996       -2.37       -2.29    1.0   97.4ms
  5   -36.74193419304       -4.51       -2.50    9.0    179ms
  6   -36.74206317447       -3.89       -2.62    1.0   95.0ms
  7   -36.74206015872   +   -5.52       -3.08    1.0    102ms
  8   -36.74207632176       -4.79       -3.38    4.0    124ms
  9   -36.74208325537       -5.16       -3.57    2.0    166ms
 10   -36.74208358234       -6.49       -4.12    1.0    103ms
 11   -36.74208339031   +   -6.72       -4.18    3.0    124ms
 12   -36.74208375191       -6.44       -4.70    2.0    112ms
 13   -36.74208376448       -7.90       -5.08    8.0    190ms
 14   -36.74208376352   +   -9.02       -5.15    3.0    149ms
 15   -36.74208376450       -9.01       -5.55    1.0    105ms
 16   -36.74208376478       -9.56       -6.13    2.0    117ms
 17   -36.74208376478      -12.21       -6.42    4.0    179ms
 18   -36.74208376479      -11.11       -6.61    3.0    117ms
 19   -36.74208376479   +  -12.11       -6.86    7.0    145ms
 20   -36.74208376479      -11.99       -7.25    2.0    140ms
 21   -36.74208376479   +    -Inf       -7.65    2.0    109ms
 22   -36.74208376479      -13.67       -7.94    3.0    148ms
 23   -36.74208376479   +    -Inf       -8.26    2.0    110ms
 24   -36.74208376479   +    -Inf       -8.70    5.0    132ms
 25   -36.74208376479   +  -14.15       -8.78    3.0    184ms
 26   -36.74208376479   +    -Inf       -9.06    2.0    112ms
 27   -36.74208376479   +    -Inf       -9.81    3.0    127ms
 28   -36.74208376479   +    -Inf       -9.64    8.0    188ms
 29   -36.74208376479      -14.15       -9.81    2.0    111ms
 30   -36.74208376479   +    -Inf      -10.31    1.0    104ms
 31   -36.74208376479   +  -13.85      -10.82    2.0    143ms
 32   -36.74208376479      -13.85      -10.90    6.0    158ms
 33   -36.74208376479   +  -13.85      -11.33    2.0    129ms
 34   -36.74208376479      -14.15      -11.50    3.0    149ms
 35   -36.74208376479   +    -Inf      -11.72    2.0    112ms
 36   -36.74208376479   +    -Inf      -12.07    2.0    108ms

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

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

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

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