Skip to article frontmatterSkip to article content
Site not loading correctly?

This may be due to an incorrect BASE_URL configuration. See the MyST Documentation for reference.

Intro to Adsorption Energies

Tutorial Overview
PropertyValue
DifficultyBeginner
Time15-30 minutes
PrerequisitesBasic Python, familiarity with ASE
GoalCalculate adsorption energies using UMA models

To introduce OCP we start with using it to calculate adsorption energies for a simple, atomic adsorbate where we specify the site we want to the adsorption energy for. Conceptually, you do this like you would do it with density functional theory. You create a slab model for the surface, place an adsorbate on it as an initial guess, run a relaxation to get the lowest energy geometry, and then compute the adsorption energy using reference states for the adsorbate.

Intro to Adsorption energies

Adsorption energies are always a reaction energy (an adsorbed species relative to some implied combination of reactants). There are many common schemes in the catalysis literature.

For example, you may want the adsorption energy of oxygen, and you might compute that from this reaction:

1/2 O2 + slab -> slab-O

DFT has known errors with the energy of a gas-phase O2 molecule, so it’s more common to compute this energy relative to a linear combination of H2O and H2. The suggested reference scheme for consistency with OC20 is a reaction

x CO + (x + y/2 - z) H2 + (z-x) H2O + w/2 N2 + * -> CxHyOzNw*

Here, x=y=w=0, z=1, so the reaction ends up as

-H2 + H2O + * -> O*

or alternatively,

H2O + * -> O* + H2

It is possible through thermodynamic cycles to compute other reactions. If we can look up rH1 below and compute rH2

H2 + 1/2 O2 -> H2O  re1 = -3.03 eV, from exp
H2O + * -> O* + H2  re2  # Get from UMA

Then, the adsorption energy for

1/2O2 + * -> O*

is just re1 + re2.

Based on https://atct.anl.gov/Thermochemical Data/version 1.118/species/?species_number=986, the formation energy of water is about -3.03 eV at standard state experimentally. You could also compute this using DFT, but you would probably get the wrong answer for this.

The first step is getting a checkpoint for the model we want to use. UMA is currently the state-of-the-art model and will provide total energy estimates at the RPBE level of theory if you use the “OC20” task.

If you find your kernel is crashing, it probably means you have exceeded the allowed amount of memory. This checkpoint works fine in this example, but it may crash your kernel if you use it in the NRR example.

This next cell will automatically download the checkpoint from huggingface and load it.

from __future__ import annotations

from fairchem.core import FAIRChemCalculator, pretrained_mlip

predictor = pretrained_mlip.get_predict_unit("uma-s-1p1")
calc = FAIRChemCalculator(predictor, task_name="oc20")
WARNING:root:device was not explicitly set, using device='cuda'.

Next we can build a slab with an adsorbate on it. Here we use the ASE module to build a Pt slab. We use the experimental lattice constant that is the default. This can introduce some small errors with DFT since the lattice constant can differ by a few percent, and it is common to use DFT lattice constants. In this example, we do not constrain any layers.

from ase.build import add_adsorbate, fcc111
from ase.optimize import BFGS
# reference energies from a linear combination of H2O/N2/CO/H2!
atomic_reference_energies = {
    "H": -3.477,
    "N": -8.083,
    "O": -7.204,
    "C": -7.282,
}

re1 = -3.03

slab = fcc111("Pt", size=(2, 2, 5), vacuum=20.0)
slab.pbc = True

adslab = slab.copy()
add_adsorbate(adslab, "O", height=1.2, position="fcc")

slab.set_calculator(calc)
opt = BFGS(slab)
opt.run(fmax=0.05, steps=100)
slab_e = slab.get_potential_energy()

adslab.set_calculator(calc)
opt = BFGS(adslab)
opt.run(fmax=0.05, steps=100)
adslab_e = adslab.get_potential_energy()

# Energy for ((H2O-H2) + * -> *O) + (H2 + 1/2O2 -> H2) leads to 1/2O2 + * -> *O!
adslab_e - slab_e - atomic_reference_energies["O"] + re1
/tmp/ipykernel_9806/3752951811.py:17: FutureWarning: Please use atoms.calc = calc
  slab.set_calculator(calc)
      Step     Time          Energy          fmax
