Creating and modelling metallic supercells

In this section we will be concerned with modelling supercells of aluminium. When dealing with periodic problems there is no unique definition of the lattice: Clearly any duplication of the lattice along an axis is also a valid repetitive unit to describe exactly the same system. This is exactly what a supercell is: An $n$-fold repetition along one of the axes of the original lattice.

The following code achieves this for aluminium:

using DFTK
using LinearAlgebra
using ASEconvert

function aluminium_setup(repeat=1; Ecut=7.0, kgrid=[2, 2, 2])
    a = 7.65339
    lattice = a * Matrix(I, 3, 3)
    Al = ElementPsp(:Al; psp=load_psp("hgh/lda/al-q3"))
    atoms     = [Al, Al, Al, Al]
    positions = [[0.0, 0.0, 0.0], [0.0, 0.5, 0.5], [0.5, 0.0, 0.5], [0.5, 0.5, 0.0]]
    unit_cell = periodic_system(lattice, atoms, positions)

    # Make supercell in ASE:
    # We convert our lattice to the conventions used in ASE, make the supercell
    # and then convert back ...
    supercell_ase = convert_ase(unit_cell) * pytuple((repeat, 1, 1))
    supercell     = pyconvert(AbstractSystem, supercell_ase)

    # Unfortunately right now the conversion to ASE drops the pseudopotential information,
    # so we need to reattach it:
    supercell = attach_psp(supercell; Al="hgh/lda/al-q3")

    # Construct an LDA model and discretise
    # Note: We disable symmetries explicitly here. Otherwise the problem sizes
    #       we are able to run on the CI are too simple to observe the numerical
    #       instabilities we want to trigger here.
    model = model_LDA(supercell; temperature=1e-3, symmetries=false)
    PlaneWaveBasis(model; Ecut, kgrid)
end;

As part of the code we are using a routine inside the ASE, the atomistic simulation environment for creating the supercell and make use of the two-way interoperability of DFTK and ASE. For more details on this aspect see the documentation on Input and output formats.

Write an example supercell structure to a file to plot it:

setup = aluminium_setup(5)
convert_ase(periodic_system(setup.model)).write("al_supercell.png")
Python: None

As we will see in this notebook the modelling of a system generally becomes harder if the system becomes larger.

  • This sounds like a trivial statement as per se the cost per SCF step increases as the system (and thus $N$) gets larger.
  • But there is more to it: If one is not careful also the number of SCF iterations increases as the system gets larger.
  • The aim of a proper computational treatment of such supercells is therefore to ensure that the number of SCF iterations remains constant when the system size increases.

For achieving the latter DFTK by default employs the LdosMixing preconditioner [HL2021] during the SCF iterations. This mixing approach is completely parameter free, but still automatically adapts to the treated system in order to efficiently prevent charge sloshing. As a result, modelling aluminium slabs indeed takes roughly the same number of SCF iterations irrespective of the supercell size:

M. F. Herbst and A. Levitt. Black-box inhomogeneous preconditioning for self-consistent field iterations in density functional theory. J. Phys. Cond. Matt 33 085503 (2021). ArXiv:2009.01665

self_consistent_field(aluminium_setup(1); tol=1e-4);
┌ Warning: Skipping atomic property pseudopotential, which is not supported in ASE.
@ ASEconvert ~/.julia/packages/ASEconvert/CNQ1A/src/ASEconvert.jl:123
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -8.298676176404                   -0.85    5.4    142ms
  2   -8.300222744283       -2.81       -1.25    1.2   75.2ms
  3   -8.300446738470       -3.65       -1.89    3.0   94.2ms
  4   -8.300461671003       -4.83       -2.77    1.0   71.8ms
  5   -8.300464266203       -5.59       -3.08    2.5    132ms
  6   -8.300464501772       -6.63       -3.24    1.4   74.0ms
  7   -8.300464578492       -7.12       -3.39    1.0   67.9ms
  8   -8.300464615307       -7.43       -3.54    1.0   68.4ms
  9   -8.300464636572       -7.67       -3.73    1.2   73.4ms
 10   -8.300464638216       -8.78       -3.77    1.5   76.7ms
 11   -8.300464642793       -8.34       -4.02    1.0   89.8ms
