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.110155136984 -0.83 0.80 7.0 25.8ms
2 -1.117160188830 -2.15 -1.81 0.80 1.0 20.1ms
3 -1.117249558762 -4.05 -2.56 0.80 1.0 20.3ms
4 -1.117250999529 -5.84 -3.20 0.80 1.0 20.4ms
5 -1.117251185694 -6.73 -4.17 0.80 1.0 20.8ms
6 -1.117251189913 -8.37 -4.62 0.80 2.0 23.5ms
7 -1.117251190131 -9.66 -5.28 0.80 1.0 21.7ms
8 -1.117251190143 -10.93 -6.12 0.80 1.0 64.8ms
9 -1.117251190143 -12.35 -6.85 0.80 2.0 24.5ms
10 -1.117251190143 -14.24 -7.45 0.80 1.0 22.4ms
11 -1.117251190143 + -15.05 -8.03 0.80 2.0 24.7ms
12 -1.117251190143 -15.05 -8.11 0.80 1.0 22.7ms
13 -1.117251190143 + -14.51 -8.83 0.80 1.0 24.5ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117251190143 -9.21 0.80 1.0 1.24s
2 -1.117251190143 + -14.81 -10.20 0.80 1.0 37.2ms
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 0 │ -1.117251190143 │ │ 0.0269318 │ 3.81s │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.117517998042 -2.46 0.80 1.0 13.7ms
2 -1.117520586052 -5.59 -3.49 0.80 1.0 19.9ms
3 -1.117520612985 -7.57 -4.23 0.80 1.0 20.2ms
4 -1.117520613449 -9.33 -5.14 0.80 1.0 20.5ms
5 -1.117520613464 -10.81 -5.92 0.80 1.0 36.0ms
6 -1.117520613465 -12.14 -6.42 0.80 1.0 21.6ms
7 -1.117520613465 -13.36 -7.39 0.80 1.0 21.8ms
8 -1.117520613465 + -15.05 -8.19 0.80 1.0 22.1ms
9 -1.117520613465 -15.35 -8.97 0.80 1.0 22.2ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118247264594 -1.69 0.80 1.0 13.5ms
2 -1.118339043933 -4.04 -2.70 0.80 1.0 20.2ms
3 -1.118340011112 -6.01 -3.44 0.80 1.0 20.4ms
4 -1.118340028553 -7.76 -4.36 0.80 1.0 20.5ms
5 -1.118340029137 -9.23 -5.12 0.80 1.0 20.8ms
6 -1.118340029164 -10.57 -5.62 0.80 1.0 26.5ms
7 -1.118340029166 -11.73 -6.68 0.80 1.0 21.8ms
8 -1.118340029166 -13.53 -7.58 0.80 1.0 22.4ms
9 -1.118340029166 + -15.05 -8.12 0.80 1.0 22.2ms
10 -1.118340029166 + -14.51 -8.55 0.80 1.0 22.0ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340029166 -9.05 0.80 1.0 13.5ms
2 -1.118340029166 + -14.75 -10.36 0.80 1.0 20.1ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118340029166 -11.04 0.80 1.0 13.4ms
2 -1.118340029166 + -14.95 -11.41 0.80 1.0 20.1ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 1 │ -1.118340029166 │ │ 0.00297055 │ 862ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118346721507 -3.09 0.80 1.0 13.5ms
2 -1.118346864758 -6.84 -4.11 0.80 1.0 20.2ms
3 -1.118346866235 -8.83 -4.86 0.80 1.0 20.3ms
4 -1.118346866261 -10.59 -5.78 0.80 1.0 20.6ms
5 -1.118346866262 -12.08 -6.55 0.80 1.0 20.9ms
6 -1.118346866262 -13.42 -7.04 0.80 1.0 30.9ms
7 -1.118346866262 + -15.05 -8.13 0.80 1.0 22.0ms
8 -1.118346866262 -14.33 -8.95 0.80 1.0 21.7ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118354724004 -2.61 0.80 1.0 13.4ms
2 -1.118356087056 -5.87 -3.62 0.80 1.0 26.6ms
3 -1.118356101131 -7.85 -4.37 0.80 1.0 20.5ms
4 -1.118356101377 -9.61 -5.29 0.80 1.0 20.5ms
5 -1.118356101385 -11.10 -6.06 0.80 1.0 20.8ms
6 -1.118356101386 -12.45 -6.55 0.80 1.0 21.1ms
7 -1.118356101386 -13.54 -7.65 0.80 1.0 21.6ms
8 -1.118356101386 + -14.35 -8.45 0.80 1.0 28.3ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356101386 -9.24 0.80 1.0 13.5ms
2 -1.118356101386 -14.37 -9.82 0.80 1.0 19.9ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356101386 -10.62 0.80 1.0 13.5ms
2 -1.118356101386 + -Inf -11.46 0.80 1.0 28.2ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 2 │ -1.118356101386 │ -4.79 │ 4.47268e-5 │ 606ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356102848 -4.94 0.80 1.0 13.6ms
2 -1.118356102877 -10.54 -5.96 0.80 1.0 19.9ms
3 -1.118356102877 -12.53 -6.70 0.80 1.0 20.1ms
4 -1.118356102877 -14.15 -7.62 0.80 1.0 29.5ms
5 -1.118356102877 -14.54 -8.39 0.80 1.0 21.4ms
6 -1.118356102877 + -14.54 -8.89 0.80 1.0 21.3ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356104771 -4.40 0.80 1.0 13.3ms
2 -1.118356105116 -9.46 -5.42 0.80 1.0 20.0ms
3 -1.118356105120 -11.45 -6.17 0.80 1.0 29.4ms
4 -1.118356105120 -13.21 -7.09 0.80 1.0 20.7ms
5 -1.118356105120 -14.15 -7.86 0.80 1.0 21.0ms
6 -1.118356105120 + -14.57 -8.35 0.80 1.0 21.2ms
7 -1.118356105120 + -14.75 -9.45 0.80 1.0 21.6ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356105120 -10.08 0.80 1.0 13.6ms
2 -1.118356105120 + -14.35 -11.02 0.80 1.0 20.2ms
n Energy log10(ΔE) log10(Δρ) α Diag Δtime
--- --------------- --------- --------- ---- ---- ------
1 -1.118356105120 -11.67 0.80 1.0 13.3ms
2 -1.118356105120 + -14.18 -12.44 0.80 1.0 20.0ms
┌─────┬─────────────────┬───────────┬─────────────┬────────┐
│ n │ Energy │ log10(ΔE) │ max(Force) │ Δtime │
├─────┼─────────────────┼───────────┼─────────────┼────────┤
│ 3 │ -1.118356105120 │ -8.43 │ 1.06082e-8 │ 534ms │
└─────┴─────────────────┴───────────┴─────────────┴────────┘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, -1.40094e-17, -1.09891e-17]u"a₀")
Atom(H, [ 1.44318, 2.13568e-17, -4.47802e-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.4863658662914383 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.