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

Feeding an AtomsBase AbstractSystem to DFTK

In this example we construct a silicon system using the ase.build.bulk routine from the atomistic simulation environment (ASE), which is exposed by ASEconvert as an AtomsBase AbstractSystem.

# Construct bulk system and convert to an AbstractSystem
using ASEconvert
system_ase = ase.build.bulk("Si")
system = pyconvert(AbstractSystem, system_ase)
FlexibleSystem(Si₂, periodic = TTT):
    bounding_box      : [       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"Å")

                       
                       
                       
                       
              Si       
                       
          Si           
                       
                       
                       
                       

To use an AbstractSystem in DFTK, we attach pseudopotentials, construct a DFT model, discretise and solve:

system = attach_psp(system; Si="hgh/lda/si-q4")

model  = model_LDA(system; 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.921700005959                   -0.69    6.0    200ms
  2   -7.926162562112       -2.35       -1.22    1.0    166ms
  3   -7.926838584544       -3.17       -2.37    1.9    194ms
  4   -7.926861197704       -4.65       -3.02    2.2    183ms
  5   -7.926861642675       -6.35       -3.36    2.0    164ms
  6   -7.926861666346       -7.63       -3.72    1.5    149ms
  7   -7.926861679188       -7.89       -4.14    1.1    160ms
  8   -7.926861681790       -8.58       -5.12    1.9    232ms
  9   -7.926861681859      -10.16       -5.18    2.9    222ms
 10   -7.926861681872      -10.88       -6.46    1.0    162ms
 11   -7.926861681873      -12.15       -6.63    3.1    251ms
 12   -7.926861681873      -14.45       -7.37    1.0    172ms
 13   -7.926861681873      -14.75       -7.70    2.5    233ms
 14   -7.926861681873   +  -15.05       -8.74    1.8    195ms

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

Reading a system using AtomsIO

using AtomsIO

# Read a file using [AtomsIO](https://github.com/mfherbst/AtomsIO.jl),
# which directly yields an AbstractSystem.
system = load_system("Si.extxyz")

# Now run the LDA calculation:
system = attach_psp(system; Si="hgh/lda/si-q4")
model  = model_LDA(system; 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.921736013788                   -0.69    6.0    326ms
  2   -7.926166973747       -2.35       -1.22    1.0    177ms
  3   -7.926838403035       -3.17       -2.37    1.9    197ms
  4   -7.926861204847       -4.64       -3.02    2.2    232ms
  5   -7.926861637489       -6.36       -3.34    1.8    193ms
  6   -7.926861664486       -7.57       -3.69    1.4    186ms
  7   -7.926861679552       -7.82       -4.18    1.1    194ms
  8   -7.926861681801       -8.65       -5.22    1.9    180ms
  9   -7.926861681861      -10.22       -5.19    3.4    238ms
 10   -7.926861681872      -10.96       -6.28    1.0    164ms
 11   -7.926861681873      -12.55       -6.77    2.6    235ms
 12   -7.926861681873      -14.21       -7.08    1.9    227ms
 13   -7.926861681873      -14.75       -8.05    1.4    175ms

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:
system = attach_psp(system; Si="hgh/lda/si-q4")
model  = model_LDA(system; 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.921731125277                   -0.69    6.1    369ms
  2   -7.926167213830       -2.35       -1.22    1.0    178ms
  3   -7.926842808776       -3.17       -2.37    1.9    231ms
  4   -7.926864673322       -4.66       -3.04    2.4    252ms
  5   -7.926865066954       -6.40       -3.41    2.0    214ms
  6   -7.926865082199       -7.82       -3.81    1.4    174ms
  7   -7.926865090101       -8.10       -4.10    1.4    184ms

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₂, periodic = TTT):
    bounding_box      : [       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₀")

                       
                       
                       
                       
              Si       
                       
          Si           
                       
                       
                       
                       

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; psp=load_psp("hgh/lda/Si-q4"))
atoms     = [Si, Si]
positions = [ones(3)/8, -ones(3)/8]

third_system = atomic_system(lattice, atoms, positions)
FlexibleSystem(Si₂, periodic = TTT):
    bounding_box      : [       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₀")

                       
                       
                       
                       
              Si       
                       
          Si