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.921738267740 -0.69 5.5 246ms
2 -7.926129632009 -2.36 -1.22 1.0 186ms
3 -7.926833824210 -3.15 -2.37 2.0 198ms
4 -7.926861299604 -4.56 -3.01 2.9 211ms
5 -7.926861656011 -6.45 -3.47 2.0 169ms
6 -7.926861675490 -7.71 -4.05 1.9 172ms
7 -7.926861677804 -8.64 -4.05 2.1 170ms
8 -7.926861677624 + -9.75 -4.04 1.0 155ms
9 -7.926861679234 -8.79 -3.77 1.0 153ms
10 -7.926861551300 + -6.89 -2.57 2.8 195ms
11 -7.926861575340 -7.62 -2.63 1.0 153ms
12 -7.926861577648 -8.64 -2.66 1.0 146ms
13 -7.926861633055 -7.26 -2.84 2.0 172ms
14 -7.926861665128 -7.49 -3.06 1.8 171ms
15 -7.926861676992 -7.93 -3.34 1.1 206ms
16 -7.926861681126 -8.38 -3.74 1.0 165ms
17 -7.926861680879 + -9.61 -3.69 1.0 147ms
18 -7.926861680975 -10.02 -3.70 1.0 166ms
19 -7.926861681431 -9.34 -3.85 1.0 149ms
20 -7.926861681860 -9.37 -6.00 1.0 155ms
21 -7.926861681873 -10.91 -6.74 3.9 224ms
22 -7.926861681873 -13.60 -6.66 2.4 180ms
23 -7.926861681873 + -15.05 -7.17 1.0 156ms
24 -7.926861681873 -14.57 -7.23 1.0 153ms
25 -7.926861681873 + -15.05 -7.16 1.0 149ms
26 -7.926861681873 + -15.05 -7.27 1.0 156ms
27 -7.926861681873 + -15.05 -7.83 1.0 150ms
28 -7.926861681873 -14.45 -9.10 1.1 160ms
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.921678831331 -0.69 5.6 203ms
2 -7.926135708443 -2.35 -1.22 1.0 174ms
3 -7.926834135171 -3.16 -2.37 2.0 168ms
4 -7.926861286565 -4.57 -3.00 2.9 210ms
5 -7.926861650008 -6.44 -3.40 2.0 178ms
6 -7.926861671821 -7.66 -3.89 1.5 154ms
7 -7.926861679205 -8.13 -4.15 1.8 282ms
8 -7.926861681782 -8.59 -4.93 1.5 154ms
9 -7.926861681851 -10.16 -5.13 2.1 1.37s
10 -7.926861681868 -10.75 -5.59 1.5 157ms
11 -7.926861681872 -11.37 -6.25 1.2 155ms
12 -7.926861681873 -12.74 -6.80 2.2 182ms
13 -7.926861681873 -14.01 -7.31 1.8 175ms
14 -7.926861681873 -15.05 -7.74 1.6 190ms
15 -7.926861681873 + -14.75 -8.24 2.1 198ms
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.921700867529 -0.69 5.5 274ms
2 -7.926135409685 -2.35 -1.22 1.0 131ms
3 -7.926837508611 -3.15 -2.37 2.0 156ms
4 -7.926864715280 -4.57 -3.02 3.0 211ms
5 -7.926865064897 -6.46 -3.45 2.0 169ms
6 -7.926865085607 -7.68 -3.99 1.8 153ms
7 -7.926865089650 -8.39 -4.09 1.9 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₂, 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₀")