catkit.gen package

Submodules

catkit.gen.adsorption module

class catkit.gen.adsorption.AdsorptionSites(slab, surface_atoms=None, tol=1e-05)[source]

Adsorption site object.

ex_sites(index, select='inner', cutoff=0)[source]

Get site indices inside or outside of a cutoff radii from a provided periodic site index. If two sites are provided, an option to return the mutually inclusive points is also available.

get_adsorption_edges(symmetric=True, periodic=True)[source]

Return the edges of adsorption sties defined as all regions with adjacent vertices.

Parameters:
  • symmetric (bool) – Return only the symmetrically reduced edges.
  • periodic (bool) – Return edges which are unique via periodicity.
Returns:

edges – All edges crossing ridge or vertices indexed by the expanded unit slab.

Return type:

ndarray (n, 2)

get_adsorption_vectors(unique=True, screen=True)[source]

Returns the vectors representing the furthest distance from the neighboring atoms.

Returns:vectors – Adsorption vectors for surface sites.
Return type:ndarray (n, 3)
get_connectivity(unique=True)[source]

Return the number of connections associated with each site.

get_coordinates(unique=True)[source]

Return the 3D coordinates associated with each site.

get_periodic_sites(screen=True)[source]

Return an index of the coordinates which are unique by periodic boundary conditions.

Parameters:screen (bool) – Return only sites inside the unit cell.
Returns:periodic_match – Indices of the coordinates which are identical by periodic boundary conditions.
Return type:ndarray (n,)
get_symmetric_sites(unique=True, screen=True)[source]

Determine the symmetrically unique adsorption sites from a list of fractional coordinates.

Parameters:
  • unique (bool) – Return only the unique symmetrically reduced sites.
  • screen (bool) – Return only sites inside the unit cell.
Returns:

sites – Dictionary of sites containing index of site

Return type:

dict of lists

get_topology(unique=True)[source]

Return the indices of adjacent surface atoms.

plot(savefile=None)[source]

Create a visualization of the sites.

class catkit.gen.adsorption.Builder(slab, surface_atoms=None, tol=1e-05)[source]

Bases: catkit.gen.adsorption.AdsorptionSites

Initial module for creation of 3D slab structures with attached adsorbates.

add_adsorbate(adsorbate, bonds=None, index=0, auto_construct=True, **kwargs)[source]

Add and adsorbate to a slab. If the auto_constructor flag is False, the atoms object provided will be attached at the active site.

Parameters:
  • adsorbate (gratoms object) – Molecule to connect to the surface.
  • bonds (int or list of 2 int) – Index of adsorbate atoms to be bonded.
  • index (int) – Index of the site or edge to use as the adsorption position. A value of -1 will return all possible structures.
  • auto_construct (bool) – Whether to automatically estimate the position of atoms in larger molecules or use the provided structure.
Returns:

slabs – Slab(s) with adsorbate attached.

Return type:

gratoms object

add_adsorbates(adsorbates, indices)[source]
catkit.gen.adsorption.get_adsorption_sites(slab, surface_atoms=None, symmetry_reduced=True, tol=1e-05)[source]

Get the adsorption sites of a slab as defined by surface symmetries of the surface atoms.

Parameters:
  • slab (Atoms object) – The slab to find adsorption sites for. Must have surface atoms defined.
  • surface_atoms (array_like (N,)) – List of the surface atom indices.
  • symmetry_reduced (bool) – Return the symmetrically unique sites only.
  • adsorption_vectors (bool) – Return the adsorption vectors.
Returns:

  • coordinates (ndarray (N, 3)) – Cartesian coordinates of activate sites.
  • connectivity (ndarray (N,)) – Number of bonds formed for a given adsorbate.
  • symmetry_index (ndarray (N,)) – Arbitrary indices of symmetric sites.

catkit.gen.adsorption.symmetry_equivalent_points(fractional_coordinates, atoms, tol=1e-05)[source]

Return the symmetrically equivalent points from a list of provided fractional coordinates.

Parameters:
  • fractional_coordinates (ndarray (N ,3)) – Fractional coordinates to find symmetrical equivalence between.
  • atoms (Atoms object) – Atoms object to use the unit cell, positions, and pbc of.
  • tol (float) – Float point precision tolerance.
