Autoadsorbate
Chemical intuition for surface science in a package.
Install / Use
/learn @basf/AutoadsorbateREADME
Table of Contents
- Installation
- AutoAdsorbate
- Making surrogate SMILES automatically
- Fully automatic - populate Surface with Fragment
Installation
The package is designed to be as lightweight as possible, to implement seamlessly into existing environments with complex dependecies. If you git clone <autoadsorbate> and just sys.path.insert(0, <path/to/autoadsorbate>), most likely it will work.
- Built on only:
aserdkit- Basic Python packages:
pandas,numpy
The package is available on PyPi:
pip install autoadsorbate
Installation from source:
git clone <autoadsorbate>
cd autoadsorbate
pip install .
AutoAdsorbate
AutoAdsorbate is a lightweight and easy-to-use Python package for generating chemically meaningful configurations of molecules and fragments on surfaces. Built with minimal dependencies and a low barrier to entry, it enables rapid setup of surface-adsorbate systems using the Surrogate-SMILES (*SMILES) representation. Ideal for researchers in catalysis, nanotech, and materials science, AutoAdsorbate streamlines dataset generation for simulations and machine learning workflows.
The challenge of generating initial structures for heterogeneous catalysis has traditionally been addressed through manual effort. This package offers an alternative, automated approach.
To effectively simulate reactive behavior at surfaces, it is crucial to establish clear definitions within our framework. The following definitions are essential for accurately characterizing the structures of interest:
-
Fragment:
- <font color='red'>Molecules</font> – species that retain their corresponding geometries even when isolated from the surface.
- <font color='red'>Reactive species</font> – species that adopt their corresponding geometries only when attached to the surface.
-
Surface:
- The surface is defined simply – <font color='red'>every atom of the slab that can be in contact with an intermediate is considered a surface atom</font>. The surface is the collection of such atoms.
- Every atom of the surface is a "top" site.
- When two "top" sites are close (close in its literal meaning), they form a "bridge" site.
- When three "top" sites are close (close in its literal meaning), they form a "3-fold" site.
- etc.
-
Active Site:
- A collection of one or more sites that can facilitate a chemical transformation is called an active site.
- A "top" site can be an active site only for Eley-Rideal transformations.
- All other transformations require that at least one intermediate binds through at least two sites. All involved sites compose an active site.
-
Intermediate:
- Intermediates are fragments bound to an active site.
Fragment
Molecules and reactive species are both initialized as the Fragment object (based on ase.Atoms). Some examples are given below.
Molecules
Before to follow this guide, you need to load the following packages:
import matplotlib.pyplot as plt
from autoadsorbate import Fragment, Surface, docs_plot_conformers, get_marked_smiles, get_drop_snapped, docs_plot_sites, _example_config, construct_smiles
from ase.visualize.plot import plot_atoms
from ase import Atoms
Let us initialize a molecule of dimethyl ether (DME):
from autoadsorbate import Fragment
f = Fragment(input = 'COC', to_initialize = 5)
import matplotlib.pyplot as plt
from autoadsorbate import docs_plot_conformers
conformer_trajectory = f.conformers
fig = docs_plot_conformers(conformer_trajectory)
plt.show()

Notice that the orientation of the fragment is arbitrary. While we could simply place these structures onto the surface of a material, it would be difficult to evaluate the quality of these initial random configurations. This uncertainty would force us to sample a large number of structures and run dynamic simulations to explore local minima and determine which configurations are the most stable.
However, in the case of DME, we can leverage chemical intuition to simplify the problem. The oxygen atom bridging the two methyl groups has two lone electron pairs. By using a simple trick—replacing one of these lone pairs with a marker atom (such as chlorine, Cl)—we can guide the placement more effectively.
Notice that we had to make two adjustments to the SMILES string. To replace the lone pair with a marker atom, we must "trick" the valence of the oxygen atom and rearrange the SMILES formula so that the marker atom appears first (for easier bookkeeping).
- COC original
- CO(Cl)C add Cl instead of the O lone pair (this is an invalid SMILES)
- C[O+](Cl)C trick to make the valence work
- Cl[O+](C)C rearrange so that the SMILES string starts with the marker first (for easy book keeping)
This can be also done with a function:
from autoadsorbate import get_marked_smiles
marked_smile = get_marked_smiles(['COC'])[0]
marked_smile
'Cl[O+](C)(C)'
These surrogate smilles can now be used to initialize a Fragment object (we can set the number of randoms conformers to be initialized):
f = Fragment(input = 'Cl[O+](C)(C)', to_initialize = 5)
len(f.conformers)
5
We can visualize these structures:
conformer_trajectory = f.conformers
fig = docs_plot_conformers(conformer_trajectory)
plt.show()

