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
Feeding an AtomsBase AbstractSystem to DFTK
In this example we construct a silicon system using the ase.build.bulk
routine from the atomistic simulation environment (ASE), which is exposed by ASEconvert as an AtomsBase AbstractSystem
.
# Construct bulk system and convert to an AbstractSystem
using ASEconvert
system_ase = ase.build.bulk("Si")
system = pyconvert(AbstractSystem, system_ase)
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
To use an AbstractSystem in DFTK, we attach pseudopotentials, construct a DFT model, discretise and solve:
system = attach_psp(system; Si="hgh/lda/si-q4")
model = model_LDA(system; 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.921680683459 -0.69 5.9
2 -7.926164798114 -2.35 -1.22 1.0 223ms
3 -7.926836807179 -3.17 -2.37 1.8 252ms
4 -7.926861511277 -4.61 -3.01 2.9 311ms
5 -7.926861631371 -6.92 -3.33 1.8 219ms
6 -7.926861666275 -7.46 -3.71 1.6 218ms
7 -7.926861680629 -7.84 -4.40 1.2 221ms
8 -7.926861681816 -8.93 -4.97 2.1 280ms
9 -7.926861681861 -10.35 -5.25 1.8 230ms
10 -7.926861681872 -10.96 -6.03 1.6 224ms
11 -7.926861681873 -12.15 -7.03 2.6 285ms
12 -7.926861681873 -13.73 -7.38 2.8 262ms
13 -7.926861681873 -15.05 -8.04 1.8 230ms
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_LDA(system; 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.921702456820 -0.69 6.0
2 -7.926161886436 -2.35 -1.22 1.0 226ms
3 -7.926837418102 -3.17 -2.37 1.9 256ms
4 -7.926861510761 -4.62 -3.02 2.6 296ms
5 -7.926861637211 -6.90 -3.37 1.6 220ms
6 -7.926861667713 -7.52 -3.75 1.6 222ms
7 -7.926861680542 -7.89 -4.37 1.4 222ms
8 -7.926861681795 -8.90 -4.91 2.2 277ms
9 -7.926861681861 -10.18 -5.27 1.6 234ms
10 -7.926861681871 -10.98 -5.85 1.6 220ms
11 -7.926861681873 -11.89 -7.04 1.9 260ms
12 -7.926861681873 -13.35 -7.44 3.4 285ms
13 -7.926861681873 + -14.57 -8.01 1.4 214ms
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_LDA(system; 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.921658530002 -0.69 5.9
2 -7.926167369828 -2.35 -1.22 1.0 230ms
3 -7.926839762432 -3.17 -2.37 1.8 259ms
4 -7.926864918720 -4.60 -3.00 2.9 327ms
5 -7.926865033180 -6.94 -3.30 1.6 233ms
6 -7.926865076159 -7.37 -3.69 1.5 224ms
7 -7.926865091813 -7.81 -4.43 1.1 210ms
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