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.110844445339 -0.82 0.80 8.0 25.0ms
2 -1.117207081986 -2.20 -1.90 0.80 1.0 17.6ms
3 -1.117250716399 -4.36 -2.72 0.80 1.0 17.7ms
4 -1.117251292116 -6.24 -3.58 0.80 1.0 18.0ms
5 -1.117251309140 -7.77 -4.12 0.80 1.0 39.2ms
6 -1.117251310263 -8.95 -4.66 0.80 1.0 18.8ms
7 -1.117251310377 -9.94 -5.57 0.80 1.0 19.1ms
8 -1.117251310380 -11.48 -5.82 0.80 2.0 21.8ms
9 -1.117251310381 -12.15 -6.83 0.80 1.0 19.9ms
10 -1.117251310381 -14.31 -7.43 0.80 1.0 24.6ms
11 -1.117251310381 -15.35 -7.99 0.80 1.0 20.1ms
12 -1.117251310381 -14.70 -8.58 0.80 1.0 20.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117251310381 -9.13 0.80 1.0 20.8ms
2 -1.117251310381 -14.57 -10.32 0.80 1.0 17.4ms
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 0 │ -1.117251310381 │ │ 0.0269317 │ 850ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117518115445 -2.46 0.80 1.0 10.9ms
2 -1.117520704451 -5.59 -3.49 0.80 1.0 17.4ms
3 -1.117520731396 -7.57 -4.23 0.80 1.0 27.0ms
4 -1.117520731854 -9.34 -5.15 0.80 1.0 17.9ms
5 -1.117520731870 -10.81 -5.93 0.80 1.0 18.2ms
6 -1.117520731870 -12.14 -6.43 0.80 1.0 23.0ms
7 -1.117520731870 -13.36 -7.44 0.80 1.0 19.6ms
8 -1.117520731870 -15.35 -8.35 0.80 1.0 19.7ms
9 -1.117520731870 + -15.05 -8.81 0.80 1.0 19.7ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118247378334 -1.69 0.80 1.0 10.9ms
2 -1.118339157269 -4.04 -2.70 0.80 1.0 17.1ms
3 -1.118340124443 -6.01 -3.44 0.80 1.0 17.5ms
4 -1.118340141884 -7.76 -4.36 0.80 1.0 26.7ms
5 -1.118340142468 -9.23 -5.12 0.80 1.0 18.6ms
6 -1.118340142495 -10.57 -5.62 0.80 1.0 18.4ms
7 -1.118340142497 -11.73 -6.68 0.80 1.0 18.9ms
8 -1.118340142497 -13.52 -7.58 0.80 1.0 59.5ms
9 -1.118340142497 + -15.05 -8.12 0.80 1.0 25.5ms
10 -1.118340142497 + -15.35 -8.55 0.80 1.0 19.6ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340142497 -9.04 0.80 1.0 19.5ms
2 -1.118340142497 -14.75 -10.36 0.80 1.0 17.1ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 1 │ -1.118340142497 │ -2.96 │ 0.00297052 │ 858ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118346834725 -3.09 0.80 1.0 10.8ms
2 -1.118346977974 -6.84 -4.11 0.80 1.0 26.0ms
3 -1.118346979451 -8.83 -4.86 0.80 1.0 17.4ms
4 -1.118346979477 -10.59 -5.78 0.80 1.0 17.6ms
5 -1.118346979478 -12.09 -6.55 0.80 1.0 17.9ms
6 -1.118346979478 -13.40 -7.04 0.80 1.0 22.6ms
7 -1.118346979478 -14.88 -8.13 0.80 1.0 18.8ms
8 -1.118346979478 -15.05 -8.95 0.80 1.0 18.9ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118354837069 -2.61 0.80 1.0 10.8ms
2 -1.118356200096 -5.87 -3.62 0.80 1.0 26.3ms
3 -1.118356214170 -7.85 -4.37 0.80 1.0 17.5ms
4 -1.118356214417 -9.61 -5.29 0.80 1.0 17.8ms
5 -1.118356214425 -11.10 -6.06 0.80 1.0 18.0ms
6 -1.118356214425 -12.45 -6.55 0.80 1.0 23.1ms
7 -1.118356214425 -13.63 -7.65 0.80 1.0 18.8ms
8 -1.118356214425 + -15.05 -8.45 0.80 1.0 18.6ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356214425 -9.24 0.80 1.0 10.8ms
2 -1.118356214425 -15.65 -9.82 0.80 1.0 17.1ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 2 │ -1.118356214425 │ -4.79 │ 4.47259e-5 │ 746ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356215887 -4.94 0.80 1.0 10.8ms
2 -1.118356215917 -10.54 -5.96 0.80 1.0 27.0ms
3 -1.118356215917 -12.52 -6.70 0.80 1.0 17.4ms
4 -1.118356215917 -14.45 -7.62 0.80 1.0 17.6ms
5 -1.118356215917 -15.18 -8.39 0.80 1.0 18.1ms
6 -1.118356215917 + -14.81 -8.89 0.80 1.0 23.1ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356217810 -4.40 0.80 1.0 10.8ms
2 -1.118356218156 -9.46 -5.42 0.80 1.0 17.0ms
3 -1.118356218159 -11.45 -6.17 0.80 1.0 17.1ms
4 -1.118356218159 -13.21 -7.09 0.80 1.0 26.4ms
5 -1.118356218159 -14.51 -7.86 0.80 1.0 17.9ms
6 -1.118356218159 + -Inf -8.35 0.80 1.0 18.2ms
7 -1.118356218159 + -15.05 -9.45 0.80 1.0 23.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356218159 -10.08 0.80 1.0 10.8ms
2 -1.118356218159 -15.65 -11.02 0.80 1.0 17.0ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 3 │ -1.118356218159 │ -8.43 │ 1.06573e-8 │ 683ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
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.01718e-10, 4.31825e-10]u"a₀")
Atom(H, [ 1.44318, 2.00193e-10, 4.34112e-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.4863655838976793 a₀
Our results (1.486 Bohr) agrees with the equivalent tutorial from ABINIT.