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.921748966874 -0.69 5.6 194ms
2 -7.926134696242 -2.36 -1.22 1.0 163ms
3 -7.926833530336 -3.16 -2.37 2.1 162ms
4 -7.926861275507 -4.56 -3.00 2.9 216ms
5 -7.926861644809 -6.43 -3.39 2.4 167ms
6 -7.926861670702 -7.59 -3.85 1.4 144ms
7 -7.926861679426 -8.06 -4.16 2.0 157ms
8 -7.926861681768 -8.63 -4.92 1.4 148ms
9 -7.926861681855 -10.06 -5.16 2.1 164ms
10 -7.926861681869 -10.84 -5.65 1.5 169ms
11 -7.926861681872 -11.51 -6.28 1.4 147ms
12 -7.926861681873 -12.77 -7.40 2.2 171ms
┌ Warning: Eigensolver not converged
│ n_iter =
│ 8-element Vector{Int64}:
│ 4
│ 3
│ 3
│ 2
│ 3
│ 3
│ 3
│ 3
└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/scf/self_consistent_field.jl:76
13 -7.926861681873 -14.10 -7.34 3.0 219ms
14 -7.926861681873 -14.75 -8.28 1.0 143ms
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.921741530760 -0.69 5.5 192ms
2 -7.926134999567 -2.36 -1.22 1.0 136ms
3 -7.926834043590 -3.16 -2.37 2.0 159ms
4 -7.926860878430 -4.57 -3.02 3.1 195ms
5 -7.926861660265 -6.11 -3.43 2.1 168ms
6 -7.926861674987 -7.83 -3.98 1.5 154ms
7 -7.926861681387 -8.19 -4.72 1.9 160ms
8 -7.926861681849 -9.33 -5.45 2.4 212ms
9 -7.926861681871 -10.68 -5.49 2.4 188ms
10 -7.926861681872 -11.73 -6.31 1.0 147ms
11 -7.926861681873 -12.81 -7.27 2.1 166ms
┌ Warning: Eigensolver not converged
│ n_iter =
│ 8-element Vector{Int64}:
│ 3
│ 3
│ 3
│ 2
│ 3
│ 4
│ 2
│ 3
└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/scf/self_consistent_field.jl:76
12 -7.926861681873 -14.10 -7.77 2.9 196ms
┌ Warning: Eigensolver not converged
│ n_iter =
│ 8-element Vector{Int64}:
│ 2
│ 3
│ 2
│ 1
│ 3
│ 1
│ 2
│ 2
└ @ DFTK ~/work/DFTK.jl/DFTK.jl/src/scf/self_consistent_field.jl:76
13 -7.926861681873 + -14.57 -8.27 2.0 166ms
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.921700727892 -0.69 5.6 220ms
2 -7.926137916871 -2.35 -1.22 1.0 139ms
3 -7.926837088806 -3.16 -2.37 2.0 147ms
4 -7.926864714687 -4.56 -2.99 3.0 202ms
5 -7.926865050668 -6.47 -3.35 2.2 164ms
6 -7.926865079490 -7.54 -3.78 1.4 140ms
7 -7.926865090979 -7.94 -4.24 1.5 164ms
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₀")