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.

using DFTK

calc = DFTKCalculator(;
    model_kwargs = (; functionals=LDA()),        # model_DFT keyword arguments
    basis_kwargs = (; kgrid=[1, 1, 1], Ecut=10)  # PlaneWaveBasis keyword arguments
)
DFTKCalculator{Nothing}(DFTK.DFTKParameters((functionals = Xc(lda_x, lda_c_pw),), (kgrid = [1, 1, 1], Ecut = 10), (callback = identity,)), nothing, Base.RefValue{Int64}(0), Base.RefValue{Int64}(0), true)

Next we set up an initial hydrogen molecule within a box of vacuum. We use the parameters of the equivalent tutorial from ABINIT, that is a simulation box of 10 bohr times 10 bohr times 10 bohr and a pseudodojo pseudopotential.

using LinearAlgebra
using LazyArtifacts

r0 = 1.4   # Initial bond length in Bohr
a  = 10.0  # Box size in Bohr

lattice = a * I(3)
H = ElementPsp(:H; psp=load_psp(artifact"pd_nc_sr_pbe_standard_0.4.1_upf/H.upf"));
atoms = [H, H]
positions = [zeros(3), lattice \ [r0, 0., 0.]]

h2_crude = periodic_system(lattice, atoms, positions)
FlexibleSystem(H₂, periodic = TTT):
    bounding_box      : [      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₀")

       .------------.  
      /|            |  
     / |            |  
    *  |            |  
    |  |            |  
    |  |            |  
    |  .------------.  
    | /            /   
    |/            /    
    H-H----------*     

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.110556589431                   -0.82   0.80    8.0   27.0ms
  2   -1.117202297379       -2.18       -1.86   0.80    1.0   17.3ms
  3   -1.117251826676       -4.31       -2.71   0.80    1.0   17.0ms
  4   -1.117252479208       -6.19       -3.55   0.80    1.0   17.5ms
  5   -1.117252499101       -7.70       -3.99   0.80    1.0   17.9ms
  6   -1.117252500639       -8.81       -4.54   0.80    1.0   18.3ms
  7   -1.117252500840       -9.70       -5.40   0.80    1.0   35.0ms
  8   -1.117252500846      -11.19       -6.24   0.80    2.0   21.7ms
  9   -1.117252500846      -12.90       -6.82   0.80    1.0   19.3ms
 10   -1.117252500846      -14.51       -6.98   0.80    2.0   21.8ms
 11   -1.117252500846      -14.37       -8.36   0.80    1.0   19.6ms
 12   -1.117252500846      -15.05       -8.29   0.80    2.0   21.9ms
 13   -1.117252500846   +  -15.35       -8.74   0.80    1.0   20.1ms
n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.117252500846                   -9.26   0.80    1.0   79.8ms
  2   -1.117252500846      -14.75      -10.31   0.80    1.0   16.9ms

Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n           Energy  log10(ΔE)  max(Force)   Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.117252500846 │           │  0.0269316 │  2.27s │
└─────┴─────────────────┴───────────┴────────────┴────────┘

n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.117519305142                   -2.46   0.80    1.0   10.9ms
  2   -1.117521893101       -5.59       -3.49   0.80    1.0   17.1ms
  3   -1.117521920034       -7.57       -4.23   0.80    1.0   17.0ms
  4   -1.117521920498       -9.33       -5.14   0.80    1.0   17.8ms
  5   -1.117521920514      -10.81       -5.92   0.80    1.0   18.1ms
  6   -1.117521920514      -12.13       -6.42   0.80    1.0   18.0ms
  7   -1.117521920514      -13.35       -7.39   0.80    1.0   19.2ms
  8   -1.117521920514   +  -14.70       -8.18   0.80    1.0   19.3ms
  9   -1.117521920514      -14.61       -8.97   0.80    1.0   19.5ms
n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118248561944                   -1.69   0.80    1.0   11.0ms
  2   -1.118340340104       -4.04       -2.70   0.80    1.0   37.6ms
  3   -1.118341307269       -6.01       -3.44   0.80    1.0   17.4ms
  4   -1.118341324709       -7.76       -4.36   0.80    1.0   17.7ms
  5   -1.118341325294       -9.23       -5.12   0.80    1.0   18.0ms
  6   -1.118341325321      -10.57       -5.62   0.80    1.0   18.4ms
  7   -1.118341325323      -11.73       -6.68   0.80    1.0   18.9ms
  8   -1.118341325323      -13.53       -7.58   0.80    1.0   19.5ms
  9   -1.118341325323      -14.95       -8.12   0.80    1.0   19.5ms
 10   -1.118341325323   +  -14.65       -8.55   0.80    1.0   19.4ms
n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118341325323                   -9.05   0.80    1.0   10.9ms
  2   -1.118341325323   +  -15.35      -10.36   0.80    1.0   17.2ms

┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n           Energy  log10(ΔE)  max(Force)   Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   1 │ -1.118341325323 │     -2.96 │ 0.00297052 │  1.86s │
└─────┴─────────────────┴───────────┴────────────┴────────┘

n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118348017537                   -3.09   0.80    1.0   15.0ms
  2   -1.118348160786       -6.84       -4.11   0.80    1.0   17.0ms
  3   -1.118348162263       -8.83       -4.86   0.80    1.0   17.1ms
  4   -1.118348162289      -10.59       -5.78   0.80    1.0   21.6ms
  5   -1.118348162289      -12.08       -6.55   0.80    1.0   25.5ms
  6   -1.118348162290      -13.41       -7.04   0.80    1.0   25.6ms
  7   -1.118348162290      -14.61       -8.13   0.80    1.0   22.2ms
  8   -1.118348162290      -14.65       -8.95   0.80    1.0   19.0ms
n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118356019869                   -2.61   0.80    1.0   23.2ms
  2   -1.118357382893       -5.87       -3.62   0.80    1.0   17.5ms
  3   -1.118357396968       -7.85       -4.37   0.80    1.0   17.4ms
  4   -1.118357397214       -9.61       -5.29   0.80    1.0   17.6ms
  5   -1.118357397222      -11.10       -6.06   0.80    1.0   17.9ms
  6   -1.118357397222      -12.45       -6.55   0.80    1.0   18.2ms
  7   -1.118357397222      -13.62       -7.65   0.80    1.0   18.8ms
  8   -1.118357397222      -14.40       -8.45   0.80    1.0   18.9ms
n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118357397222                   -9.24   0.80    1.0   10.9ms
  2   -1.118357397222   +  -14.57       -9.82   0.80    1.0   17.1ms

┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n           Energy  log10(ΔE)  max(Force)   Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   2 │ -1.118357397222 │     -4.79 │ 4.47258e-5 │  1.80s │
└─────┴─────────────────┴───────────┴────────────┴────────┘

n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118357398684                   -4.94   0.80    1.0   11.0ms
  2   -1.118357398714      -10.54       -5.96   0.80    1.0   17.0ms
  3   -1.118357398714      -12.52       -6.70   0.80    1.0   17.2ms
  4   -1.118357398714      -14.45       -7.62   0.80    1.0   17.6ms
  5   -1.118357398714      -14.51       -8.39   0.80    1.0   18.1ms
  6   -1.118357398714   +    -Inf       -8.89   0.80    1.0   18.3ms
n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118357400607                   -4.40   0.80    1.0   10.9ms
  2   -1.118357400953       -9.46       -5.42   0.80    1.0   17.2ms
  3   -1.118357400956      -11.45       -6.17   0.80    1.0   17.1ms
  4   -1.118357400956      -13.22       -7.09   0.80    1.0   17.5ms
  5   -1.118357400956      -14.65       -7.86   0.80    1.0   17.9ms
  6   -1.118357400956      -15.35       -8.35   0.80    1.0   32.6ms
  7   -1.118357400956   +    -Inf       -9.45   0.80    1.0   18.9ms
n     Energy            log10(ΔE)   log10(Δρ)   α      Diag   Δtime
---   ---------------   ---------   ---------   ----   ----   ------
  1   -1.118357400956                  -10.08   0.80    1.0   11.1ms
  2   -1.118357400956   +  -15.18      -11.02   0.80    1.0   17.1ms

┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n           Energy  log10(ΔE)  max(Force)   Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   3 │ -1.118357400956 │     -8.43 │ 1.15243e-8 │  1.73s │
└─────┴─────────────────┴───────────┴────────────┴────────┘

Structure after optimisation (note that the atom has wrapped around)

results.system
FlexibleSystem(H₂, periodic = TTT):
    bounding_box      : [      10        0        0;
                                0       10        0;
                                0        0       10]u"a₀"

    Atom(H,  [-0.0431826, 1.28704e-10, -2.04063e-10]u"a₀")
    Atom(H,  [ 1.44318, 1.27281e-10, -2.02326e-10]u"a₀")

       .------------.  
      /|            |  
     / |            |  
    * H|         H  |  
    |  |            |  
    |  |            |  
    |  .------------.  
    | /            /   
    |/            /    
    *------------*     

Compute final bond length:

rmin = norm(position(results.system[1]) - position(results.system[2]))
println("Optimal bond length: ", rmin)
Optimal bond length: 1.4863652613886027 a₀

Our results (1.486 Bohr) agrees with the equivalent tutorial from ABINIT.