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₂
    pseudopot. family    : PseudoFamily("dojo.nc.sr.lda.v0_4_1.standard.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.921724718075                   -0.69    5.5    203ms
  2   -7.926132416715       -2.36       -1.22    1.0    186ms
  3   -7.926834148028       -3.15       -2.37    2.0    167ms
  4   -7.926861305475       -4.57       -3.02    3.1    207ms
  5   -7.926861655820       -6.46       -3.47    2.2    194ms
  6   -7.926861675616       -7.70       -4.05    1.8    159ms
  7   -7.926861678481       -8.54       -4.09    1.6    165ms
  8   -7.926861680452       -8.71       -4.25    1.0    145ms
  9   -7.926861681799       -8.87       -4.59    1.0    235ms
 10   -7.926861681859      -10.22       -4.78    1.1    147ms
 11   -7.926861681868      -11.03       -5.44    1.1    973ms
 12   -7.926861681872      -11.41       -6.08    1.8    161ms
 13   -7.926861681873      -12.57       -5.86    2.2    176ms
 14   -7.926861681873      -13.40       -6.29    1.0    147ms
 15   -7.926861681873      -13.97       -7.51    1.0    145ms
 16   -7.926861681873      -15.05       -7.84    2.5    184ms
 17   -7.926861681873   +    -Inf       -7.29    1.2    150ms
 18   -7.926861681873      -15.05       -8.91    1.4    154ms

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.921693905046                   -0.69    5.6    290ms
  2   -7.926130331243       -2.35       -1.22    1.0    147ms
  3   -7.926834064790       -3.15       -2.37    2.0    166ms
  4   -7.926861313773       -4.56       -3.02    3.0    207ms
  5   -7.926861656641       -6.46       -3.49    1.8    202ms
  6   -7.926861676443       -7.70       -4.11    2.1    169ms
  7   -7.926861678045       -8.80       -4.08    1.9    168ms
  8   -7.926861678671       -9.20       -4.07    1.0    147ms
  9   -7.926861681702       -8.52       -4.11    1.0    154ms
 10   -7.926861681833       -9.88       -4.49    1.0    149ms
 11   -7.926861681841      -11.08       -4.46    1.0    149ms
 12   -7.926861681850      -11.07       -4.45    1.0    154ms
 13   -7.926861681871      -10.67       -5.51    1.0    148ms
 14   -7.926861681872      -12.15       -5.75    1.9    169ms
 15   -7.926861681871   +  -11.83       -5.26    1.4    156ms
 16   -7.926861681872      -11.74       -5.63    1.0    155ms
 17   -7.926861681872   +  -13.50       -5.45    1.0    149ms
 18   -7.926861681873      -12.68       -6.71    1.0    157ms
 19   -7.926861681873   +  -15.05       -6.56    2.0    173ms
 20   -7.926861681873      -13.91       -7.50    1.0    154ms
 21   -7.926861681873      -14.75       -7.58    2.4    182ms
 22   -7.926861681873   +  -14.75       -7.83    1.1    151ms
 23   -7.926861681873   +  -15.05       -8.07    1.1    157ms

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.921742373804                   -0.69    5.6    262ms
  2   -7.926135799980       -2.36       -1.22    1.0    133ms
  3   -7.926836164265       -3.15       -2.37    2.0    181ms
  4   -7.926864716894       -4.54       -2.98    2.8    204ms
  5   -7.926865044185       -6.49       -3.34    2.0    162ms
  6   -7.926865077701       -7.47       -3.75    1.5    150ms
  7   -7.926865090890       -7.88       -4.24    1.2    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₀")