Using a different Optimization.jl compatible solver

In this example we perform the simplistic optimisation the bond length of a Hydrogen molecule using a trust region quasi-Newton method from NLopt.

We create a calculator employing the density-functional toolkit to compute energies and forces at using the LDA density functional.

using DFTK
using PseudoPotentialData

pseudopotentials = PseudoFamily("dojo.nc.sr.lda.v0_4_1.oncvpsp3.standard.upf")
model_kwargs = (; functionals=LDA(), pseudopotentials)
basis_kwargs = (; kgrid=(1, 1, 1), Ecut=20.0)
calc = DFTKCalculator(; model_kwargs, basis_kwargs)
DFTKCalculator(functionals=Xc(lda_x, lda_c_pw), pseudopotentials=PseudoFamily("dojo.nc.sr.lda.v0_4_1.oncvpsp3.standard.upf"), Ecut=20.0, kgrid=(1, 1, 1))

and we build the hydrogen molecular system, where we attach pseudopotential information for DFTK:

using AtomsBuilder
using Unitful
using UnitfulAtomic

bounding_box = [[10.0, 0.0, 0.0], [0.0, 10.0, 0.0], [0.0, 0.0, 10.0]]u"Å"
system = periodic_system([:H => [0, 0, 1.]u"bohr",
                          :H => [0, 0, 3.]u"bohr"],
                         bounding_box)
nothing

We now run GeometryOptimization.minimize_energy!, but notably pass the NLopt.LD_TNEWTON solver from NLopt as the third argument to employ this solver. Extra keyword argument to NLopt can be added, e.g. here the maxevel=100, which limits the solver to 100 function evaluations:

using GeometryOptimization
using OptimizationNLopt
solver = NLopt.LD_TNEWTON()

results = minimize_energy!(system, calc, solver;
                           tol_forces=1e-4u"eV/Å", verbosity=1,
                           maxeval=100)
nothing
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.105713124409 │           │  0.0820222 │  23.7s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.105713128077 │           │  0.0820221 │  29.7s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.105713125945 │           │  0.0820222 │  35.6s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.105713123945 │           │  0.0820222 │  41.5s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.102361689098 │           │   0.224993 │  52.8s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.131248584688 │           │  0.0436285 │  64.8s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.131248586639 │           │  0.0436284 │  70.8s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.134305430805 │           │  0.0275855 │  81.4s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.134305432039 │           │  0.0275855 │  87.7s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135280022438 │           │ 0.00235839 │  96.8s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135280022544 │           │ 0.00235837 │   103s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135288284108 │           │  5.4289e-5 │   111s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135288284110 │           │ 5.42729e-5 │   117s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135288288523 │           │  4.2555e-6 │   125s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135288288523 │           │ 4.24147e-6 │   130s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135288288523 │           │ 4.25592e-6 │   136s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135288288523 │           │ 4.25016e-6 │   142s │
Geometry optimisation convergence (in atomic units)
┌─────┬─────────────────┬───────────┬────────────┬────────┐
│   n │          Energy │ log10(ΔE) │ max(Force) │  Δtime │
├─────┼─────────────────┼───────────┼────────────┼────────┤
│   0 │ -1.135288288551 │           │ 3.51837e-7 │   149s │
┌ Warning: NLopt failed to converge: FORCED_STOP
@ OptimizationNLopt ~/.julia/packages/OptimizationNLopt/LsFvb/src/OptimizationNLopt.jl:209

The final hydrogen bond length is:

using LinearAlgebra
norm(position(results.system[1]) - position(results.system[2]))
1.4524747684369537 a₀