AtomsBase integration

AtomsBase.jl is a common interface for representing atomic structures in Julia. DFTK directly supports using such structures to run a calculation as is demonstrated here.

using DFTK
using AtomsBuilder

Feeding an AtomsBase AbstractSystem to DFTK

In this example we construct a bulk silicon system using the bulk function from AtomsBuilder. This function uses tabulated data to set up a reasonable starting geometry and lattice for bulk silicon.

system = bulk(:Si)
FlexibleSystem(Si₂, periodicity = TTT):
    cell_vectors      : [       0    2.715    2.715;
                            2.715        0    2.715;
                            2.715    2.715        0]u"Å"

    Atom(Si, [       0,        0,        0]u"Å")
    Atom(Si, [  1.3575,   1.3575,   1.3575]u"Å")

By default the atoms of an AbstractSystem employ the bare Coulomb potential. To employ pseudpotential models (which is almost always advisable for plane-wave DFT) one employs the pseudopotential keyword argument in model constructors such as model_DFT. For example we can employ a PseudoFamily object from the PseudoPotentialData package. See its documentation for more information on the available pseudopotential families and how to select them.

using PseudoPotentialData  # defines PseudoFamily

pd_lda_family = PseudoFamily("dojo.nc.sr.lda.v0_4_1.standard.upf")
model = model_DFT(system; functionals=LDA(), temperature=1e-3,
                  pseudopotentials=pd_lda_family)
