Image Dependent Pair Potential for improved interpolation of NEB initial guess

Reference: S. Smidstrup, A. Pedersen, K. Stokbro and H. Jonsson, Improved initial guess for minimum energy path calculations, J. Chem. Phys. 140, 214106 (2014).

Use of the NEB method is dependent upon generating an initial guess for the images lying between the initial and final states. The most simple approach is to use linear interpolation of the atomic coordinates. However, this can be problematic as the quality of the interpolated path can ofter be far from the real one. The implication being that a lot of time is spent in the NEB routine optimising the shape of the path, before the transition state is homed-in upon.

The image dependent pair potential is a method that has been developed to provide an improvement to the initial guess for the NEB path. The IDPP method uses the bond distance between the atoms involved in the transition state to create target structures for the images, rather than interpolating the atomic positions. By defining an objective function in terms of the distances between atoms, the NEB algorithm is used with this image dependent pair potential (IDPP) to create the initial guess for the full NEB calculation.

Note: The examples below utilise the EMT calculator for illustrative purposes, the results should not be over interpreted.

Example 1: Ethane

This example illustrates the use of the IDPP interpolation scheme to generate an initial guess for rotation of a methyl group around the CC bond.

Using the standard linear interpolation approach, as in the following example, we can see that 46 iterations are required to find the transition state.

from ase.build import molecule
from ase.neb import NEB
from ase.calculators.emt import EMT
from ase.optimize.fire import FIRE as QuasiNewton

# Optimise molecule
initial = molecule('C2H6')
initial.set_calculator(EMT())
relax = QuasiNewton(initial)
relax.run(fmax=0.05)

# Create final state
final = initial.copy()
final.positions[2:5] = initial.positions[[3, 4, 2]]

# Generate blank images
images = [initial]

for i in range(9):
    images.append(initial.copy())

for image in images:
    image.set_calculator(EMT())
   
images.append(final)

# Run IDPP interpolation
neb = NEB(images)
neb.interpolate()

# Run NEB calculation
qn = QuasiNewton(neb, trajectory='ethane_linear.traj',
                 logfile='ethane_linear.log')
qn.run(fmax=0.05)

However if we modify our script slightly and use the IDPP method to find the initial guess, we can see that the number of iterations required to find the transition state is reduced to 14.

from ase.build import molecule
from ase.neb import NEB
from ase.calculators.emt import EMT
from ase.optimize.fire import FIRE as QuasiNewton
from ase.visualize import view

#Optimise molecule
initial = molecule('C2H6')
initial.set_calculator(EMT())
relax = QuasiNewton(initial)
relax.run(fmax=0.05)
view(initial)

#Create final state
final = initial.copy()
final.positions[2:5] = initial.positions[[3, 4, 2]]

#Generate blank images
images = [initial]

for i in range(9):
    images.append(initial.copy())

for image in images:
    image.set_calculator(EMT())
   
images.append(final)

#Run IDPP interpolation
neb = NEB(images)
neb.interpolate('idpp')

#Run NEB calculation
qn = QuasiNewton(neb, trajectory='ethane_idpp.traj',  logfile='ethane_idpp.log')
qn.run(fmax=0.05)

Clearly, if one was using a full DFT calculator one can potentially gain a significant time improvement.

Example 2: N diffusion over a step edge

Often we are interested in generating an initial guess for a surface reaction. This example illustrates how we can optimise our initial and final state structures before using the IDPP interpolation to generate our initial guess for the NEB calculation:

import numpy as np
from ase import Atoms
from ase.constraints import FixAtoms
from ase.calculators.emt import EMT
from ase.neb import NEB
from ase.visualize import view
from ase.optimize.fire import FIRE as QuasiNewton
from ase.lattice.cubic import FaceCenteredCubic

#set the number of images you want
nimages = 5

#some algebra to determine surface normal and the plane of the surface
d3=[2,1,1]
a1=np.array([0,1,1])
d1=np.cross(a1,d3)
a2=np.array([0,-1,1])
d2=np.cross(a2,d3)

#create your slab
slab  =FaceCenteredCubic(directions=[d1,d2,d3],
                         size=(2,1,2),
                         symbol=('Pt'),
                         latticeconstant=3.9)