BFGS:    0 18:57:31     -104.695193        0.709591
BFGS:    1 18:57:31     -104.753187        0.607748
BFGS:    2 18:57:31     -104.906057        0.369388
BFGS:    3 18:57:31     -104.938753        0.439757
BFGS:    4 18:57:32     -105.016144        0.464013
BFGS:    5 18:57:32     -105.076561        0.356061
BFGS:    6 18:57:32     -105.112621        0.189415
BFGS:    7 18:57:32     -105.126758        0.045372
/tmp/ipykernel_9806/3752951811.py:22: FutureWarning: Please use atoms.calc = calc
  adslab.set_calculator(calc)
      Step     Time          Energy          fmax
BFGS:    0 18:57:32     -110.048788        1.757515
BFGS:    1 18:57:32     -110.230720        0.986276
BFGS:    2 18:57:32     -110.379353        0.731504
BFGS:    3 18:57:32     -110.430289        0.807338
BFGS:    4 18:57:33     -110.543670        0.692009
BFGS:    5 18:57:33     -110.617475        0.492773
BFGS:    6 18:57:33     -110.673457        0.666983
BFGS:    7 18:57:33     -110.721678        0.710522
BFGS:    8 18:57:33     -110.756902        0.443872
BFGS:    9 18:57:33     -110.769294        0.214552
BFGS:   10 18:57:34     -110.772452        0.091728
BFGS:   11 18:57:34     -110.772964        0.062466
BFGS:   12 18:57:34     -110.773259        0.046193
-1.4725015361139024

It is good practice to look at your geometries to make sure they are what you expect.

import matplotlib.pyplot as plt
from ase.visualize.plot import plot_atoms

fig, axs = plt.subplots(1, 2)
plot_atoms(slab, axs[0])
plot_atoms(slab, axs[1], rotation=("-90x"))
axs[0].set_axis_off()
axs[1].set_axis_off()
<Figure size 640x480 with 2 Axes>
import matplotlib.pyplot as plt
from ase.visualize.plot import plot_atoms

fig, axs = plt.subplots(1, 2)
plot_atoms(adslab, axs[0])
plot_atoms(adslab, axs[1], rotation=("-90x"))
axs[0].set_axis_off()
axs[1].set_axis_off()
<Figure size 640x480 with 2 Axes>

How did we do? We need a reference point. In the paper below, there is an atomic adsorption energy for O on Pt(111) of about -4.264 eV. This is for the reaction O + * -> O*. To convert this to the dissociative adsorption energy, we have to add the reaction:

1/2 O2 -> O   D = 2.58 eV (expt)

to get a comparable energy of about -1.68 eV. There is about ~0.2 eV difference (we predicted -1.47 eV above, and the reference comparison is -1.68 eV) to account for. The biggest difference is likely due to the differences in exchange-correlation functional. The reference data used the PBE functional, and eSCN was trained on RPBE data. To additional places where there are differences include:

  1. Difference in lattice constant

  2. The reference energy used for the experiment references. These can differ by up to 0.5 eV from comparable DFT calculations.

  3. How many layers are relaxed in the calculation

Some of these differences tend to be systematic, and you can calibrate and correct these, especially if you can augment these with your own DFT calculations.

See convergence study for some additional studies of factors that influence this number.

Exercises

  1. Explore the effect of the lattice constant on the adsorption energy.

  2. Try different sites, including the bridge and top sites. Compare the energies, and inspect the resulting geometries.

Xu, Z., & Kitchin, J. R. (2014). Probing the coverage dependence of site and adsorbate configurational correlations on (111) surfaces of late transition metals. J. Phys. Chem. C, 118(44), 25597–25602. Xu & Kitchin (2014)

Supporting information.

These are atomic adsorption energies:

O + * -> O*

We have to do some work to get comparable numbers from OCP

H2 + 1/2 O2 -> H2O  re1 = -3.03 eV
H2O + * -> O* + H2  re2   # Get from UMA
O -> 1/2 O2         re3 = -2.58 eV

Then, the adsorption energy for

O + * -> O*

is just re1 + re2 + re3.

Here we just look at the fcc site on Pt. First, we get the data stored in the paper.

Next we get the structures and compute their energies. Some subtle points are that we have to account for stoichiometry, and normalize the adsorption energy by the number of oxygens.

First we get a reference energy from the paper (PBE, 0.25 ML O on Pt(111)).

import json

with open("energies.json") as f:
    edata = json.load(f)

with open("structures.json") as f:
    sdata = json.load(f)

edata["Pt"]["O"]["fcc"]["0.25"]
-4.263842000000002

Next, we load data from the SI to get the geometry to start from.

with open("structures.json") as f:
    s = json.load(f)

sfcc = s["Pt"]["O"]["fcc"]["0.25"]

Next, we construct the atomic geometry, run the geometry optimization, and compute the energy.

