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.921746284713 -0.69 5.5 223ms
2 -7.926135949678 -2.36 -1.22 1.0 230ms
3 -7.926833178059 -3.16 -2.37 2.0 716ms
4 -7.926861183215 -4.55 -2.99 3.0 230ms
5 -7.926861648831 -6.33 -3.37 2.0 165ms
6 -7.926861670831 -7.66 -3.84 1.5 151ms
7 -7.926861678954 -8.09 -4.11 1.9 172ms
8 -7.926861681778 -8.55 -4.97 1.4 157ms
9 -7.926861681857 -10.10 -5.18 2.4 171ms
10 -7.926861681868 -10.96 -5.57 1.2 176ms
11 -7.926861681872 -11.36 -6.32 1.1 149ms
12 -7.926861681873 -12.75 -7.02 2.2 173ms
13 -7.926861681873 -14.35 -7.43 2.1 169ms
14 -7.926861681873 + -15.05 -8.44 1.9 157ms
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.921746828707 -0.69 5.5 198ms
2 -7.926136461420 -2.36 -1.22 1.0 139ms
3 -7.926834397869 -3.16 -2.37 2.0 155ms
4 -7.926861121656 -4.57 -3.01 3.1 193ms
5 -7.926861661353 -6.27 -3.44 1.9 159ms
6 -7.926861674962 -7.87 -4.00 1.6 147ms
7 -7.926861678993 -8.39 -4.11 1.9 153ms
8 -7.926861681468 -8.61 -4.41 1.1 157ms
9 -7.926861681854 -9.41 -5.02 1.0 136ms
10 -7.926861681868 -10.85 -5.49 2.0 156ms
11 -7.926861681872 -11.39 -5.70 1.4 146ms
12 -7.926861681873 -12.61 -6.35 1.0 138ms
13 -7.926861681873 -13.27 -7.42 2.2 165ms
14 -7.926861681873 -14.15 -6.99 3.1 191ms
15 -7.926861681873 + -Inf -7.74 1.0 140ms
16 -7.926861681873 + -14.57 -8.40 1.0 143ms
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.921729493697 -0.69 5.5 229ms
2 -7.926139248420 -2.36 -1.22 1.0 137ms
3 -7.926837775195 -3.16 -2.37 2.0 147ms
4 -7.926864711965 -4.57 -3.01 3.0 195ms
5 -7.926865061246 -6.46 -3.42 1.9 150ms
6 -7.926865084192 -7.64 -3.92 1.6 139ms
7 -7.926865090172 -8.22 -4.13 1.5 131ms
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₀")