Model(lda_x+lda_c_pw, 3D):
    lattice (in Bohr)    : [0         , 5.13061   , 5.13061   ]
                           [5.13061   , 0         , 5.13061   ]
                           [5.13061   , 5.13061   , 0         ]
    unit cell volume     : 270.11 Bohr³

    atoms                : Si₂
    atom potentials      : ElementPsp(Si, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Si.upf")
                           ElementPsp(Si, "/home/runner/.julia/artifacts/326db5c901e2681584ec5c06fc17f6c96e516ff9/Si.upf")

    num. electrons       : 8
    spin polarization    : none
    temperature          : 0.001 Ha
    smearing             : DFTK.Smearing.FermiDirac()

    terms                : Kinetic()
                           AtomicLocal()
                           AtomicNonlocal()
                           Ewald(nothing)
                           PspCorrection()
                           Hartree()
                           Xc(lda_x, lda_c_pw)
                           Entropy()

Alternatively the pseudopotentials object also accepts a Dict{Symbol,String}, which provides for each element symbol the filename or identifier of the pseudopotential to be employed, e.g.

path_to_pspfile = PseudoFamily("cp2k.nc.sr.lda.v0_1.semicore.gth")[:Si]
model = model_DFT(system; functionals=LDA(), temperature=1e-3,
                  pseudopotentials=Dict(:Si => path_to_pspfile))
Model(lda_x+lda_c_pw, 3D):
    lattice (in Bohr)    : [0         , 5.13061   , 5.13061   ]
                           [5.13061   , 0         , 5.13061   ]
                           [5.13061   , 5.13061   , 0         ]
    unit cell volume     : 270.11 Bohr³

    atoms                : Si₂
    atom potentials      : ElementPsp(Si, "/home/runner/.julia/artifacts/966fd9cdcd7dbaba6dc2bf43ee50dd81e63e8837/Si.gth")
                           ElementPsp(Si, "/home/runner/.julia/artifacts/966fd9cdcd7dbaba6dc2bf43ee50dd81e63e8837/Si.gth")

    num. electrons       : 8
    spin polarization    : none
    temperature          : 0.001 Ha
    smearing             : DFTK.Smearing.FermiDirac()

    terms                : Kinetic()
                           AtomicLocal()
                           AtomicNonlocal()
                           Ewald(nothing)
                           PspCorrection()
                           Hartree()
                           Xc(lda_x, lda_c_pw)
                           Entropy()

We can then discretise such a model and solve:

basis  = PlaneWaveBasis(model; Ecut=15, kgrid=[4, 4, 4])
scfres = self_consistent_field(basis, tol=1e-8);
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -7.921749810259                   -0.69    5.6    212ms
  2   -7.926135165825       -2.36       -1.22    1.0    155ms
  3   -7.926833634746       -3.16       -2.37    2.0    150ms
  4   -7.926861268133       -4.56       -3.00    3.1    195ms
  5   -7.926861646513       -6.42       -3.39    1.9    228ms
  6   -7.926861670708       -7.62       -3.86    1.5    138ms
  7   -7.926861679426       -8.06       -4.17    1.8    775ms
  8   -7.926861681773       -8.63       -4.98    1.1    135ms
  9   -7.926861681863      -10.05       -5.27    2.9    168ms
 10   -7.926861681872      -11.06       -6.02    1.2    134ms
 11   -7.926861681873      -11.97       -6.46    2.1    160ms
 12   -7.926861681873      -13.47       -7.38    1.8    148ms
 13   -7.926861681873      -14.45       -7.84    3.0    175ms
 14   -7.926861681873      -15.05       -8.61    1.9    146ms

If we did not want to use AtomsBuilder we could of course use any other package which yields an AbstractSystem object. This includes:

Reading a system using AtomsIO

Read a file using AtomsIO, which directly yields an AbstractSystem.

using AtomsIO
system = load_system("Si.extxyz");

Run the LDA calculation:

pseudopotentials = PseudoFamily("cp2k.nc.sr.lda.v0_1.semicore.gth")
model  = model_DFT(system; pseudopotentials, functionals=LDA(), temperature=1e-3)
basis  = PlaneWaveBasis(model; Ecut=15, kgrid=[4, 4, 4])
scfres = self_consistent_field(basis, tol=1e-8);
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -7.921746927738                   -0.69    5.6    181ms
  2   -7.926132185743       -2.36       -1.22    1.0    147ms
  3   -7.926833320918       -3.15       -2.37    2.1    151ms
  4   -7.926861305142       -4.55       -2.99    2.8    194ms
  5   -7.926861642135       -6.47       -3.37    2.2    162ms
  6   -7.926861670211       -7.55       -3.82    1.6    139ms
  7   -7.926861679929       -8.01       -4.23    1.6    144ms
  8   -7.926861681788       -8.73       -5.00    1.5    143ms
  9   -7.926861681863      -10.13       -5.28    2.6    166ms
 10   -7.926861681871      -11.09       -5.79    1.2    141ms
 11   -7.926861681873      -11.77       -6.50    1.9    154ms
 12   -7.926861681873      -13.01       -7.63    2.2    167ms
 13   -7.926861681873      -14.45       -8.01    3.9    188ms

The same could be achieved using ExtXYZ by system = Atoms(read_frame("Si.extxyz")), since the ExtXYZ.Atoms object is directly AtomsBase-compatible.

Directly setting up a system in AtomsBase

using AtomsBase
using Unitful
using UnitfulAtomic

# Construct a system in the AtomsBase world
a = 10.26u"bohr"  # Silicon lattice constant
lattice = a / 2 * [[0, 1, 1.],  # Lattice as vector of vectors
                   [1, 0, 1.],
                   [1, 1, 0.]]
atoms  = [:Si => ones(3)/8, :Si => -ones(3)/8]
system = periodic_system(atoms, lattice; fractional=true)

# Now run the LDA calculation:
pseudopotentials = PseudoFamily("cp2k.nc.sr.lda.v0_1.semicore.gth")
model  = model_DFT(system; pseudopotentials, functionals=LDA(), temperature=1e-3)
basis  = PlaneWaveBasis(model; Ecut=15, kgrid=[4, 4, 4])
scfres = self_consistent_field(basis, tol=1e-4);
n     Energy            log10(ΔE)   log10(Δρ)   Diag   Δtime
---   ---------------   ---------   ---------   ----   ------
  1   -7.921725854050                   -0.69    5.6    227ms
  2   -7.926136221338       -2.36       -1.22    1.0    123ms
  3   -7.926836955122       -3.15       -2.37    2.0    150ms
  4   -7.926864480735       -4.56       -2.99    3.1    193ms
  5   -7.926865063410       -6.23       -3.38    1.9    151ms
  6   -7.926865082719       -7.71       -3.86    1.4    129ms
  7   -7.926865089120       -8.19       -4.03    1.6    136ms

Obtaining an AbstractSystem from DFTK data

At any point we can also get back the DFTK model as an AtomsBase-compatible AbstractSystem:

second_system = atomic_system(model)
FlexibleSystem(Si₂, periodicity = TTT):
    cell_vectors      : [       0     5.13     5.13;
                             5.13        0     5.13;
                             5.13     5.13        0]u"a₀"

    Atom(Si, [  1.2825,   1.2825,   1.2825]u"a₀")
    Atom(Si, [ -1.2825,  -1.2825,  -1.2825]u"a₀")

Similarly DFTK offers a method to the atomic_system and periodic_system functions (from AtomsBase), which enable a seamless conversion of the usual data structures for setting up DFTK calculations into an AbstractSystem:

lattice = 5.431u"Å" / 2 * [[0 1 1.];
                           [1 0 1.];
                           [1 1 0.]];
Si = ElementPsp(:Si, pseudopotentials)
atoms     = [Si, Si]
positions = [ones(3)/8, -ones(3)/8]

third_system = atomic_system(lattice, atoms, positions)
FlexibleSystem(Si₂, periodicity = TTT):
    cell_vectors      : [       0  5.13155  5.13155;
                          5.13155        0  5.13155;
                          5.13155  5.13155        0]u"a₀"

    Atom(Si, [ 1.28289,  1.28289,  1.28289]u"a₀")
    Atom(Si, [-1.28289, -1.28289, -1.28289]u"a₀")