Geometry optimization
We use DFTK and the GeometryOptimization package to find the minimal-energy bond length of the $H_2$ molecule. First we set up an appropriate DFTKCalculator
(see AtomsCalculators integration), for which we use the LDA model just like in the Tutorial for silicon in combination with a pseudodojo pseudopotential (see Pseudopotentials).
using DFTK
using PseudoPotentialData
pseudopotentials = PseudoFamily("dojo.nc.sr.pbe.v0_4_1.standard.upf")
calc = DFTKCalculator(;
model_kwargs = (; functionals=LDA(), pseudopotentials), # model_DFT keyword arguments
basis_kwargs = (; kgrid=[1, 1, 1], Ecut=10) # PlaneWaveBasis keyword arguments
)
DFTKCalculator(functionals=Xc(lda_x, lda_c_pw), pseudopotentials=PseudoFamily("dojo.nc.sr.pbe.v0_4_1.standard.upf"), Ecut=10, kgrid=[1, 1, 1])
Next we set up an initial hydrogen molecule within a box of vacuum. We use the parameters of the equivalent tutorial from ABINIT and DFTK's AtomsBase integration to setup the hydrogen molecule. We employ a simulation box of 10 bohr times 10 bohr times 10 bohr and a pseudodojo pseudopotential.
using LinearAlgebra
using Unitful
using UnitfulAtomic
r0 = 1.4 # Initial bond length in Bohr
a = 10.0 # Box size in Bohr
cell_vectors = [[a, 0, 0]u"bohr", [0, a, 0]u"bohr", [0, 0, a]u"bohr"]
h2_crude = periodic_system([:H => [0, 0, 0.]u"bohr",
:H => [r0, 0, 0]u"bohr"],
cell_vectors)
FlexibleSystem(H₂, periodicity = TTT):
cell_vectors : [ 10 0 0;
0 10 0;
0 0 10]u"a₀"
Atom(H, [ 0, 0, 0]u"a₀")
Atom(H, [ 1.4, 0, 0]u"a₀")
Finally we call minimize_energy!
to start the geometry optimisation. We use verbosity=2
to get some insight into the minimisation. With verbosity=1
only a summarising table would be printed and with verbosity=0
(default) the minimisation would be quiet.
using GeometryOptimization
results = minimize_energy!(h2_crude, calc; tol_forces=2e-6, verbosity=2)
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.110543828153 -0.82 0.80 7.0 22.5ms
2 -1.117203592624 -2.18 -1.86 0.80 1.0 41.6ms
3 -1.117250624712 -4.33 -2.71 0.80 1.0 17.5ms
4 -1.117251177947 -6.26 -3.61 0.80 1.0 17.4ms
5 -1.117251188908 -7.96 -4.02 0.80 1.0 17.7ms
6 -1.117251190117 -8.92 -5.06 0.80 1.0 18.1ms
7 -1.117251190142 -10.61 -5.83 0.80 1.0 32.4ms
8 -1.117251190143 -11.91 -6.33 0.80 2.0 21.1ms
9 -1.117251190143 + -13.37 -6.23 0.80 1.0 19.7ms
10 -1.117251190143 -13.11 -7.44 0.80 1.0 20.1ms
11 -1.117251190143 -14.33 -7.74 0.80 2.0 21.9ms
12 -1.117251190143 + -14.81 -8.22 0.80 1.0 23.0ms
13 -1.117251190143 + -15.05 -8.76 0.80 1.0 20.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117251190143 -9.13 0.80 1.0 94.1ms
2 -1.117251190143 -14.61 -10.23 0.80 1.0 24.4ms
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 0 │ -1.117251190143 │ │ 0.0269318 │ 1.28s │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117517998019 -2.46 0.80 1.0 10.5ms
2 -1.117520586052 -5.59 -3.49 0.80 1.0 17.1ms
3 -1.117520612985 -7.57 -4.23 0.80 1.0 23.0ms
4 -1.117520613449 -9.33 -5.14 0.80 1.0 18.2ms
5 -1.117520613464 -10.81 -5.92 0.80 1.0 18.1ms
6 -1.117520613465 -12.14 -6.42 0.80 1.0 18.3ms
7 -1.117520613465 -13.33 -7.39 0.80 1.0 22.2ms
8 -1.117520613465 + -14.88 -8.19 0.80 1.0 19.4ms
9 -1.117520613465 -14.75 -8.96 0.80 1.0 19.3ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118247264594 -1.69 0.80 1.0 10.7ms
2 -1.118339043935 -4.04 -2.70 0.80 1.0 16.8ms
3 -1.118340011115 -6.01 -3.44 0.80 1.0 17.0ms
4 -1.118340028555 -7.76 -4.36 0.80 1.0 669ms
5 -1.118340029140 -9.23 -5.12 0.80 1.0 18.4ms
6 -1.118340029167 -10.57 -5.62 0.80 1.0 18.4ms
7 -1.118340029168 -11.73 -6.68 0.80 1.0 19.2ms
8 -1.118340029168 -13.57 -7.58 0.80 1.0 19.7ms
9 -1.118340029168 + -14.88 -8.12 0.80 1.0 19.6ms
10 -1.118340029168 -14.45 -8.55 0.80 1.0 19.6ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340029168 -9.05 0.80 1.0 11.0ms
2 -1.118340029168 + -14.70 -10.36 0.80 1.0 17.5ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340029168 -11.04 0.80 1.0 10.9ms
2 -1.118340029168 -14.88 -11.41 0.80 1.0 17.4ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 1 │ -1.118340029168 │ │ 0.00297055 │ 1.40s │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118346721508 -3.09 0.80 1.0 10.8ms
2 -1.118346864760 -6.84 -4.11 0.80 1.0 17.6ms
3 -1.118346866237 -8.83 -4.86 0.80 1.0 17.7ms
4 -1.118346866262 -10.59 -5.78 0.80 1.0 18.0ms
5 -1.118346866263 -12.09 -6.55 0.80 1.0 18.5ms
6 -1.118346866263 -13.39 -7.04 0.80 1.0 18.7ms
7 -1.118346866263 -14.51 -8.13 0.80 1.0 19.4ms
8 -1.118346866263 + -14.95 -8.95 0.80 1.0 19.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118354724004 -2.61 0.80 1.0 10.6ms
2 -1.118356087056 -5.87 -3.62 0.80 1.0 24.8ms
3 -1.118356101131 -7.85 -4.37 0.80 1.0 17.4ms
4 -1.118356101377 -9.61 -5.29 0.80 1.0 18.0ms
5 -1.118356101385 -11.10 -6.06 0.80 1.0 18.4ms
6 -1.118356101386 -12.45 -6.55 0.80 1.0 18.5ms
7 -1.118356101386 -13.63 -7.65 0.80 1.0 19.1ms
8 -1.118356101386 -14.45 -8.45 0.80 1.0 18.9ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356101386 -9.24 0.80 1.0 10.7ms
2 -1.118356101386 -14.42 -9.82 0.80 1.0 17.2ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356101386 -10.62 0.80 1.0 10.5ms
2 -1.118356101386 + -14.70 -11.46 0.80 1.0 16.9ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 2 │ -1.118356101386 │ -4.79 │ 4.47268e-5 │ 577ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356102848 -4.94 0.80 1.0 17.7ms
2 -1.118356102877 -10.54 -5.96 0.80 1.0 17.3ms
3 -1.118356102877 -12.52 -6.70 0.80 1.0 17.2ms
4 -1.118356102877 -14.22 -7.62 0.80 1.0 17.6ms
5 -1.118356102877 + -14.95 -8.40 0.80 1.0 17.9ms
6 -1.118356102877 + -14.75 -8.89 0.80 1.0 18.1ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356104771 -4.40 0.80 1.0 10.6ms
2 -1.118356105116 -9.46 -5.42 0.80 1.0 16.9ms
3 -1.118356105120 -11.45 -6.17 0.80 1.0 17.2ms
4 -1.118356105120 -13.19 -7.09 0.80 1.0 17.4ms
5 -1.118356105120 + -14.88 -7.86 0.80 1.0 22.8ms
6 -1.118356105120 + -14.88 -8.35 0.80 1.0 18.5ms
7 -1.118356105120 -14.88 -9.45 0.80 1.0 18.8ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356105120 -10.08 0.80 1.0 10.7ms
2 -1.118356105120 -14.88 -11.02 0.80 1.0 17.1ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356105120 -11.67 0.80 1.0 10.5ms
2 -1.118356105120 + -15.05 -12.44 0.80 1.0 16.8ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 3 │ -1.118356105120 │ -8.43 │ 1.13986e-8 │ 525ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
Structure after optimisation (note that the atom has wrapped around)
results.system
FlexibleSystem(H₂, periodicity = TTT):
cell_vectors : [ 10 0 0;
0 10 0;
0 0 10]u"a₀"
Atom(H, [-0.0431829, -1.85996e-11, -2.09476e-11]u"a₀")
Atom(H, [ 1.44318, -1.83106e-11, -2.09361e-11]u"a₀")
Compute final bond length:
rmin = norm(position(results.system[1]) - position(results.system[2]))
println("Optimal bond length: ", rmin)
Optimal bond length: 1.4863658633393886 a₀
Our results (1.486 Bohr) agrees with the equivalent tutorial from ABINIT.
Recent versions of GeometryOptimization support cell shape optimisations as well by passing variablecell=true
to minimize_energy!
. See the GeometryOptimization documentation for an example.