re3 = -2.58  # O -> 1/2 O2         re3 = -2.58 eV

from ase import Atoms

adslab = Atoms(sfcc["symbols"], positions=sfcc["pos"], cell=sfcc["cell"], pbc=True)

# Grab just the metal surface atoms
slab = adslab[adslab.arrays["numbers"] == adslab.arrays["numbers"][0]]
adsorbates = adslab[~(adslab.arrays["numbers"] == adslab.arrays["numbers"][0])]
slab.set_calculator(calc)
opt = BFGS(slab)
opt.run(fmax=0.05, steps=100)

adslab.set_calculator(calc)
opt = BFGS(adslab)

opt.run(fmax=0.05, steps=100)
re2 = (
    adslab.get_potential_energy()
    - slab.get_potential_energy()
    - sum([atomic_reference_energies[x] for x in adsorbates.get_chemical_symbols()])
)

nO = 0
for atom in adslab:
    if atom.symbol == "O":
        nO += 1
        re2 += re1 + re3

print(re2 / nO)
/tmp/ipykernel_9806/647904475.py:10: FutureWarning: Please use atoms.calc = calc
  slab.set_calculator(calc)
      Step     Time          Energy          fmax
BFGS:    0 18:57:37      -82.862062        1.000262
BFGS:    1 18:57:38      -82.918113        0.752920
BFGS:    2 18:57:38      -83.011694        0.335860
BFGS:    3 18:57:38      -83.016496        0.303120
BFGS:    4 18:57:38      -83.024496        0.224962
BFGS:    5 18:57:38      -83.030036        0.152567
BFGS:    6 18:57:38      -83.033361        0.088625
BFGS:    7 18:57:38      -83.034468        0.065747
BFGS:    8 18:57:39      -83.035209        0.068690
BFGS:    9 18:57:39      -83.035562        0.046754
      Step     Time          Energy          fmax
BFGS:    0 18:57:39      -88.770343        0.312349
/tmp/ipykernel_9806/647904475.py:14: FutureWarning: Please use atoms.calc = calc
  adslab.set_calculator(calc)
BFGS:    1 18:57:39      -88.773630        0.271703
BFGS:    2 18:57:39      -88.784473        0.104050
BFGS:    3 18:57:40      -88.786014        0.101458
BFGS:    4 18:57:40      -88.788553        0.106390
BFGS:    5 18:57:40      -88.790503        0.089089
BFGS:    6 18:57:40      -88.791972        0.096440
BFGS:    7 18:57:41      -88.792659        0.099418
BFGS:    8 18:57:41      -88.793179        0.077363
BFGS:    9 18:57:41      -88.793677        0.036091
-4.164113856060191

Site correlations

This cell reproduces a portion of a figure in the paper. We compare oxygen adsorption energies in the fcc and hcp sites across metals and coverages. These adsorption energies are highly correlated with each other because the adsorption sites are so similar.

At higher coverages, the agreement is not as good. This is likely because the model is extrapolating and needs to be fine-tuned.

import time

from tqdm import tqdm

t0 = time.time()

data = {"fcc": [], "hcp": []}

refdata = {"fcc": [], "hcp": []}


for metal in ["Cu", "Ag", "Pd", "Pt", "Rh", "Ir"]:
    print(metal)
    for site in ["fcc", "hcp"]:
        for adsorbate in ["O"]:
            for coverage in tqdm(["0.25"]):

                entry = s[metal][adsorbate][site][coverage]

                adslab = Atoms(
                    entry["symbols"],
                    positions=entry["pos"],
                    cell=entry["cell"],
                    pbc=True,
                )

                # Grab just the metal surface atoms
                adsorbates = adslab[
                    ~(adslab.arrays["numbers"] == adslab.arrays["numbers"][0])
                ]

                slab = adslab[adslab.arrays["numbers"] == adslab.arrays["numbers"][0]]
                slab.set_calculator(calc)
                opt = BFGS(slab)
                opt.run(fmax=0.05, steps=100)

                adslab.set_calculator(calc)
                opt = BFGS(adslab)
                opt.run(fmax=0.05, steps=100)

                re2 = (
                    adslab.get_potential_energy()
                    - slab.get_potential_energy()
                    - sum(
                        [
                            atomic_reference_energies[x]
                            for x in adsorbates.get_chemical_symbols()
                        ]
                    )
                )

                nO = 0
                for atom in adslab:
                    if atom.symbol == "O":
                        nO += 1
                        re2 += re1 + re3

                data[site] += [re2 / nO]
                refdata[site] += [edata[metal][adsorbate][site][coverage]]

