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 ASEconvert
using DFTK
using LazyArtifacts

ase_Al    = ase.build.bulk("Al"; cubic=true) * pytuple((4, 1, 1))
system_Al = attach_psp(pyconvert(AbstractSystem, ase_Al);
                       Al=artifact"pd_nc_sr_pbe_standard_0.4.1_upf/Al.upf")
FlexibleSystem(Al₁₆, periodic = TTT):
    bounding_box      : [    16.2        0        0;
                                0     4.05        0;
                                0        0     4.05]u"Å"

      .---------------------------------------.  
     /|                                       |  
    * |                                       |  
    |Al        Al        Al        Al         |  
    | |                                       |  
    | .--Al--------Al--------Al--------Al-----.  
    |/    Al        Al        Al        Al   /   
    Al--------Al--------Al--------Al--------*    

and we discretise:

model_Al = model_LDA(system_Al; temperature=1e-3, symmetries=false)
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   -35.98148248935                   -0.86   13.0    341ms
  2   -35.96393481880   +   -1.76       -1.59    1.0    113ms
  3   -5.210587155496   +    1.49       -0.33    6.0    191ms
  4   -35.91660264314        1.49       -1.27    8.0    203ms
  5   -35.95690016469       -1.39       -1.66    2.0    118ms
  6   -35.27740118167   +   -0.17       -1.14    5.0    162ms
  7   -35.95353871405       -0.17       -1.69    4.0    140ms
  8   -35.98932173216       -1.45       -2.39    2.0    106ms
  9   -35.98895834444   +   -3.44       -2.17    4.0    172ms
 10   -35.98837972258   +   -3.24       -2.17    2.0    121ms
 11   -35.98965508036       -2.89       -2.61    2.0    121ms
 12   -35.98973544367       -4.09       -2.83    1.0   92.3ms
 13   -35.98975996143       -4.61       -3.19    1.0    117ms
 14   -35.98975923204   +   -6.14       -3.32    2.0    126ms
 15   -35.98723179396   +   -2.60       -2.37    4.0    157ms
 16   -35.98976458536       -2.60       -3.70    5.0    167ms
 17   -35.98976613087       -5.81       -4.10    2.0   97.0ms
 18   -35.98975481386   +   -4.95       -3.51    4.0    157ms
 19   -35.98976641764       -4.94       -4.24    4.0    139ms
 20   -35.98976662133       -6.69       -4.35    3.0    123ms
 21   -35.98976676579       -6.84       -4.76    2.0    104ms
 22   -35.98976676793       -8.67       -4.66    2.0    125ms
 23   -35.98976678078       -7.89       -5.00    2.0    126ms
 24   -35.98976678269       -8.72       -5.10    1.0   92.3ms
 25   -35.98976678421       -8.82       -5.45    1.0   91.9ms
 26   -35.98976678436       -9.80       -5.82    2.0    109ms
 27   -35.98976678257   +   -8.75       -5.41    4.0    153ms
 28   -35.98976678395       -8.86       -5.67    4.0    144ms
 29   -35.98976678452       -9.24       -6.71    2.0    109ms
 30   -35.98976678451   +  -11.00       -6.47    5.0    218ms
 31   -35.98976678452      -10.97       -6.82    3.0    127ms
 32   -35.98976678452   +  -11.94       -6.88    2.0    110ms
 33   -35.98976678452   +  -11.92       -6.71    3.0    127ms
 34   -35.98976678452      -11.52       -7.17    2.0    105ms
 35   -35.98976678452   +  -13.45       -7.29    3.0    159ms
 36   -35.98976678452      -12.70       -7.45    2.0    109ms
 37   -35.98976678452      -13.45       -7.76    4.0    150ms
 38   -35.98976678452   +    -Inf       -8.02    2.0    127ms
 39   -35.98976678452   +  -13.67       -7.79    2.0    113ms
 40   -35.98976678452      -13.67       -8.49    3.0    126ms
 41   -35.98976678452   +    -Inf       -8.28    2.0    127ms
 42   -35.98976678452   +  -13.85       -7.93    4.0    146ms
 43   -35.98976678452      -14.15       -9.03    4.0    153ms
 44   -35.98976678452   +  -14.15       -8.97    2.0    126ms
 45   -35.98976678452      -13.85       -9.15    2.0    149ms
 46   -35.98976678452   +    -Inf       -8.97    3.0    129ms
 47   -35.98976678452   +  -14.15       -9.65    3.0    126ms
 48   -35.98976678452      -14.15      -10.23    2.0    127ms
 49   -35.98976678452   +    -Inf      -10.10    3.0    136ms
 50   -35.98976678452   +    -Inf      -10.35    2.0   98.1ms
 51   -35.98976678452   +    -Inf      -10.70    2.0    114ms
 52   -35.98976678452   +    -Inf      -10.48    3.0    117ms
 53   -35.98976678452   +  -14.15      -10.65    3.0    126ms
 54   -35.98976678452   +  -14.15      -10.92    2.0    127ms
 55   -35.98976678452      -14.15      -11.23    3.0    134ms
 56   -35.98976678452      -14.15      -11.63    2.0    109ms
 57   -35.98976678452   +    -Inf      -11.63    2.0    127ms
 58   -35.98976678452   +  -13.85      -11.89    2.0    109ms
 59   -35.98976678452   +    -Inf      -11.79    3.0    135ms
 60   -35.98976678452      -13.85      -11.88    3.0    146ms
 61   -35.98976678452      -14.15      -12.05    2.0    110ms

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   -35.98153549458                   -0.86   12.0    355ms
  2   -35.98750838449       -2.22       -1.34    1.0   88.3ms
  3   -35.98606916490   +   -2.84       -1.61    3.0    128ms
  4   -35.98952584638       -2.46       -2.23    1.0   89.5ms
  5   -35.98948597679   +   -4.40       -2.42    7.0    156ms
  6   -35.98973160713       -3.61       -2.56    8.0    136ms
  7   -35.98971741463   +   -4.85       -2.94    1.0    117ms
  8   -35.98975262448       -4.45       -3.02    3.0    132ms
  9   -35.98975434462       -5.76       -3.36    1.0   90.9ms
 10   -35.98976641689       -4.92       -4.04    2.0    104ms
 11   -35.98976655456       -6.86       -4.22    3.0    137ms
 12   -35.98976675633       -6.70       -4.53    3.0    110ms
 13   -35.98976676859       -7.91       -4.82    3.0    110ms
 14   -35.98976678427       -7.80       -5.42    7.0    147ms
 15   -35.98976678446       -9.72       -5.64    3.0    148ms
 16   -35.98976678435   +   -9.94       -5.82    2.0    174ms
 17   -35.98976678444      -10.03       -6.03    7.0    144ms
 18   -35.98976678452      -10.09       -6.38    3.0    137ms
 19   -35.98976678452      -11.55       -6.67    1.0   94.4ms
 20   -35.98976678452      -13.37       -6.97    5.0    124ms
 21   -35.98976678452      -12.41       -7.18    3.0    130ms
 22   -35.98976678452      -13.00       -7.38    2.0    131ms
 23   -35.98976678452      -13.11       -7.95    2.0    104ms
 24   -35.98976678452   +  -14.15       -8.05    4.0    194ms
 25   -35.98976678452      -14.15       -8.33    1.0   94.7ms
 26   -35.98976678452   +    -Inf       -8.90    6.0    129ms
 27   -35.98976678452   +  -14.15       -9.48    4.0    146ms
 28   -35.98976678452   +    -Inf       -9.83    6.0    140ms
 29   -35.98976678452   +  -14.15      -10.24    3.0    141ms
 30   -35.98976678452   +    -Inf      -10.37    3.0    114ms
 31   -35.98976678452      -13.85      -11.02    8.0    180ms
 32   -35.98976678452   +    -Inf      -11.25    3.0    140ms
 33   -35.98976678452   +  -14.15      -11.51    3.0    139ms
 34   -35.98976678452   +  -14.15      -11.65    2.0    103ms
 35   -35.98976678452      -13.85      -12.11    2.0    102ms

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

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

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

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