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.110565300443 -0.82 0.80 8.0 24.1ms
2 -1.117199497860 -2.18 -1.86 0.80 1.0 39.1ms
3 -1.117250617723 -4.29 -2.70 0.80 1.0 16.7ms
4 -1.117251283418 -6.18 -3.52 0.80 1.0 16.9ms
5 -1.117251308505 -7.60 -4.02 0.80 1.0 17.3ms
6 -1.117251310274 -8.75 -4.80 0.80 1.0 17.7ms
7 -1.117251310377 -9.99 -5.63 0.80 1.0 32.9ms
8 -1.117251310381 -11.40 -6.36 0.80 2.0 20.7ms
9 -1.117251310381 -13.07 -7.16 0.80 1.0 18.7ms
10 -1.117251310381 -14.54 -7.55 0.80 2.0 20.8ms
11 -1.117251310381 -14.88 -8.44 0.80 1.0 18.9ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117251310381 -8.79 0.80 1.0 10.6ms
2 -1.117251310381 + -14.75 -9.75 0.80 1.0 16.7ms
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 0 │ -1.117251310381 │ │ 0.0269317 │ 792ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117518116330 -2.46 0.80 1.0 10.5ms
2 -1.117520704450 -5.59 -3.49 0.80 1.0 31.2ms
3 -1.117520731384 -7.57 -4.23 0.80 1.0 16.7ms
4 -1.117520731847 -9.33 -5.14 0.80 1.0 17.0ms
5 -1.117520731863 -10.81 -5.92 0.80 1.0 17.3ms
6 -1.117520731863 -12.14 -6.42 0.80 1.0 22.1ms
7 -1.117520731864 -13.37 -7.40 0.80 1.0 18.2ms
8 -1.117520731864 + -15.35 -8.20 0.80 1.0 18.8ms
9 -1.117520731864 -15.05 -8.92 0.80 1.0 44.5ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118247378331 -1.69 0.80 1.0 10.5ms
2 -1.118339157348 -4.04 -2.70 0.80 1.0 16.5ms
3 -1.118340124524 -6.01 -3.44 0.80 1.0 16.6ms
4 -1.118340141964 -7.76 -4.36 0.80 1.0 26.0ms
5 -1.118340142549 -9.23 -5.12 0.80 1.0 17.3ms
6 -1.118340142576 -10.57 -5.62 0.80 1.0 17.6ms
7 -1.118340142578 -11.73 -6.68 0.80 1.0 18.1ms
8 -1.118340142578 -13.53 -7.58 0.80 1.0 22.9ms
9 -1.118340142578 + -14.70 -8.12 0.80 1.0 18.6ms
10 -1.118340142578 -14.75 -8.55 0.80 1.0 18.6ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340142578 -9.05 0.80 1.0 10.5ms
2 -1.118340142578 -15.05 -10.36 0.80 1.0 16.4ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 1 │ -1.118340142578 │ -2.96 │ 0.00297051 │ 833ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118346834774 -3.09 0.80 1.0 10.4ms
2 -1.118346978023 -6.84 -4.11 0.80 1.0 30.7ms
3 -1.118346979500 -8.83 -4.86 0.80 1.0 16.6ms
4 -1.118346979525 -10.59 -5.78 0.80 1.0 17.2ms
5 -1.118346979526 -12.08 -6.55 0.80 1.0 17.3ms
6 -1.118346979526 -13.44 -7.04 0.80 1.0 21.8ms
7 -1.118346979526 -14.21 -8.13 0.80 1.0 18.2ms
8 -1.118346979526 + -14.88 -8.95 0.80 1.0 18.1ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118354837076 -2.61 0.80 1.0 18.9ms
2 -1.118356200096 -5.87 -3.62 0.80 1.0 16.4ms
3 -1.118356214171 -7.85 -4.37 0.80 1.0 16.6ms
4 -1.118356214417 -9.61 -5.29 0.80 1.0 16.9ms
5 -1.118356214425 -11.10 -6.06 0.80 1.0 21.4ms
6 -1.118356214425 -12.45 -6.55 0.80 1.0 17.6ms
7 -1.118356214426 -13.66 -7.65 0.80 1.0 18.8ms
8 -1.118356214426 -14.65 -8.45 0.80 1.0 18.4ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356214426 -9.24 0.80 1.0 10.6ms
2 -1.118356214426 -15.05 -9.82 0.80 1.0 16.5ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 2 │ -1.118356214426 │ -4.79 │ 4.47257e-5 │ 749ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356215887 -4.94 0.80 1.0 10.4ms
2 -1.118356215917 -10.54 -5.96 0.80 1.0 29.8ms
3 -1.118356215917 -12.52 -6.70 0.80 1.0 16.7ms
4 -1.118356215917 -14.33 -7.62 0.80 1.0 17.0ms
5 -1.118356215917 -15.05 -8.39 0.80 1.0 17.3ms
6 -1.118356215917 + -Inf -8.89 0.80 1.0 21.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.4ms
3 -1.118356218159 -11.45 -6.17 0.80 1.0 25.4ms
4 -1.118356218159 -13.17 -7.09 0.80 1.0 17.6ms
5 -1.118356218159 + -14.37 -7.86 0.80 1.0 17.5ms
6 -1.118356218159 -14.65 -8.35 0.80 1.0 17.7ms
7 -1.118356218159 -14.75 -9.45 0.80 1.0 22.2ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356218159 -10.08 0.80 1.0 10.6ms
2 -1.118356218159 -14.75 -11.02 0.80 1.0 16.4ms
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│ 3 │ -1.118356218159 │ -8.43 │ 1.07179e-8 │ 679ms │
└─────┴─────────────────┴───────────┴────────────┴────────┘
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, -1.20765e-09, -5.6322e-10]u"a₀")
Atom(H, [ 1.44318, -1.22061e-09, -5.85087e-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.4863655836710123 a₀
Our results (1.486 Bohr) agrees with the equivalent tutorial from ABINIT.