Returns:

symmetry_match – Indices of fractional coordinates which are unique or matching.

Return type:

ndarray (N,)

catkit.gen.geometry module

catkit.gen.molecules module

catkit.gen.molecules.bin_hydrogen(hydrogens=1, bins=1)[source]

Recursive function for determining distributions of hydrogens across bins.

catkit.gen.molecules.get_3D_positions(atoms, bond_index=None)[source]

Return an estimation of the 3D structure of a Gratoms object based on its graph.

WARNING: This function operates on the atoms object in-place.

Parameters:
  • atoms (Gratoms object) – Structure with connectivity matrix to provide a 3D structure.
  • bond_index (int) – Index of the atoms to consider as the origin of a surface bonding site.
Returns:

atoms – Structure with updated 3D positions.

Return type:

Gratoms object

catkit.gen.molecules.get_topologies(symbols, saturate=False)[source]

Return the possible topologies of a given chemical species.

Parameters:
  • symbols (str) – Atomic symbols to construct the topologies from.
  • saturate (bool) – Saturate the molecule with hydrogen based on the default.radicals set.
Returns:

molecules – Gratoms objects with unique connectivity matrix attached. No 3D positions will be provided for these structures.

Return type:

list (N,)

catkit.gen.molecules.hydrogenate(atoms, bins, copy=True)[source]

Add hydrogens to a gratoms object via provided bins

catkit.gen.pathways module

class catkit.gen.pathways.ReactionNetwork(db_name='reaction-network.db', base_valence=None, nbond_limits=None)[source]

A class for accessing a temporary SQLite database. This function works as a context manager and should be used as follows:

with ReactionNetwork() as rn:
(Perform operation here)

This syntax will automatically construct the temporary database, or access an existing one. Upon exiting the indentation, the changes to the database will be automatically committed.

create_table()[source]

Create the SQLite database table framework.

load_3d_structures(ids=None)[source]

Return Gratoms objects from the ReactionNetwork database.

Parameters:ids (int or list of int) – Identifier of the molecule in the database. If None, return all structure.
Returns:images – All Gratoms objects in the database.
Return type:list
load_molecules(ids=None, binned=False)[source]

Load 2D molecule graphs from the database.

Parameters:binned (bool) – Return the molecules in sub-dictionaries of their corresponding composition and bonding tags.
Returns:molecules – All molecules present in the database.
Return type:dict
load_pathways(broken_bonds=False)[source]

Save enumerated pathways the ReactionNetwork database. More than two reactants or two products is not supported.

Parameters:broken_bonds (bool) – Return the index information of which bond was broken. Only supported for elementary steps.
Returns:pathways – All pathways present in the database.
Return type:list

Return the enumeration of molecules which can be produced from the specified atoms.

Parameters:
  • element_pool (dict) – Atomic symbols keys paired with the maximum number of that atom.
  • load_molecules (bool) – Load any existing molecules from the database.
  • multiple_bond_search (bool) – Allow atoms to form bonds with other atoms in the molecule.

Search for reaction mechanisms from a pre-populated database of molecules. By default, only single bond addition pathways are enumerated (Also called elementary steps).

Parameters:
  • reconfiguration (bool) – Search for reconfiguration paths. Reconfiguration paths are all those where only the bond order is changed. R1 –> P1.
  • substitution (bool) – Search for substitution paths. Substitution paths are all those where one bond is broken and one bond is formed simultaneously. R1 + R2 –> P1 + P2.
plot_reaction_network(file_name=None)[source]

Plot the reaction network present in the database.

save_3d_structure(gratoms, overwrite=False)[source]

Save Cartesian coordinates into the ReactionNetwork database.

Parameters:
  • gratoms (Atoms object) – Structure containing Cartesian coordinates to be saved.
  • overwrite (bool) – Allow the database to overwrite a matching index.
save_molecules(molecules)[source]

Save enumerated molecules to the ReactionNetwork database.

Parameters:molecules (dict) – Molecules to be saved to the database.
save_pathways(pathways, broken_bonds=None)[source]

Save enumerated pathways the ReactionNetwork database. More than two reactants or two products is not supported.

Parameters:
  • pathways (list) – Sorted pathways in the form [R1, R2, P1, P2].
  • broken_bonds (list) – Comma separated strings of index associated with the two atoms whos bond is broken. List order must match pathways.

