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.111137571173 -0.83 0.80 9.0 25.4ms
2 -1.117207538046 -2.22 -1.99 0.80 1.0 16.6ms
3 -1.117250561322 -4.37 -2.75 0.80 1.0 16.6ms
4 -1.117251259679 -6.16 -3.65 0.80 1.0 17.0ms
5 -1.117251310267 -7.30 -4.46 0.80 2.0 44.6ms
6 -1.117251310380 -9.95 -5.79 0.80 1.0 17.9ms
7 -1.117251310381 -11.84 -6.35 0.80 1.0 18.3ms
8 -1.117251310381 -13.11 -7.30 0.80 1.0 18.4ms
9 -1.117251310381 + -15.65 -8.02 0.80 2.0 20.3ms
10 -1.117251310381 -15.05 -8.50 0.80 1.0 18.2ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117251310381 -9.05 0.80 1.0 10.4ms
2 -1.117251310381 -14.88 -10.30 0.80 1.0 16.4ms
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 0 │ -1.117251310381 │ │ 0.0269317 │ 611ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117518116472 -2.46 0.80 1.0 10.5ms
2 -1.117520704474 -5.59 -3.49 0.80 1.0 16.6ms
3 -1.117520731408 -7.57 -4.23 0.80 1.0 16.7ms
4 -1.117520731872 -9.33 -5.14 0.80 1.0 17.0ms
5 -1.117520731888 -10.81 -5.92 0.80 1.0 17.6ms
6 -1.117520731888 -12.14 -6.42 0.80 1.0 17.8ms
7 -1.117520731888 -13.30 -7.39 0.80 1.0 18.3ms
8 -1.117520731888 + -14.65 -8.19 0.80 1.0 18.9ms
9 -1.117520731888 -15.35 -8.93 0.80 1.0 24.8ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118247378322 -1.69 0.80 1.0 10.4ms
2 -1.118339157060 -4.04 -2.70 0.80 1.0 16.5ms
3 -1.118340124233 -6.01 -3.44 0.80 1.0 16.7ms
4 -1.118340141673 -7.76 -4.36 0.80 1.0 17.2ms
5 -1.118340142258 -9.23 -5.12 0.80 1.0 17.4ms
6 -1.118340142285 -10.57 -5.62 0.80 1.0 24.5ms
7 -1.118340142286 -11.73 -6.68 0.80 1.0 18.3ms
8 -1.118340142286 -13.56 -7.58 0.80 1.0 18.8ms
9 -1.118340142286 -14.70 -8.12 0.80 1.0 18.9ms
10 -1.118340142286 + -15.05 -8.55 0.80 1.0 18.9ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340142286 -9.05 0.80 1.0 10.6ms
2 -1.118340142286 + -14.88 -10.36 0.80 1.0 16.6ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 1 │ -1.118340142286 │ -2.96 │ 0.00297054 │ 609ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118346834596 -3.09 0.80 1.0 10.4ms
2 -1.118346977847 -6.84 -4.11 0.80 1.0 23.5ms
3 -1.118346979324 -8.83 -4.86 0.80 1.0 16.7ms
4 -1.118346979350 -10.59 -5.78 0.80 1.0 17.1ms
5 -1.118346979350 -12.08 -6.55 0.80 1.0 17.7ms
6 -1.118346979350 -13.42 -7.04 0.80 1.0 22.8ms
7 -1.118346979350 -14.57 -8.13 0.80 1.0 18.4ms
8 -1.118346979350 + -14.61 -8.95 0.80 1.0 18.2ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118354837050 -2.61 0.80 1.0 10.4ms
2 -1.118356200096 -5.87 -3.62 0.80 1.0 16.6ms
3 -1.118356214170 -7.85 -4.37 0.80 1.0 16.7ms
4 -1.118356214417 -9.61 -5.29 0.80 1.0 23.3ms
5 -1.118356214425 -11.10 -6.06 0.80 1.0 17.6ms
6 -1.118356214425 -12.45 -6.55 0.80 1.0 17.7ms
7 -1.118356214425 -13.65 -7.65 0.80 1.0 18.1ms
8 -1.118356214425 -14.81 -8.45 0.80 1.0 23.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356214425 -9.24 0.80 1.0 10.4ms
2 -1.118356214425 -14.75 -9.82 0.80 1.0 16.6ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 2 │ -1.118356214425 │ -4.79 │ 4.47265e-5 │ 506ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356215887 -4.94 0.80 1.0 10.6ms
2 -1.118356215917 -10.54 -5.96 0.80 1.0 16.6ms
3 -1.118356215917 -12.52 -6.70 0.80 1.0 16.7ms
4 -1.118356215917 -14.37 -7.62 0.80 1.0 23.7ms
5 -1.118356215917 -15.18 -8.39 0.80 1.0 17.7ms
6 -1.118356215917 -15.35 -8.89 0.80 1.0 17.8ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356217810 -4.40 0.80 1.0 10.4ms
2 -1.118356218156 -9.46 -5.42 0.80 1.0 16.7ms
3 -1.118356218159 -11.45 -6.17 0.80 1.0 16.9ms
4 -1.118356218159 -13.19 -7.09 0.80 1.0 17.0ms
5 -1.118356218159 -14.65 -7.86 0.80 1.0 22.5ms
6 -1.118356218159 + -14.88 -8.35 0.80 1.0 17.7ms
7 -1.118356218159 + -14.51 -9.45 0.80 1.0 18.3ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356218159 -10.08 0.80 1.0 10.5ms
2 -1.118356218159 + -15.35 -11.02 0.80 1.0 16.6ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 3 │ -1.118356218159 │ -8.43 │ 1.06793e-8 │ 452ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
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.0431828, -2.8462e-10, -7.23186e-10]u"a₀")
Atom(H, [ 1.44318, -2.82715e-10, -6.97579e-10]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.4863655838152152 a₀
Our results (1.486 Bohr) agrees with the equivalent tutorial from ABINIT.