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

al_supercell = bulk(:Al; cubic=true) * (4, 1, 1)
system_Al = attach_psp(al_supercell;
                       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_DFT(system_Al; functionals=LDA(), 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.98030221165                   -0.86   12.0    395ms
  2   -35.94563716634   +   -1.46       -1.53    1.0   91.9ms
  3   +7.968857763187   +    1.64       -0.25    7.0    228ms
  4   -35.90090928684        1.64       -1.17    9.0    233ms
  5   -35.68700922690   +   -0.67       -1.28    3.0    143ms
  6   -35.30007109112   +   -0.41       -1.14    4.0    155ms
  7   -35.97387483223       -0.17       -1.82    4.0    156ms
  8   -35.98812496510       -1.85       -2.11    2.0    118ms
  9   -35.98537913505   +   -2.56       -2.09    2.0    130ms
 10   -35.98783421416       -2.61       -2.12    2.0    128ms
 11   -35.98942964347       -2.80       -2.40    2.0    115ms
 12   -35.98946130373       -4.50       -2.76    1.0    102ms
 13   -35.98962048581       -3.80       -3.08    1.0   97.5ms
 14   -35.98962699496       -5.19       -3.25    4.0    145ms
 15   -35.98915640365   +   -3.33       -2.73    3.0    146ms
 16   -35.98958053758       -3.37       -3.08    4.0    199ms
 17   -35.98956923684   +   -4.95       -3.10    3.0    140ms
 18   -35.98917368079   +   -3.40       -2.73    3.0    147ms
 19   -35.98962819307       -3.34       -3.63    4.0    156ms
 20   -35.98963172448       -5.45       -4.14    2.0    131ms
 21   -35.98963187538       -6.82       -4.20    2.0    134ms
 22   -35.98963195790       -7.08       -4.56    1.0    101ms
 23   -35.98963197843       -7.69       -4.79    3.0    163ms
 24   -35.98963198574       -8.14       -5.16    2.0    137ms
 25   -35.98963198489   +   -9.07       -5.31    2.0    108ms
 26   -35.98963198540       -9.29       -5.47    2.0    104ms
 27   -35.98963198582       -9.37       -5.75    1.0    101ms
 28   -35.98963193318   +   -7.28       -4.70    5.0    201ms
 29   -35.98963198604       -7.28       -6.01    4.0    189ms
 30   -35.98963198567   +   -9.43       -5.73    3.0    154ms
 31   -35.98963198611       -9.35       -6.35    4.0    162ms
 32   -35.98963198612      -11.02       -6.54    2.0    138ms
 33   -35.98963198613      -11.21       -6.89    2.0    103ms
 34   -35.98963198613      -12.46       -7.13    2.0    109ms
 35   -35.98963198613      -13.37       -7.16    3.0    144ms
 36   -35.98963198613      -12.85       -7.28    3.0    126ms
 37   -35.98963198613      -12.77       -7.61    2.0    107ms
 38   -35.98963198613      -14.15       -8.08    1.0    101ms
 39   -35.98963198613   +    -Inf       -7.91    3.0    167ms
 40   -35.98963198613      -13.85       -8.27    5.0    141ms
 41   -35.98963198613   +  -13.67       -7.84    4.0    164ms
 42   -35.98963198613      -13.85       -8.56    3.0    147ms
 43   -35.98963198613      -13.85       -8.44    3.0    158ms
 44   -35.98963198613   +  -13.85       -9.09    2.0    115ms
 45   -35.98963198613      -13.85       -8.91    4.0    143ms
 46   -35.98963198613   +  -14.15       -9.38    2.0    123ms
 47   -35.98963198613   +  -14.15       -9.50    3.0    163ms
 48   -35.98963198613   +    -Inf       -9.84    1.0   97.6ms
 49   -35.98963198613      -14.15      -10.17    1.0    101ms
 50   -35.98963198613   +    -Inf      -10.16    2.0    138ms
 51   -35.98963198613   +  -14.15      -10.42    2.0    120ms
 52   -35.98963198613   +    -Inf      -10.72    2.0    107ms
 53   -35.98963198613      -13.85      -10.48    3.0    208ms
 54   -35.98963198613   +  -13.85      -10.71    4.0    798ms
 55   -35.98963198613      -14.15      -10.58    2.0    115ms
 56   -35.98963198613   +    -Inf      -10.38    3.0    162ms
 57   -35.98963198613   +  -14.15      -11.25    3.0    137ms
 58   -35.98963198613      -13.85      -11.62    3.0    149ms
 59   -35.98963198613   +  -13.85      -11.57    3.0    143ms
 60   -35.98963198613      -13.85      -12.07    2.0    104ms

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.98280387627                   -0.86   10.0    377ms
  2   -35.98799533752       -2.28       -1.34    1.0   93.4ms
  3   -35.98856402881       -3.25       -1.94    3.0    129ms
  4   -35.98919142259       -3.20       -2.04    5.0    154ms
  5   -35.98959435790       -3.39       -2.78    1.0    101ms
  6   -35.98958547570   +   -5.05       -2.68   10.0    189ms
  7   -35.98962454826       -4.41       -3.00    1.0   96.2ms
  8   -35.98963104025       -5.19       -3.43    1.0    103ms
  9   -35.98963091812   +   -6.91       -3.37    3.0    125ms
 10   -35.98963169043       -6.11       -3.72    1.0    102ms
 11   -35.98963191268       -6.65       -4.35    2.0    106ms
 12   -35.98963194109       -7.55       -4.59    5.0    170ms
 13   -35.98963198479       -7.36       -4.91    6.0    145ms
 14   -35.98963198158   +   -8.49       -5.17    3.0    120ms
 15   -35.98963198608       -8.35       -5.60    3.0    147ms
 16   -35.98963198612      -10.39       -6.09    4.0    125ms
 17   -35.98963198613      -11.03       -6.38    4.0    158ms
 18   -35.98963198613      -11.80       -6.88    7.0    145ms
 19   -35.98963198613      -12.89       -7.01    3.0    150ms
 20   -35.98963198613      -12.87       -7.33    1.0    109ms
 21   -35.98963198613      -14.15       -7.82    5.0    142ms
 22   -35.98963198613      -13.85       -8.16    4.0    156ms
 23   -35.98963198613   +    -Inf       -8.37    6.0    140ms
 24   -35.98963198613      -14.15       -8.53    3.0    150ms
 25   -35.98963198613   +  -14.15       -9.10    1.0    104ms
 26   -35.98963198613   +  -14.15       -9.45    3.0    149ms
 27   -35.98963198613      -14.15       -9.74    6.0    157ms
 28   -35.98963198613   +  -14.15      -10.14    4.0    135ms
 29   -35.98963198613   +    -Inf      -10.35    3.0    148ms
 30   -35.98963198613      -14.15      -10.66    2.0    111ms
 31   -35.98963198613      -14.15      -11.17    4.0    148ms
 32   -35.98963198613   +  -13.85      -11.51    6.0    169ms
 33   -35.98963198613      -13.85      -11.60    3.0    120ms
 34   -35.98963198613   +  -13.67      -11.90    2.0    111ms
 35   -35.98963198613      -13.67      -12.47    4.0    161ms

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

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

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

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