f"Elapsed time = {time.time() - t0} seconds"
Cu
  0%|          | 0/1 [00:00<?, ?it/s]
/tmp/ipykernel_9806/1356342052.py:33: FutureWarning: Please use atoms.calc = calc
  slab.set_calculator(calc)
      Step     Time          Energy          fmax
BFGS:    0 18:57:41      -48.910948        0.646351
BFGS:    1 18:57:42      -48.934031        0.539303
BFGS:    2 18:57:42      -49.000043        0.280103
BFGS:    3 18:57:42      -49.002374        0.251321
BFGS:    4 18:57:42      -49.008810        0.158909
BFGS:    5 18:57:42      -49.013985        0.134401
BFGS:    6 18:57:42      -49.016868        0.071526
BFGS:    7 18:57:42      -49.017633        0.051269
BFGS:    8 18:57:43      -49.018011        0.047564
/tmp/ipykernel_9806/1356342052.py:37: FutureWarning: Please use atoms.calc = calc
  adslab.set_calculator(calc)
      Step     Time          Energy          fmax
BFGS:    0 18:57:43      -55.204221        0.339047
BFGS:    1 18:57:43      -55.206662        0.267938
BFGS:    2 18:57:43      -55.213711        0.160041
BFGS:    3 18:57:43      -55.215874        0.166184
BFGS:    4 18:57:44      -55.219675        0.116680
BFGS:    5 18:57:44      -55.221814        0.079417
BFGS:    6 18:57:44      -55.223360        0.075349
BFGS:    7 18:57:44      -55.224364        0.090287
BFGS:    8 18:57:44      -55.225517        0.098894
BFGS:    9 18:57:45      -55.226530        0.068555
BFGS:   10 18:57:45      -55.226977        0.042996
100%|██████████| 1/1 [00:03<00:00,  3.65s/it]
100%|██████████| 1/1 [00:03<00:00,  3.65s/it]

  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:57:45      -48.936503        0.552695
BFGS:    1 18:57:45      -48.954084        0.470533
BFGS:    2 18:57:45      -49.008577        0.210998
BFGS:    3 18:57:46      -49.009768        0.197685
BFGS:    4 18:57:46      -49.017746        0.041126
      Step     Time          Energy          fmax
BFGS:    0 18:57:46      -55.119070        0.305598
BFGS:    1 18:57:46      -55.120918        0.239019
BFGS:    2 18:57:46      -55.125705        0.129148
BFGS:    3 18:57:47      -55.127340        0.147608
BFGS:    4 18:57:47      -55.130398        0.128614
BFGS:    5 18:57:47      -55.132308        0.103136
BFGS:    6 18:57:47      -55.133382        0.051684
BFGS:    7 18:57:47      -55.133851        0.056058
BFGS:    8 18:57:48      -55.134358        0.071220
BFGS:    9 18:57:48      -55.135043        0.069432
BFGS:   10 18:57:48      -55.135574        0.039921
100%|██████████| 1/1 [00:03<00:00,  3.15s/it]
100%|██████████| 1/1 [00:03<00:00,  3.15s/it]

Ag
  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:57:48      -33.034314        0.634820
BFGS:    1 18:57:48      -33.053470        0.551048
BFGS:    2 18:57:49      -33.120905        0.172021
BFGS:    3 18:57:49      -33.122778        0.157330
BFGS:    4 18:57:49      -33.124181        0.143386
BFGS:    5 18:57:49      -33.126944        0.112458
BFGS:    6 18:57:49      -33.130771        0.110654
BFGS:    7 18:57:50      -33.133992        0.067529
BFGS:    8 18:57:50      -33.134835        0.034887
      Step     Time          Energy          fmax
BFGS:    0 18:57:50      -38.173509        0.098015
BFGS:    1 18:57:50      -38.174435        0.094365
BFGS:    2 18:57:51      -38.184909        0.084033
BFGS:    3 18:57:51      -38.185732        0.085434
BFGS:    4 18:57:51      -38.189035        0.086084
BFGS:    5 18:57:51      -38.190952        0.068838
BFGS:    6 18:57:52      -38.192256        0.087453
BFGS:    7 18:57:52      -38.193188        0.076912
BFGS:    8 18:57:52      -38.193675        0.042571
100%|██████████| 1/1 [00:04<00:00,  4.32s/it]
100%|██████████| 1/1 [00:04<00:00,  4.32s/it]

  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:57:53      -33.056750        0.557366