catkit.gen.route module

catkit.gen.route.get_heppel_sellers(nu, terminal)[source]

Returns an array of linearly independent reaction routes as described by Heppel-Sellers reaction route enumeration.

Parameters:
  • nu (ndarray (n, m)) – The stoichiometric matrix of n species and m mechanisms.
  • terminal (ndarray (j,)) – Indices of the m species to be considered as terminal
Returns:

sigma – Linearly independent set of Heppel-Sellers reaction routes.

Return type:

ndarray (m, k)

catkit.gen.route.get_reaction_routes(nu, sigma, empty_routes=True, independent_only=False)[source]

Returns an array of reaction routes. Returns all full reaction routes by default.

Parameters:
  • nu (ndarray (n, m)) – The stoichiometric matrix of n species and m mechanisms.
  • sigma (ndarray (m, j)) – A linearly independent set of reaction routes.
  • empty_routes (bool) – Return the empty routes along with the full routes.
  • independent_only (bool) – Return only a linearly independent set of full reaction routes. Can take less time.
Returns:

  • FR (ndarray (m, k)) – Enumerated full reaction routes.
  • ER (ndarray (m, l)) – Enumerated empty reaction routes.

catkit.gen.route.get_response_reactions(epsilon, selection=None, species=False)[source]

Returns an array of possible response reaction routes for a given chemical formula array.

Parameters:
  • epsilon (ndarray (n, m)) – Chemical formula array of n elements by m molecular species.
  • selection (ndarray (j,)) – Indices of the m species to be considered as terminal
  • species (bool) – Return the indices of the chemical species used.
Returns:

  • RER (ndarray (m, k)) – Possible response reaction routes.
  • index (ndarray (j, k)) – Indices of the k terminal species use to produce the l response reactions.

catkit.gen.surface module

class catkit.gen.surface.SlabGenerator(bulk, miller_index, layers, vacuum=None, fixed=None, layer_type='ang', attach_graph=True, standardize_bulk=False, primitive=True, tol=1e-08)[source]

Bases: object

Class for generation of slab unit cells from bulk unit cells.

Many surface operations rely upon / are made easier through the bulk basis cell they are created from. The SlabGenerator class is designed to house these operations.

Return the miller indices associated with the users requested values. Follows the following steps:

  • Convert Miller-Bravais notation into standard Miller index.
  • (optional) Ensure the bulk cell is in its standard form.
  • Convert the indices to the cell for the primitive lattice.
  • Reduce the indices by their greatest common divisor.
Parameters:
  • bulk (Atoms object) – Bulk system to be converted into slab.
  • miller_index (list (3,) or (4,)) – Miller index to construct surface from. If length 4, Miller-Bravais notation is assumed.
  • layers (int) – Number of layers to include in the slab. A slab layer is defined as a unique z-coordinate.
  • vacuum (float) – Angstroms of vacuum to apply to the slab.
  • fixed (int) – Number of slab layers to constrain.
  • layer_type ('angs', 'trim', 'stoich', or 'sym') –

    Determines how to perform slab layering.

    ’angs’: Layers denotes the thickness of the slab in Angstroms. ‘trim’: The slab will be trimmed to a number of layers equal to the exact number of unique z-coordinates. Useful for precision control. ‘stoich’ : Constraints any slab generated to have the same stoichiometric ratio as the provided bulk. ‘sym’ : Return a slab which is inversion symmetric. i.e. The same on both sides.

  • attach_graph (bool) – Attach the connectivity graph generated from the bulk structure. This is only necessary for fingerprinting and setting it to False can save time. Surface atoms will be found regardless.
  • standardize_bulk (bool) – Covert the bulk input to its standard form before and produce the cleave from it. This is highly recommended as Miller indices are not defined for non-standard cells.
  • tol (float) – Tolerance for floating point rounding errors.
adsorption_sites(slab, **kwargs)[source]

Helper function to return the adsorption sites of the provided slab.

Parameters:slab (atoms object) – The slab to find adsorption sites for. Assumes you are using the same basis.
Returns:output – Coordinates and connectivity of the adsorption sites. The symmetry indices can also be returned.
Return type:tuple (n, n) | (n, n, n)
align_crystal(bulk, miller_index)[source]

