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₂, 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           
                       
                       
                       
                       

By default the atoms of an AbstractSystem employ the bare Coulomb potential. To make calculations feasible for plane-wave DFT we thus attach pseudopotential information, before passing the system to construct a DFT model, discretise and solve:

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

model  = model_DFT(system; 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.921726746289                   -0.69    6.0    256ms
  2   -7.926164308317       -2.35       -1.22    1.0    155ms
  3   -7.926839237591       -3.17       -2.37    1.9    176ms
  4   -7.926861199790       -4.66       -3.04    2.1    189ms
  5   -7.926861652648       -6.34       -3.41    2.0    171ms
  6   -7.926861671100       -7.73       -3.82    1.6    206ms
  7   -7.926861678324       -8.14       -4.05    1.1    149ms
  8   -7.926861681780       -8.46       -4.88    1.4    154ms
  9   -7.926861681858      -10.11       -5.18    2.8    183ms
 10   -7.926861681872      -10.86       -6.19    1.5    158ms
 11   -7.926861681873      -11.98       -6.43    2.9    198ms
 12   -7.926861681873      -14.21       -6.28    1.0    150ms
 13   -7.926861681873      -14.35       -7.35    1.0    148ms
 14   -7.926861681873   +  -15.05       -7.47    2.5    184ms
 15   -7.926861681873      -14.35       -8.56    1.0    152ms

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_DFT(system; 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.921712243822                   -0.69    6.0    232ms
  2   -7.926163011018       -2.35       -1.22    1.0    157ms
  3   -7.926838940325       -3.17       -2.37    2.0    182ms
  4   -7.926861252050       -4.65       -3.03    2.0    186ms
  5   -7.926861648028       -6.40       -3.39    1.9    169ms
  6   -7.926861668698       -7.68       -3.76    1.6    170ms
  7   -7.926861679070       -7.98       -4.12    1.1    149ms
  8   -7.926861681795       -8.56       -5.11    1.8    185ms
  9   -7.926861681855      -10.22       -5.14    3.0    198ms
 10   -7.926861681872      -10.77       -6.29    1.0    149ms
 11   -7.926861681873      -12.13       -6.54    2.8    257ms
 12   -7.926861681873      -14.05       -6.94    1.0    173ms
 13   -7.926861681873   +  -15.05       -7.61    1.6    186ms
 14   -7.926861681873      -14.57       -8.25    2.9    272ms

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_DFT(system; 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.921715946320                   -0.69    6.0    409ms
  2   -7.926169061291       -2.35       -1.22    1.0    190ms
  3   -7.926842414158       -3.17       -2.37    1.9    250ms
  4   -7.926864609915       -4.65       -3.02    2.1    260ms
  5   -7.926865049611       -6.36       -3.34    1.9    230ms
  6   -7.926865075274       -7.59       -3.68    1.5    190ms
  7   -7.926865090567       -7.82       -4.16    1.0    175ms

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