BFGS:    1 18:57:53      -33.071228        0.489412
BFGS:    2 18:57:53      -33.126396        0.123825
BFGS:    3 18:57:53      -33.127308        0.112960
BFGS:    4 18:57:53      -33.128516        0.097373
BFGS:    5 18:57:54      -33.130261        0.079537
BFGS:    6 18:57:54      -33.133158        0.076350
BFGS:    7 18:57:54      -33.134741        0.040732
      Step     Time          Energy          fmax
BFGS:    0 18:57:55      -38.102025        0.092717
BFGS:    1 18:57:55      -38.102856        0.088836
BFGS:    2 18:57:55      -38.111088        0.134656
BFGS:    3 18:57:55      -38.112012        0.120324
BFGS:    4 18:57:56      -38.115289        0.063997
BFGS:    5 18:57:56      -38.116155        0.056402
BFGS:    6 18:57:56      -38.116710        0.058194
BFGS:    7 18:57:56      -38.117581        0.073183
BFGS:    8 18:57:56      -38.118451        0.057896
100%|██████████| 1/1 [00:04<00:00,  4.32s/it]
100%|██████████| 1/1 [00:04<00:00,  4.33s/it]

BFGS:    9 18:57:57      -38.118834        0.026042
Pd
  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:57:57      -70.165062        0.631446
BFGS:    1 18:57:57      -70.189910        0.508128
BFGS:    2 18:57:57      -70.241649        0.191934
BFGS:    3 18:57:58      -70.243346        0.182565
BFGS:    4 18:57:58      -70.248576        0.144558
BFGS:    5 18:57:58      -70.251930        0.115933
BFGS:    6 18:57:58      -70.254308        0.074837
BFGS:    7 18:57:58      -70.255204        0.059446
BFGS:    8 18:57:59      -70.255944        0.041108
      Step     Time          Energy          fmax
BFGS:    0 18:57:59      -76.130241        0.213656
BFGS:    1 18:57:59      -76.133242        0.190305
BFGS:    2 18:58:00      -76.148717        0.169645
BFGS:    3 18:58:00      -76.150812        0.152254
BFGS:    4 18:58:00      -76.154312        0.130504
BFGS:    5 18:58:00      -76.157058        0.106138
BFGS:    6 18:58:01      -76.159302        0.088115
BFGS:    7 18:58:01      -76.160431        0.080499
BFGS:    8 18:58:01      -76.161130        0.064502
BFGS:    9 18:58:02      -76.161614        0.042989
100%|██████████| 1/1 [00:05<00:00,  5.11s/it]
100%|██████████| 1/1 [00:05<00:00,  5.12s/it]

  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:58:02      -70.197711        0.454636
BFGS:    1 18:58:02      -70.211774        0.373893
BFGS:    2 18:58:03      -70.245977        0.178178
BFGS:    3 18:58:03      -70.247162        0.164572
BFGS:    4 18:58:03      -70.253348        0.078706
BFGS:    5 18:58:03      -70.254511        0.075507
BFGS:    6 18:58:04      -70.255271        0.058129
BFGS:    7 18:58:04      -70.255844        0.036952
      Step     Time          Energy          fmax
BFGS:    0 18:58:04      -75.892387        0.169440
BFGS:    1 18:58:04      -75.895234        0.150799
BFGS:    2 18:58:05      -75.906632        0.144106
BFGS:    3 18:58:05      -75.908424        0.135675
BFGS:    4 18:58:05      -75.912427        0.129022
BFGS:    5 18:58:06      -75.915495        0.108690
BFGS:    6 18:58:06      -75.917927        0.075286
BFGS:    7 18:58:06      -75.918955        0.067276
BFGS:    8 18:58:06      -75.919479        0.044570
100%|██████████| 1/1 [00:04<00:00,  4.79s/it]
100%|██████████| 1/1 [00:04<00:00,  4.79s/it]

Pt
  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:58:07      -82.862062        1.000262
BFGS:    1 18:58:07      -82.918113        0.752919
BFGS:    2 18:58:07      -83.011694        0.335860
BFGS:    3 18:58:08      -83.016496        0.303120
BFGS:    4 18:58:08      -83.024497        0.224963
BFGS:    5 18:58:08      -83.030035        0.152566
BFGS:    6 18:58:08      -83.033360        0.088625
BFGS:    7 18:58:09      -83.034469        0.065747
BFGS:    8 18:58:09      -83.035209        0.068689
BFGS:    9 18:58:09      -83.035564        0.046754
      Step     Time          Energy          fmax