#add some vacuum to your slab
uc = slab.get_cell()
print(uc)
uc[2] += [0,0,10]  #there are ten layers of vacuum
uc = slab.set_cell(uc,scale_atoms=False)
#view the slab to make sure it is how you expect
view(slab)

#some positions needed to place the atom in the correct place
x1 = 1.379
x2 = 4.137
x3 = 2.759
y1 = 0.0
y2 = 2.238
z1 = 7.165
z2 = 6.439


#Add the adatom to the list of atoms and set constraints of surface atoms.
slab += Atoms('N', [ ((x2+x1)/2,y1,z1+1.5)])
mask = [atom.symbol == 'Pt' for atom in slab]
slab.set_constraint(FixAtoms(mask=mask))

#optimise the initial state
# Atom below step
initial = slab.copy()
initial.set_calculator(EMT())
relax = QuasiNewton(initial)
relax.run(fmax=0.05)
view(initial)

#optimise the initial state
# Atom above step
slab[-1].position = (x3,y2+1,z2+3.5)
final = slab.copy()
final.set_calculator(EMT())
relax = QuasiNewton(final)
relax.run(fmax=0.05)
view(final)

#create a list of images for interpolation
images = [initial]
for i in range(nimages):
    images.append(initial.copy())

for image in images:
    image.set_calculator(EMT())

images.append(final)
view(images)

#carry out idpp interpolation
neb = NEB(images)
neb.interpolate('idpp')

#Run NEB calculation
qn = QuasiNewton(neb, trajectory='N_diffusion.traj',  logfile='N_diffusion.log')
qn.run(fmax=0.05)


To again illustrate the potential speedup, the following script which uses the linear interpolation takes 23 iterations to find a transition state, compared to 12 using the IDPP interpolation.

import numpy as np
from ase import Atoms
from ase.constraints import FixAtoms
from ase.calculators.emt import EMT
from ase.neb import NEB
from ase.visualize import view
from ase.optimize.fire import FIRE as QuasiNewton
from ase.lattice.cubic import FaceCenteredCubic

#set the number of images you want
nimages = 5

#some algebra to determine surface normal and the plane of the surface
d3=[2,1,1]
a1=np.array([0,1,1])
d1=np.cross(a1,d3)
a2=np.array([0,-1,1])
d2=np.cross(a2,d3)

#create your slab
slab  =FaceCenteredCubic(directions=[d1,d2,d3],
                         size=(2,1,2),
                         symbol=('Pt'),
                         latticeconstant=3.9)

#add some vacuum to your slab
uc = slab.get_cell()
print(uc)
uc[2] += [0,0,10]  #there are ten layers of vacuum
uc = slab.set_cell(uc,scale_atoms=False)
#view the slab to make sure it is how you expect
view(slab)

#some positions needed to place the atom in the correct place
x1 = 1.379
x2 = 4.137
x3 = 2.759
y1 = 0.0
y2 = 2.238
z1 = 7.165
z2 = 6.439


#Add the adatom to the list of atoms and set constraints of surface atoms.
slab += Atoms('N', [ ((x2+x1)/2,y1,z1+1.5)])
mask = [atom.symbol == 'Pt' for atom in slab]
slab.set_constraint(FixAtoms(mask=mask))

#optimise the initial state
# Atom below step
initial = slab.copy()
initial.set_calculator(EMT())
relax = QuasiNewton(initial)
relax.run(fmax=0.05)
view(initial)

#optimise the initial state
# Atom above step
slab[-1].position = (x3,y2+1,z2+3.5)
final = slab.copy()
final.set_calculator(EMT())
relax = QuasiNewton(final)
relax.run(fmax=0.05)
view(final)

#create a list of images for interpolation
images = [initial]
for i in range(nimages):
    images.append(initial.copy())

for image in images:
    image.set_calculator(EMT())

images.append(final)
view(images)

#carry out idpp interpolation
neb = NEB(images)
#neb.interpolate('idpp')
neb.interpolate()

#Run NEB calculation
qn = QuasiNewton(neb, trajectory='N_diffusion_lin.traj',  logfile='N_diffusion_lin.log')
qn.run(fmax=0.05)