Now we can use the marker atom to orient our molecule:
from autoadsorbate import docs_plot_sites
oriented_conformer_trajectory = [f.get_conformer(i) for i, _ in enumerate(f.conformers)]
fig = docs_plot_conformers(oriented_conformer_trajectory)
plt.show()

We can also easily remove the marker:
clean_conformer_trajectory = [atoms[1:] for atoms in oriented_conformer_trajectory]
fig = docs_plot_conformers(clean_conformer_trajectory)
plt.show()

Reactive species
Methoxy
f = Fragment(input = 'ClOC', to_initialize = 5)
oriented_conformer_trajectory = [f.get_conformer(i) for i, _ in enumerate(f.conformers)]
fig = docs_plot_conformers(oriented_conformer_trajectory)
plt.show()

Methyl
f = Fragment(input = 'ClC', to_initialize = 5)
oriented_conformer_trajectory = [f.get_conformer(i) for i, _ in enumerate(f.conformers)]
fig = docs_plot_conformers(oriented_conformer_trajectory)
plt.show()

Frangments with more than one binding mode (e.g. 1,2-PDO)
bound through single site:
f = Fragment(input = 'Cl[OH+]CC(O)C', to_initialize = 5)
oriented_conformer_trajectory = [f.get_conformer(i) for i, _ in enumerate(f.conformers)]
fig = docs_plot_conformers(oriented_conformer_trajectory)
plt.show()

Coordinated withboth hydroxil:
f = Fragment(input = 'S1S[OH+]CC([OH+]1)C', to_initialize = 5)
oriented_conformer_trajectory = [f.get_conformer(i) for i, _ in enumerate(f.conformers)]
fig = docs_plot_conformers(oriented_conformer_trajectory)
plt.show()

Generating fragments from your own atomic structure:
from ase.io import read
f = Fragment(input = read('my_slab_with_adsorbate.traj'), to_initialize = 5)
oriented_conformer_trajectory = [f.get_conformer(i) for i, _ in enumerate(f.conformers)]
Surface
Defining the surface of a slab may seem like a simple task, but different approaches can yield varying results depending on the context. When considering catalytic sites, we can define these as surface regions capable of binding a fragment. By using reasonable steric criteria—essentially asking, "Is there enough space for a molecule to bind to that site?"—we can identify all possible binding sites on the slab's surface. These sites can be classified as top, bridge, or multi-fold, depending on how many atoms surround the site.
As an example: First, we need to define a slab (any ase.Atoms object). A slab is an arrangement of atoms that represents the boundary between a material and another phase, such as gas, fluid, or another material. We can either read an existing slab, or a new slab:
from ase.build import fcc111
slab = fcc111('Cu', (4,4,4), periodic=True, vacuum=10)
Now we can initalize the Surface object which associates the constructed slab (ase.Atoms) with additional information required for placing Fragments. We can view which atoms are in the surface:
s = Surface(slab)
plot_atoms(s.view_surface(return_atoms=True))
Visualizing surface Cu atoms as Zn

We have access to all the sites info as a pandas dataframe:
s.site_df.head()
<div>
<style scoped>
.dataframe tbody tr th:only-of-type {
vertical-align: middle;
}
.dataframe tbody tr th {
vertical-align: top;
}
.dataframe thead th {
text-align: right;
}
</style>
<table border="1" class="dataframe">
<thead>
<tr style="text-align: right;">
<th></th>
<th>coordinates</th>
<th>connectivity</th>
<th>topology</th>
<th>n_vector</th>
<th>h_vector</th>
<th>site_formula</th>
</tr>
</thead>
<tbody>
<tr>
<th>0</th>
<td>[0.0, 0.0, 16.252703415323644]</td>
<td>1</td>
<td>[48]</td>
<td>[-0.004670396521231514, -0.00314499039640268