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.110642617738 -0.82 0.80 9.0 29.0ms
2 -1.117202471555 -2.18 -1.88 0.80 1.0 19.5ms
3 -1.117250460642 -4.32 -2.71 0.80 1.0 19.5ms
4 -1.117251161071 -6.15 -3.50 0.80 1.0 20.1ms
5 -1.117251188387 -7.56 -4.06 0.80 1.0 134ms
6 -1.117251189968 -8.80 -4.56 0.80 1.0 21.0ms
7 -1.117251190136 -9.77 -5.40 0.80 1.0 21.4ms
8 -1.117251190143 -11.16 -6.31 0.80 1.0 21.9ms
9 -1.117251190143 -13.05 -7.06 0.80 1.0 21.7ms
10 -1.117251190143 -14.24 -7.79 0.80 2.0 24.0ms
11 -1.117251190143 -14.81 -8.91 0.80 1.0 22.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117251190143 -9.06 0.80 1.0 787ms
2 -1.117251190143 -14.88 -10.12 0.80 1.0 62.1ms
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 0 │ -1.117251190143 │ │ 0.0269318 │ 3.23s │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117517997923 -2.46 0.80 1.0 13.1ms
2 -1.117520586054 -5.59 -3.49 0.80 1.0 19.5ms
3 -1.117520612988 -7.57 -4.23 0.80 1.0 19.6ms
4 -1.117520613449 -9.34 -5.14 0.80 1.0 19.9ms
5 -1.117520613464 -10.81 -5.93 0.80 1.0 20.4ms
6 -1.117520613465 -12.16 -6.43 0.80 1.0 20.9ms
7 -1.117520613465 -13.34 -7.43 0.80 1.0 21.0ms
8 -1.117520613465 + -15.18 -8.28 0.80 1.0 35.3ms
9 -1.117520613465 -15.65 -8.91 0.80 1.0 22.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118247264597 -1.69 0.80 1.0 13.0ms
2 -1.118339043932 -4.04 -2.70 0.80 1.0 19.4ms
3 -1.118340011111 -6.01 -3.44 0.80 1.0 19.7ms
4 -1.118340028552 -7.76 -4.36 0.80 1.0 20.0ms
5 -1.118340029136 -9.23 -5.12 0.80 1.0 20.4ms
6 -1.118340029163 -10.57 -5.62 0.80 1.0 31.8ms
7 -1.118340029165 -11.73 -6.68 0.80 1.0 21.3ms
8 -1.118340029165 -13.57 -7.58 0.80 1.0 21.6ms
9 -1.118340029165 + -14.95 -8.12 0.80 1.0 21.8ms
10 -1.118340029165 -14.75 -8.55 0.80 1.0 21.7ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340029165 -9.05 0.80 1.0 20.8ms
2 -1.118340029165 + -14.75 -10.36 0.80 1.0 19.8ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340029165 -11.05 0.80 1.0 13.1ms
2 -1.118340029165 -14.40 -11.42 0.80 1.0 19.4ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 1 │ -1.118340029165 │ │ 0.00297055 │ 824ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118346721506 -3.09 0.80 1.0 21.4ms
2 -1.118346864758 -6.84 -4.11 0.80 1.0 19.4ms
3 -1.118346866235 -8.83 -4.86 0.80 1.0 19.9ms
4 -1.118346866261 -10.59 -5.78 0.80 1.0 20.1ms
5 -1.118346866261 -12.08 -6.55 0.80 1.0 20.6ms
6 -1.118346866261 -13.43 -7.04 0.80 1.0 26.9ms
7 -1.118346866261 -14.95 -8.13 0.80 1.0 21.3ms
8 -1.118346866261 -14.75 -8.95 0.80 1.0 21.5ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118354724004 -2.61 0.80 1.0 13.0ms
2 -1.118356087056 -5.87 -3.62 0.80 1.0 25.9ms
3 -1.118356101131 -7.85 -4.37 0.80 1.0 19.7ms
4 -1.118356101377 -9.61 -5.29 0.80 1.0 19.9ms
5 -1.118356101385 -11.10 -6.06 0.80 1.0 20.2ms
6 -1.118356101386 -12.45 -6.55 0.80 1.0 20.7ms
7 -1.118356101386 -13.60 -7.65 0.80 1.0 26.0ms
8 -1.118356101386 -14.54 -8.45 0.80 1.0 21.3ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356101386 -9.24 0.80 1.0 13.0ms
2 -1.118356101386 + -14.65 -9.82 0.80 1.0 19.4ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356101386 -10.62 0.80 1.0 13.2ms
2 -1.118356101386 -14.65 -11.46 0.80 1.0 19.4ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 2 │ -1.118356101386 │ -4.79 │ 4.47268e-5 │ 595ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356102848 -4.94 0.80 1.0 13.0ms
2 -1.118356102877 -10.54 -5.96 0.80 1.0 19.5ms
3 -1.118356102877 -12.52 -6.70 0.80 1.0 27.7ms
4 -1.118356102877 -14.07 -7.62 0.80 1.0 20.0ms
5 -1.118356102877 -14.95 -8.40 0.80 1.0 20.4ms
6 -1.118356102877 + -14.16 -8.89 0.80 1.0 20.7ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356104771 -4.40 0.80 1.0 21.1ms
2 -1.118356105116 -9.46 -5.42 0.80 1.0 19.8ms
3 -1.118356105120 -11.45 -6.17 0.80 1.0 19.7ms
4 -1.118356105120 -13.21 -7.09 0.80 1.0 19.9ms
5 -1.118356105120 -14.70 -7.86 0.80 1.0 20.5ms
6 -1.118356105120 + -14.61 -8.35 0.80 1.0 25.2ms
7 -1.118356105120 -14.57 -9.45 0.80 1.0 21.5ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356105120 -10.08 0.80 1.0 13.1ms
2 -1.118356105120 + -14.75 -11.02 0.80 1.0 19.5ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356105120 -11.67 0.80 1.0 20.0ms
2 -1.118356105120 -14.33 -12.44 0.80 1.0 19.5ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 3 │ -1.118356105120 │ -8.43 │ 1.19554e-8 │ 527ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘Structure after optimisation (note that the atom has wrapped around)
results.systemFlexibleSystem(H₂, periodicity = TTT):
cell_vectors : [ 10 0 0;
0 10 0;
0 0 10]u"a₀"
Atom(H, [-0.0431829, -7.21732e-18, -1.584e-17]u"a₀")
Atom(H, [ 1.44318, -1.04581e-18, -4.4952e-18]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.486365861260072 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.