BFGS:    0 18:58:09      -88.770343        0.312348
BFGS:    1 18:58:10      -88.773630        0.271702
BFGS:    2 18:58:10      -88.784473        0.104050
BFGS:    3 18:58:10      -88.786014        0.101457
BFGS:    4 18:58:10      -88.788552        0.106391
BFGS:    5 18:58:11      -88.790502        0.089089
BFGS:    6 18:58:11      -88.791972        0.096441
BFGS:    7 18:58:11      -88.792660        0.099423
BFGS:    8 18:58:11      -88.793179        0.077357
BFGS:    9 18:58:12      -88.793678        0.036095
100%|██████████| 1/1 [00:05<00:00,  5.05s/it]
100%|██████████| 1/1 [00:05<00:00,  5.05s/it]

  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:58:12      -82.946995        0.675644
BFGS:    1 18:58:12      -82.973149        0.546704
BFGS:    2 18:58:12      -83.026506        0.203005
BFGS:    3 18:58:12      -83.028070        0.183991
BFGS:    4 18:58:13      -83.033697        0.083372
BFGS:    5 18:58:13      -83.034396        0.061085
BFGS:    6 18:58:13      -83.035121        0.048391
      Step     Time          Energy          fmax
BFGS:    0 18:58:14      -88.350562        0.234982
BFGS:    1 18:58:14      -88.353605        0.203898
BFGS:    2 18:58:14      -88.362976        0.117094
BFGS:    3 18:58:14      -88.364342        0.121252
BFGS:    4 18:58:15      -88.368499        0.090134
BFGS:    5 18:58:15      -88.370482        0.073282
BFGS:    6 18:58:15      -88.371785        0.075223
BFGS:    7 18:58:15      -88.372318        0.079578
BFGS:    8 18:58:16      -88.372731        0.061423
BFGS:    9 18:58:16      -88.373020        0.027205
100%|██████████| 1/1 [00:04<00:00,  4.30s/it]
100%|██████████| 1/1 [00:04<00:00,  4.30s/it]

Rh
  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:58:16     -100.192566        0.688933
BFGS:    1 18:58:17     -100.218240        0.600213
BFGS:    2 18:58:17     -100.286834        0.174739
BFGS:    3 18:58:17     -100.289767        0.129382
BFGS:    4 18:58:17     -100.298129        0.066639
BFGS:    5 18:58:18     -100.299264        0.053478
BFGS:    6 18:58:18     -100.299777        0.045082
      Step     Time          Energy          fmax
BFGS:    0 18:58:18     -106.921794        0.203149
BFGS:    1 18:58:18     -106.924540        0.168907
BFGS:    2 18:58:19     -106.930976        0.052441
BFGS:    3 18:58:19     -106.931584        0.045559
100%|██████████| 1/1 [00:03<00:00,  3.20s/it]
100%|██████████| 1/1 [00:03<00:00,  3.20s/it]

  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:58:19     -100.171741        0.714774
BFGS:    1 18:58:20     -100.203866        0.626030
BFGS:    2 18:58:20     -100.287845        0.263192
BFGS:    3 18:58:20     -100.290543        0.195991
BFGS:    4 18:58:21     -100.297165        0.080377
BFGS:    5 18:58:21     -100.299032        0.056460
BFGS:    6 18:58:21     -100.299772        0.043871
      Step     Time          Energy          fmax
BFGS:    0 18:58:21     -106.876839        0.168839
BFGS:    1 18:58:22     -106.879420        0.139489
BFGS:    2 18:58:22     -106.885549        0.071255
BFGS:    3 18:58:22     -106.886238        0.064285
BFGS:    4 18:58:22     -106.886580        0.054980
100%|██████████| 1/1 [00:03<00:00,  3.48s/it]
100%|██████████| 1/1 [00:03<00:00,  3.48s/it]

BFGS:    5 18:58:23     -106.886842        0.043052
Ir
  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:58:23     -124.235610        1.189926
BFGS:    1 18:58:23     -124.310841        0.937183
BFGS:    2 18:58:23     -124.422189        0.140826
BFGS:    3 18:58:24     -124.425270        0.133780
BFGS:    4 18:58:24     -124.429598        0.060643
BFGS:    5 18:58:24     -124.430258        0.038672
      Step     Time          Energy          fmax