Return an aligned unit cell from bulk unit cell. This alignment rotates the a and b basis vectors to be parallel to the Miller index.

Parameters:
  • bulk (Atoms object) – Bulk system to be standardized.
  • miller_index (list (3,)) – Miller indices to align with the basis vectors.
Returns:

new_bulk – Standardized bulk unit cell.

Return type:

Gratoms object

get_slab(size=1, iterm=0)[source]

Generate a slab from the bulk structure. This function is meant specifically for selection of an individual termination or enumeration through surfaces of various size.

This function will orthogonalize the c basis vector and align it with the z-axis which breaks bulk symmetry along the z-axis.

Parameters:
  • size (int, array_like (2,) or (2, 2)) – Size of the unit cell to create as described in set_size().
  • iterm (int) – A termination index in reference to the list of possible terminations.
Returns:

slab – The modified basis slab produced based on the layer specifications given.

Return type:

Gratoms object

get_slab_basis(iterm=0, maxn=20)[source]

Return a list of all terminations which have been properly shifted and with an appropriate number of layers added. This function is mainly for performance, to prevent looping over other operations which are not related the size of the slab.

This step also contains periodically constrained orthogonalization of the c basis. This implementation only works if the a anb b basis vectors are properly aligned with the x and y axis. This is strictly to assist the correct identification of surface atoms.

Only produces the terminations requested as a lazy evaluator.

Parameters:
  • iterm (int) – Index of the slab termination to return.
  • maxn (int) – The maximum integer component to search for a more orthogonal bulk.
Returns:

ibasis – Prepared, ith basis.

Return type:

Gratoms object

get_unique_terminations()[source]

Determine the fractional coordinate shift that will result in a unique surface termination. This is not required if bulk standardization has been performed, since all available z shifts will result in a unique termination for a primitive cell.

Returns:unique_shift – Fractional coordinate shifts which will result in unique terminations.
Return type:array (n,)
make_symmetric(slab)[source]

Returns a symmetric slab. Note, this will trim the slab potentially resulting in loss of stoichiometry.

set_size(slab, size)[source]

Set the size of a slab based one of three methods.

1. An integer value performs a search of valid matrix operations to perform on the ab-basis vectors to return a set which with a minimal sum of distances and an angle closest to 90 degrees.

2. An array_like of length 2 will multiply the existing basis vectors by that amount.

3. An array of shape (2, 2) will be interpreted as matrix notation to multiply the ab-basis vectors by.

Parameters:
  • slab (Atoms object) – Slab to be made into the requested size.
  • size (int, array_like (2,) or (2, 2)) – Size of the unit cell to create as described above.
Returns:

supercell – Supercell of the requested size.

Return type:

Gratoms object

catkit.gen.surface.convert_miller_index(miller_index, atoms1, atoms2)[source]

Return a converted miller index between two atoms objects.

catkit.gen.surface.generate_indices(max_index)[source]

Return an array of miller indices enumerated up to values plus or minus some maximum. Filters out lists with greatest common divisors greater than one. Only positive values need to be considered for the first index.

Parameters:max_index (int) – Maximum number that will be considered for a given surface.
Returns:unique_index – Unique miller indices
Return type:ndarray (n, 3)
catkit.gen.surface.get_degenerate_indices(bulk, miller_index)[source]

Return the miller indices which are degenerate to a given miller index for a particular bulk structure.

Parameters:
  • bulk (Atoms object) – Bulk structure to get the degenerate miller indices.
  • miller_index (array_like (3,)) – Miller index to get the degenerate indices for.
Returns:

degenerate_indices – Degenerate miller indices to the provided index.

Return type:

array (N, 3)

catkit.gen.surface.get_unique_indices(bulk, max_index)[source]

Returns an array of miller indices which will produce unique surface terminations based on a provided bulk structure.

Parameters:
  • bulk (Atoms object) – Bulk structure to get the unique miller indices.
  • max_index (int) – Maximum number that will be considered for a given surface.
Returns:

unique_millers – Symmetrically distinct miller indices for a given bulk.

Return type:

ndarray (n, 3)

catkit.gen.surface.transform_ab(slab, matrix, tol=1e-05)[source]

