NEB calculations parallelized over images

The Surface diffusion energy barriers using the Nudged Elastic Band (NEB) method example can be used with GPAW like this:

from ase.io import read
from ase.constraints import FixAtoms
from ase.mep import NEB
from ase.optimize import BFGS
from gpaw.mpi import rank, size

from gpaw import GPAW

initial = read('initial.traj')
final = read('final.traj')

constraint = FixAtoms(mask=[atom.tag > 1 for atom in initial])

n = size // 3      # number of cpu's per image
j = 1 + rank // n  # my image number
assert 3 * n == size

images = [initial]

for i in range(3):
    ranks = range(i * n, (i + 1) * n)
    image = initial.copy()

    if rank in ranks:

        calc = GPAW(mode='fd',
                    h=0.3,
                    kpts=(2, 2, 1),
                    txt=f'neb{j}.txt',
                    communicator=ranks)

        image.calc = calc

    image.set_constraint(constraint)
    images.append(image)

images.append(final)

neb = NEB(images, parallel=True, climb=True)
neb.interpolate()

qn = BFGS(neb, logfile='qn.log', trajectory='neb.traj')
qn.run(fmax=0.05)

If we run the job on 12 cpu’s:

$ mpiexec -np 12 gpaw python neb.py

then each of the three internal images will be parallelized over 4 cpu’s.

The results are read with:

$ ase gui neb.traj@-5:

The energy barrier is found to be ~0.3 eV for LDA.