BFGS:    0 18:58:25     -130.605623        0.371646
BFGS:    1 18:58:25     -130.616146        0.265890
BFGS:    2 18:58:25     -130.630316        0.110763
BFGS:    3 18:58:25     -130.633615        0.094356
BFGS:    4 18:58:26     -130.634276        0.075693
BFGS:    5 18:58:26     -130.634990        0.065689
BFGS:    6 18:58:26     -130.636200        0.083046
BFGS:    7 18:58:26     -130.636656        0.075129
BFGS:    8 18:58:27     -130.637048        0.047434
100%|██████████| 1/1 [00:04<00:00,  4.11s/it]
100%|██████████| 1/1 [00:04<00:00,  4.11s/it]

  0%|          | 0/1 [00:00<?, ?it/s]
      Step     Time          Energy          fmax
BFGS:    0 18:58:27     -124.227593        1.165685
BFGS:    1 18:58:27     -124.311476        0.915619
BFGS:    2 18:58:27     -124.422556        0.165138
BFGS:    3 18:58:27     -124.425477        0.200268
BFGS:    4 18:58:28     -124.428763        0.168862
BFGS:    5 18:58:28     -124.430253        0.088404
BFGS:    6 18:58:28     -124.430560        0.041796
      Step     Time          Energy          fmax
BFGS:    0 18:58:28     -130.496158        0.405482
BFGS:    1 18:58:28     -130.507202        0.266991
BFGS:    2 18:58:28     -130.518615        0.126414
BFGS:    3 18:58:29     -130.522313        0.080268
BFGS:    4 18:58:29     -130.523179        0.066553
BFGS:    5 18:58:29     -130.523452        0.056611
BFGS:    6 18:58:29     -130.524511        0.055659
BFGS:    7 18:58:30     -130.524677        0.042692
100%|██████████| 1/1 [00:03<00:00,  3.05s/it]
100%|██████████| 1/1 [00:03<00:00,  3.05s/it]

'Elapsed time = 48.602822065353394 seconds'

First, we compare the computed data and reference data. There is a systematic difference of about 0.5 eV due to the difference between RPBE and PBE functionals, and other subtle differences like lattice constant differences and reference energy differences. This is pretty typical, and an expected deviation.

plt.plot(refdata["fcc"], data["fcc"], "r.", label="fcc")
plt.plot(refdata["hcp"], data["hcp"], "b.", label="hcp")
plt.plot([-5.5, -3.5], [-5.5, -3.5], "k-")
plt.xlabel("Ref. data (DFT)")
plt.ylabel("UMA-OC20 prediction");
<Figure size 640x480 with 1 Axes>

Next we compare the correlation between the hcp and fcc sites. Here we see the same trends. The data falls below the parity line because the hcp sites tend to be a little weaker binding than the fcc sites.

plt.plot(refdata["hcp"], refdata["fcc"], "r.")
plt.plot(data["hcp"], data["fcc"], ".")
plt.plot([-6, -1], [-6, -1], "k-")
plt.xlabel("$H_{ads, hcp}$")
plt.ylabel("$H_{ads, fcc}$")
plt.legend(["DFT (PBE)", "UMA-OC20"]);
<Figure size 640x480 with 1 Axes>

Exercises

  1. You can also explore a few other adsorbates: C, H, N.

  2. Explore the higher coverages. The deviations from the reference data are expected to be higher, but relative differences tend to be better. You probably need fine tuning to improve this performance. This data set doesn’t have forces though, so it isn’t practical to do it here.

Next steps

In the next step, we consider some more complex adsorbates in nitrogen reduction, and how we can leverage OCP to automate the search for the most stable adsorbate geometry. See the next step.

Convergence study

In the adsorption energies section we discussed some possible reasons we might see a discrepancy. Here we investigate some factors that impact the computed energies.

In this section, the energies refer to the reaction 1/2 O2 -> O*.

Effects of number of layers

Slab thickness could be a factor. Here we relax the whole slab, and see by about 4 layers the energy is converged to ~0.02 eV.

for nlayers in [3, 4, 5, 6, 7, 8]:
    slab = fcc111("Pt", size=(2, 2, nlayers), vacuum=10.0)

    slab.pbc = True
    slab.set_calculator(calc)
    opt_slab = BFGS(slab, logfile=None)
    opt_slab.run(fmax=0.05, steps=100)
    slab_e = slab.get_potential_energy()

    adslab = slab.copy()
    add_adsorbate(adslab, "O", height=1.2, position="fcc")

    adslab.pbc = True
    adslab.set_calculator(calc)
    opt_adslab = BFGS(adslab, logfile=None)
    opt_adslab.run(fmax=0.05, steps=100)
    adslab_e = adslab.get_potential_energy()

    print(
        f"nlayers = {nlayers}: {adslab_e - slab_e - atomic_reference_energies['O'] + re1:1.2f} eV"
    )