Transform the slab basis vectors parallel to the z-plane by matrix notation. This can result in changing the slabs cell size. This can also result in very unusual slab dimensions, use with caution.

Parameters:
  • slab (Atoms object) – The slab to be transformed.
  • matrix (array_like (2, 2)) – The matrix notation transformation of the a and b basis vectors.
  • tol (float) – Float point precision tolerance.
Returns:

slab – Slab after transformation.

Return type:

Atoms object

catkit.gen.utils module

catkit.gen.utils.get_voronoi_neighbors(atoms, cutoff=5.0)[source]

Return the connectivity matrix from the Voronoi method. Multi-bonding occurs through periodic boundary conditions.

Parameters:
  • atoms (atoms object) – Atoms object with the periodic boundary conditions and unit cell information to use.
  • cutoff (float) – Radius of maximum atomic bond distance to consider.
Returns:

connectivity – Number of edges formed between atoms in a system.

Return type:

ndarray (n, n)

catkit.gen.utils.get_cutoff_neighbors(atoms, cutoff=None, atol=1e-08)[source]

Return the connectivity matrix from a simple radial cutoff. Multi-bonding occurs through periodic boundary conditions.

Parameters:
  • atoms (atoms object) – Atoms object with the periodic boundary conditions and unit cell information to use.
  • cutoff (float) – Cutoff radius to use when determining neighbors.
  • atol (float) – Absolute tolerance to use when computing distances.
Returns:

connectivity – Number of edges formed between atoms in a system.

Return type:

ndarray (n, n)

catkit.gen.utils.get_integer_enumeration(N=3, span=[0, 2])[source]

Return the enumerated array of a span of integer values. These enumerations are limited to the length N.

For the default span of [0, 2], the enumeration equates to the corners of an N-dimensional hypercube.

Parameters:
  • N (int) – Length of enumerated lists to produce.
  • span (list | slice) – The range of integers to be considered for enumeration.
Returns:

enumeration – Enumeration of the requested integers.

Return type:

ndarray (M, N)

catkit.gen.utils.trilaterate(centers, r, zvector=None)[source]

Find the intersection of two or three spheres. In the case of two sphere intersection, the z-coordinate is assumed to be an intersection of a plane whose normal is aligned with the points and perpendicular to the positive z-coordinate.

If more than three spheres are supplied, the centroid of the points is returned (no radii used).

Parameters:
  • centers (list | ndarray (n, 3)) – Cartesian coordinates representing the center of each sphere
  • r (list | ndarray (n,)) – The radii of the spheres.
  • zvector (ndarray (3,)) – The vector associated with the upward direction for under-specified coordinations (1 and 2).
Returns:

intersection – The point where all spheres/planes intersect.

Return type:

ndarray (3,)

catkit.gen.utils.get_unique_xy(xyz_coords, cutoff=0.1)[source]

Return the unique coordinates of an atoms object for the requrested atoms indices. Z-coordinates are projected to maximum z-coordinate by default.

Parameters:
  • xyz_coords (ndarray (n, 3)) – Cartesian coordinates to identify unique xy positions from.
  • cutoff (float) – Distance in Angstrons to consider xy-coordinate unique within.
Returns:

xy_pos – Unique xy coordinates projected onto a maximal z coordinate.

Return type:

ndarray (m, 3)

catkit.gen.utils.expand_cell(atoms, cutoff=None, padding=None)[source]

Return Cartesian coordinates atoms within a supercell which contains repetitions of the unit cell which contains at least one neighboring atom.

Parameters:
  • atoms (Atoms object) – Atoms with the periodic boundary conditions and unit cell information to use.
  • cutoff (float) – Radius of maximum atomic bond distance to consider.
  • padding (ndarray (3,)) – Padding of repetition of the unit cell in the x, y, z directions. e.g. [1, 0, 1].
Returns:

  • index (ndarray (N,)) – Indices associated with the original unit cell positions.
  • coords (ndarray (N, 3)) – Cartesian coordinates associated with positions in the supercell.
  • offsets (ndarray (M, 3)) – Integer offsets of each unit cell.

catkit.gen.utils.matching_sites(position, comparators, tol=1e-08)[source]

Get the indices of all points in a comparator list that are equal to a given position (with a tolerance), taking into account periodic boundary conditions (adaptation from Pymatgen).