self_consistent_field(aluminium_setup(2); tol=1e-4);
┌ Warning: Skipping atomic property pseudopotential, which is not supported in ASE.
@ ASEconvert ~/.julia/packages/ASEconvert/CNQ1A/src/ASEconvert.jl:123
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -16.64326548602                   -0.70    5.0    351ms
  2   -16.64614328443       -2.54       -1.15    1.0    246ms
  3   -16.67920977948       -1.48       -1.88    1.8    209ms
  4   -16.67927207156       -4.21       -2.69    4.1    248ms
  5   -16.67928574731       -4.86       -3.01    4.8    296ms
  6   -16.67928615047       -6.39       -3.33    2.1    250ms
  7   -16.67928620702       -7.25       -3.73    1.5    180ms
  8   -16.67928621939       -7.91       -4.48    2.4    207ms
self_consistent_field(aluminium_setup(4); tol=1e-4);
┌ Warning: Skipping atomic property pseudopotential, which is not supported in ASE.
@ ASEconvert ~/.julia/packages/ASEconvert/CNQ1A/src/ASEconvert.jl:123
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -33.32864612071                   -0.56    8.0    1.36s
  2   -33.33473045315       -2.22       -1.00    1.5    692ms
  3   -33.33598398590       -2.90       -1.72    5.4    1.11s
  4   -33.33611200965       -3.89       -2.63    3.4    783ms
  5   -33.33675006503       -3.20       -2.19    5.0    1.03s
  6   -33.33678386179       -4.47       -2.23    4.4    768ms
  7   -33.33694160471       -3.80       -3.25    1.4    645ms
  8   -33.33694369306       -5.68       -3.59    3.2    857ms
  9   -33.33694372275       -7.53       -3.62    2.6    767ms
 10   -33.33694357539   +   -6.83       -3.67    1.1    593ms
 11   -33.33694371526       -6.85       -3.93    1.0    619ms
 12   -33.33694377469       -7.23       -4.40    1.9    633ms

When switching off explicitly the LdosMixing, by selecting mixing=SimpleMixing(), the performance of number of required SCF steps starts to increase as we increase the size of the modelled problem:

self_consistent_field(aluminium_setup(1); tol=1e-4, mixing=SimpleMixing());
┌ Warning: Skipping atomic property pseudopotential, which is not supported in ASE.
@ ASEconvert ~/.julia/packages/ASEconvert/CNQ1A/src/ASEconvert.jl:123
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -8.298609023355                   -0.85    5.6    132ms
  2   -8.300276646083       -2.78       -1.59    1.0   60.2ms
  3   -8.300388403971       -3.95       -2.33    2.0   73.1ms
  4   -8.300318466567   +   -4.16       -2.17    3.9    106ms
  5   -8.300464145001       -3.84       -3.52    1.0   73.5ms
  6   -8.300464598581       -6.34       -3.86    3.5    116ms
  7   -8.300464638328       -7.40       -4.26    1.9   85.0ms
self_consistent_field(aluminium_setup(4); tol=1e-4, mixing=SimpleMixing());
┌ Warning: Skipping atomic property pseudopotential, which is not supported in ASE.
@ ASEconvert ~/.julia/packages/ASEconvert/CNQ1A/src/ASEconvert.jl:123
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -33.32658792619                   -0.56    7.0    1.24s
  2   -33.32449273693   +   -2.68       -1.27    1.0    550ms
  3   -21.75949852728   +    1.06       -0.54    6.1    1.13s
  4   -33.27322301781        1.06       -1.56    6.9    1.31s
  5   -33.08572891802   +   -0.73       -1.24    3.1    874ms
  6   -30.35270564076   +    0.44       -0.83    5.0    1.04s
  7   -33.33580603561        0.47       -2.22    4.2    955ms
  8   -33.32965993016   +   -2.21       -1.97    3.2    895ms
  9   -33.33520549005       -2.26       -2.25    2.5    789ms
 10   -33.33600385085       -3.10       -2.54    1.5    612ms
 11   -33.33681709172       -3.09       -2.92    1.9    680ms
 12   -33.33688367793       -4.18       -3.10    3.1    841ms
 13   -33.33692722470       -4.36       -3.34    2.9    747ms
 14   -33.33693902402       -4.93       -3.62    1.9    721ms
 15   -33.33693618208   +   -5.55       -3.62    3.6    965ms
 16   -33.33694337595       -5.14       -4.19    1.8    631ms

For completion let us note that the more traditional mixing=KerkerMixing() approach would also help in this particular setting to obtain a constant number of SCF iterations for an increasing system size (try it!). In contrast to LdosMixing, however, KerkerMixing is only suitable to model bulk metallic system (like the case we are considering here). When modelling metallic surfaces or mixtures of metals and insulators, KerkerMixing fails, while LdosMixing still works well. See the Modelling a gallium arsenide surface example or [HL2021] for details. Due to the general applicability of LdosMixing this method is the default mixing approach in DFTK.