/tmp/ipykernel_9806/338101817.py:5: FutureWarning: Please use atoms.calc = calc
  slab.set_calculator(calc)
/tmp/ipykernel_9806/338101817.py:14: FutureWarning: Please use atoms.calc = calc
  adslab.set_calculator(calc)
nlayers = 3: -1.58 eV
nlayers = 4: -1.47 eV
nlayers = 5: -1.47 eV
nlayers = 6: -1.46 eV
nlayers = 7: -1.47 eV
nlayers = 8: -1.47 eV

Effects of relaxation

It is common to only relax a few layers, and constrain lower layers to bulk coordinates. We do that here. We only relax the adsorbate and the top layer.

This has a small effect (0.1 eV).

from ase.constraints import FixAtoms

for nlayers in [3, 4, 5, 6, 7, 8]:
    slab = fcc111("Pt", size=(2, 2, nlayers), vacuum=10.0)

    slab.set_constraint(FixAtoms(mask=[atom.tag > 1 for atom in slab]))
    slab.pbc = True
    slab.set_calculator(calc)
    opt_slab = BFGS(slab, logfile=None)
    opt_slab.run(fmax=0.05, steps=100)
    slab_e = slab.get_potential_energy()

    adslab = slab.copy()
    add_adsorbate(adslab, "O", height=1.2, position="fcc")

    adslab.set_constraint(FixAtoms(mask=[atom.tag > 1 for atom in adslab]))
    adslab.pbc = True
    adslab.set_calculator(calc)
    opt_adslab = BFGS(adslab, logfile=None)
    opt_adslab.run(fmax=0.05, steps=100)
    adslab_e = adslab.get_potential_energy()

    print(
        f"nlayers = {nlayers}: {adslab_e - slab_e - atomic_reference_energies['O'] + re1:1.2f} eV"
    )
/tmp/ipykernel_9806/1426773950.py:8: FutureWarning: Please use atoms.calc = calc
  slab.set_calculator(calc)
/tmp/ipykernel_9806/1426773950.py:18: FutureWarning: Please use atoms.calc = calc
  adslab.set_calculator(calc)
nlayers = 3: -1.47 eV
nlayers = 4: -1.35 eV
nlayers = 5: -1.35 eV
nlayers = 6: -1.34 eV
nlayers = 7: -1.35 eV
nlayers = 8: -1.35 eV

Unit cell size

Coverage effects are quite noticeable with oxygen. Here we consider larger unit cells. This effect is large, and the results don’t look right, usually adsorption energies get more favorable at lower coverage, not less. This suggests fine-tuning could be important even at low coverages.

for size in [1, 2, 3, 4, 5]:

    slab = fcc111("Pt", size=(size, size, 5), vacuum=10.0)

    slab.set_constraint(FixAtoms(mask=[atom.tag > 1 for atom in slab]))
    slab.pbc = True
    slab.set_calculator(calc)
    opt_slab = BFGS(slab, logfile=None)
    opt_slab.run(fmax=0.05, steps=100)
    slab_e = slab.get_potential_energy()

    adslab = slab.copy()
    add_adsorbate(adslab, "O", height=1.2, position="fcc")

    adslab.set_constraint(FixAtoms(mask=[atom.tag > 1 for atom in adslab]))
    adslab.pbc = True
    adslab.set_calculator(calc)
    opt_adslab = BFGS(adslab, logfile=None)
    opt_adslab.run(fmax=0.05, steps=100)
    adslab_e = adslab.get_potential_energy()

    print(
        f"({size}x{size}): {adslab_e - slab_e - atomic_reference_energies['O'] + re1:1.2f} eV"
    )
/tmp/ipykernel_9806/3371624330.py:7: FutureWarning: Please use atoms.calc = calc
  slab.set_calculator(calc)
/tmp/ipykernel_9806/3371624330.py:17: FutureWarning: Please use atoms.calc = calc
  adslab.set_calculator(calc)
(1x1): -0.12 eV
(2x2): -1.35 eV
(3x3): -1.38 eV
(4x4): -1.41 eV
(5x5): -1.42 eV

Summary

As with DFT, you should take care to see how these kinds of decisions affect your results, and determine if they would change any interpretations or not.

References
  1. Xu, Z., & Kitchin, J. R. (2014). Probing the Coverage Dependence of Site and Adsorbate Configurational Correlations on (111) Surfaces of Late Transition Metals. The Journal of Physical Chemistry C, 118(44), 25597–25602. 10.1021/jp508805h