This will only accept a fractional coordinate scheme.

Parameters:
  • position (list (3,)) – Fractional coordinate to compare to list.
  • comparators (list (3, n)) – Fractional coordinates to compare against.
  • tol (float) – Absolute tolerance.
Returns:

match – Indices of matches.

Return type:

list (n,)

catkit.gen.utils.get_basis_vectors(coordinates)[source]

Return a set of basis vectors for a given array of 3D coordinates.

Parameters:coordinates (array_like (3, 3) | (2, 3)) – Cartesian coordinates to determine the basis of. If only 2 positions are given 3rd is chosen as the positive y-axis.
Returns:basis_vectors – Automatically generated basis vectors from the given positions.
Return type:ndarray (3, 3)
catkit.gen.utils.connectivity_to_edges(connectivity, indices=None)[source]

Convert a Numpy connectivity matrix into a list of NetworkX compatible edges.

catkit.gen.utils.isomorphic_molecules(graph0, graph1)[source]

Check whether two molecule graphs are isomorphic.

catkit.gen.utils.matching_coordinates(position, comparators, tol=1e-08)[source]

Get the indices of all points in a comparator list that are equal to a given position (with a tolerance), taking into account periodic boundary conditions (adaptation from Pymatgen).

This will only accept a Cartesian coordinate scheme. TODO: merge this with matching_sites.

Parameters:
  • position (list (3,)) – Fractional coordinate to compare to list.
  • comparators (list (3, N)) – Fractional coordinates to compare against.
  • tol (float) – Absolute tolerance.
Returns:

match – Indices of matches.

Return type:

list (N,)

catkit.gen.utils.get_unique_coordinates(atoms, axis=2, tag=False, tol=0.001)[source]

Return unique coordinate values of a given atoms object for a specified axis.

Parameters:
  • atoms (object) – Atoms object to search for unique values along.
  • axis (int (0, 1, or 2)) – Look for unique values along the x, y, or z axis.
  • tag (bool) – Assign ASE-like tags to each layer of the slab.
  • tol (float) – The tolerance to search for unique values within.
Returns:

values – Array of unique positions in fractional coordinates.

Return type:

ndarray (n,)

catkit.gen.utils.get_reciprocal_vectors(atoms)[source]

Return the reciprocal lattice vectors to a atoms unit cell.

catkit.gen.utils.plane_normal(xyz)[source]

Return the surface normal vector to a plane of best fit.

Parameters:xyz (ndarray (n, 3)) – 3D points to fit plane to.
Returns:vec – Unit vector normal to the plane of best fit.
Return type:ndarray (1, 3)
catkit.gen.utils.running_mean(array, N=5)[source]

Calculate the running mean of array for N instances.

Parameters:
  • array (array_like | ndarray (N,)) – Array of values to have a average taken from.
  • N (int) – Number of values to take an average with.
Returns:

running_mean – Mean value of the running average.

Return type:

ndarray (N + 1,)

catkit.gen.utils.to_gratoms(atoms, edges=None)[source]

Convert and atom object to a gratoms object.

catkit.gen.utils.get_atomic_numbers(formula, return_count=False)[source]

Return the atomic numbers associated with a chemical formula.

Parameters:
  • formula (string) – A chemical formula to parse into atomic numbers.
  • return_count (bool) – Return the count of each element in the formula.
Returns:

  • numbers (ndarray (n,)) – Element numbers in associated species.
  • counts (ndarray (n,)) – Count of each element in a species.

catkit.gen.utils.get_reference_energies(species, energies)[source]

Get reference energies for the elements in a set of molecules.

Parameters:
  • species (list (n,)) – Chemical formulas for each molecular species.
  • energies (list (n,)) – Total energies associated with each species.
Returns:

  • elements (ndarray (n,)) – Atomic elements associated with all species.
  • references (ndarray (n,)) – Reference energies associated with each element.

catkit.gen.utils.parse_slice(slice_name)[source]

Return a correctly parsed slice from input of varying types.

catkit.gen.utils.ext_gcd(a, b)[source]

Extension of greatest common divisor.

catkit.gen.utils.list_gcd(values)[source]

Return the greatest common divisor of a list of values.

Module contents

Catalysis Generator.