carmm.build package

Subpackages

Submodules

carmm.build.adsorbate_placer module

class carmm.build.adsorbate_placer.RotationBox(atoms_ads, atoms_site, ads_idx, site_idx, bond_length, neighb_idx=0, lps=1, lp_idx=1, cutoff_mult=1)

Bases: object

Object intended to store functions for rotation

Args:
atoms_ads: Atoms object

Contains positions of adsorbate.

atoms_site: Atoms object

Contains positions of the adsorption site.

ads_idx: integer

Index of the adsorbate atom bonding adsorbate and adsorption site.

site_idx: integer

Index of the site atom bonding adsorbate and adsorption site.

bond_length: float

Length of the desired bond between the site and adsorbate.

neighb_idx: integer

Index of neighbour used to form principle axes of rotation.

lps: integer

Number of lone pairs on the site atom.

lp_idx: integer

Index of the potential lone pair site (i.e., if two lone pair sites are present, controls which site the atom adsorbs to).

cutoff_mult: float

Constant factor multiplies species dependent cut-off radii used in carmm.analysis.neighbours to find the first nearest neighbour to the adsorption site (i.e., larger value potentially finds more neighbours).

Returns:
RotationBox: RotationBox Object

Object used to adsorb a given adsorbate to a given site and rotate as desired using internal functions.

find_adsorbate_normal(atoms, index)

Returns a normalised vector specifying the direction of the desired bonding atom of the adsorbate, and the mean of the cartesian positions of the adsorbate molecule.

Parameters:
  • atoms – Atoms object Contains the atomic position of the adsorbate.

  • index – integer The atomic index of the bonding atom.

Returns:

numpy array, (3).

Contains the normalised direction of bonding atom->center of positions about (0,0,0).

Return type:

site_normal

find_adsorbate_rotation_axes()

Find the x, y and z rotation axes, defined by the direction of the bond to the adsorbate (site_normal) and the direction of the bond to a given pre-existing neighbour atom. If only one around the adsorbate site is found, the next atom its neighbours are used.

Parameters:
  • atoms_site – Atoms object Contains positions of the adsorption site.

  • site_idx – integer Index of the site atom bonding adsorbate and adsorption site.

  • neighb_idx – integer Index of neighbour used to form principle axes of rotation.

  • cutoff_mult – float Multiplier for the cutoff radii of neighbouring atoms (1 = natural cutoff radii).

Returns:

numpy array, (3)

x rotation axis (perpendicular to site bonds to the adsorbate and another neighbour atom).

y_axis, numpy array, (3)

y rotation axis (perpendicular given x and y rotation axes).

z_axis, numpy array, (3)

z rotation axis (defined by the site normal).

Return type:

x_axis

find_generic_normal(atoms, index1, index2)

Returns a normalised vector specifying the direction of a generic bond.

Parameters:
  • atoms – Atoms object Contains the atomic position of the adsorbate.

  • index1 – integer The atomic index of atom 1.

  • index2 – integer The atomic index of atom 2.

Returns:

numpy array, (3).

Contains the normalised direction of atom1 to atom2. about (0,0,0).

Return type:

site_normal

find_site_normal(atoms, index, cutoff_mult)

Returns a normalised vector specifying the direction of the adsorbate-adsorption site bond.

Parameters:
  • atoms – Atoms object Contains the atomic position of the site.

  • index – integer The atomic index of the desired site atom.

  • cutoff_mult – float Multiplier for the cutoff radii of neighbouring atoms (1 = natural cutoff radii).

Returns:

numpy array, (3).

Contains the normalised direction of the ads+site bond about (0,0,0).

Return type:

site_normal

normal_rotation_matrix(theta, ax)

Returns a (3x3) rotation matrix for a rotation about a given axis by a theta value.

Parameters:
  • theta – float Angle of rotation in radians.

  • ax – numpy array, size = (3) Direction of the rotation axis, centered on [0, 0, 0]

Returns:

numpy array, size = (3x3)

Rotation matrix to be applied to atomic positions.

Return type:

rotation_matrix

place_adsorbate()

Place and rotates the adsorbate along a bond length defined by vector. Molecule rotated such that the center of positions for the molecule points along the bonding direction. Equivalent to rotate_and_place_adsorbate(…,rotation=[0,0,0]).

Parameters:
  • atoms_ads – Atoms object Contains positions of adsorbate.

  • atoms_site – Atoms object Contains positions of the adsorption site.

  • ads_idx – integer Index of the adsorbate atom bonding adsorbate and adsorption site.

  • site_idx – integer Index of the site atom bonding adsorbate and adsorption site.

  • bond_length – float Length of the bond for adsorption.

Returns:

numpy array, (:,3)

Array of positions of the adsorbed species and site.

atom_ads: numpy array, (:,3)

Array of positions for the rotated adsorbed species only.

Return type:

ads_and_site

rotate(rotation)
rotate_ads2site_vec(atoms, site_vec, mol_vec)

Performs a rotation of the molecule onto the site normal. ROTATION MATRIX TAKEN FROM https://gist.github.com/kevinmoran/b45980723e53edeb8a5a43c49f134724

Parameters:
  • atoms – Atoms object Contains positions of adsorbate.

  • atoms_idx – integer Index of the adsorbate atom bonding adsorbate and adsorption site.

  • site_vec – numpy array, size=(3) vector defining the direction

  • mol_vec – numpy array, size=(3) vector defining direction from the bonding atom to the mean of the atomic positions.

Returns:

numpy array, size=(Nx3)

Position of the atoms rotated so site_vec and mol_vec align.

Return type:

rotated_positions

rotate_and_place_adsorbate(rotation=[0, 0, 0])

Place and rotates the adsorbate along a bond length defined by vector - vector calculated assuming VSEPR model with one remaining lone pair. Rotation takes places along the bonding direction. See documentation for rotation variable for directions of rotation.

Parameters:
  • atoms_ads – Atoms object Contains positions of adsorbate.

  • atoms_site – Atoms object Contains positions of the adsorption site.

  • bond_length – float Length of the bond for adsorption.

  • ads_idx – integer Index of the adsorbate atom bonding adsorbate and adsorption site.

  • site_idx – integer Index of the site atom bonding adsorbate and adsorption site.

  • neighb_idx – integer Index of neighbour used to form principle axes of rotation.

  • rotation – list Angles in degrees of rotation for absorbate. z - along ads and site bond, x - perpendicular to neighbour and bond, y - perpendicular to x and z.

Returns:

numpy array, (:,3)

Array of positions of the adsorbed species and site.

atom_ads: numpy array, (:,3)

Array of positions for the rotated adsorbed species only.

Return type:

ads_and_site

carmm.build.adsorbate_placer_gui module

carmm.build.adsorbate_placer_gui.get_rotations()
carmm.build.adsorbate_placer_gui.main_ase_gui(gui)
carmm.build.adsorbate_placer_gui.print_position()
carmm.build.adsorbate_placer_gui.rotationbox_gui(RotationObj, output='output.xyz')

A graphical view of the adsorbate-site system, allowing you to visually rotate the adsorbate in place and produce a structure file with that configuration. Begin the visualisation by pressing “Start ASE”, rotate the adsorbate using the x, y and z sliders, and write the structure to the desired filename using the “Print Position” button.

Parameters:
  • RotationObj – RotationBox object The desired adsorbate-site system as a RotationBox() object from build/adsorbate_placer.py

  • output – string The desired filename to print the structure to (include desired file extension)

Returns:

On press of “Print Position” button, file containing Atoms object

carmm.build.adsorbate_placer_gui.start_ase_gui()

carmm.build.alloy module

carmm.build.cutout module

A collection of scripts that can generate and manipulate Chemshell fragment and cluster objects. These can be created from ASE Atoms objects and can use either a large, non-periodic structure to cut out a spherical cluster or can make a cluster from a periodic model.

carmm.build.cutout.cutout_sphere(atoms, centre, distance_cutoff=5.0)

Returns a spherical cutout of a structure TODO: This doesn’t work with periodic boundary conditions. It’d be nice to get that working and update the regression to test this.

Parameters:

atoms: Atoms object

Input structure to cutout from

centre: Integer

Index of central atom in cutout

distance_cutoff: Float

Distance outside of which atoms are removed

carmm.build.cutout.transpose(periodic, cluster, start, stop, centre_periodic, centre_cluster, file_name)

Returns a ChemShell cluster representation of a periodic model

Parameters: periodic: Atoms object

Input periodic model of interest

cluster: Atoms object

Input cluster model of interest

start: Integer

Starting index of atoms to copy into cluster

stop: Integer

End index of atoms to copy into cluster

centre_periodic: Integer

Centre of periodic model

centre_cluster: Integer

Centre of cluster model

name: String

Name of new model to save

carmm.build.facets module

carmm.build.facets.generate(bulk_model, layers=2, facets=[(1, 1, 1)], supercell=(1, 1, 1), vacuum=10, save=False)

Function to create variety of different facets

Parameters:

bulk_model: Atoms object or string

specify the unit cell of the bulk geometry, or element from which to construct default bulk

layers: int

specify the repeating layers along the direction defined in ‘axis’ of the ‘center’ function

facets: list of array of ints

specify the list of facets required as output

supercell: array of ints

specify the periodicity of the supercell; (x,y,z). # can this just be an int? - No, has to be an array #

vacuum: float

specifies the vacuum addition to the slab in the direction of ‘axis-’ of the ‘center’ function

carmm.build.slab_consistent_bulk_generator module

carmm.build.slab_consistent_bulk_generator.bulk_identifier(slab, cutoff_distance=10.0)

— 08/01/2024 testing ASE version 3.24.0 —

Slab models can be generated using various libraries like ASE, pymatgen, atomman, etc. However, it becomes important to check how these slab models are generated within various code bases.

Pymatgen for example generates desired surface facets by transforming the basis vectors of the conventional bulk unit cell such that the c-vector of new bulk is normal to the desired surface plane. The basis transformation however changes the convergence behavior of the slab models with respect to the thickness of the slab. Hence, a new bulk unit cell is necessary, generated for each individual slab model, to make sure that the surface energies are calculated using a correct bulk reference.

This function helps to identify the repeating unit in the z-direction which ultimately represents the bulk atoms in a slab structure. This new bulk structure will be consistent with the slab model and the computed bulk energy for this structure will help us to obtain surface energies that are convergent with increasing slab thickness.

NOTE: Please use this functionality only if you observe divergent surface energy behaviour using the already available total energy of the conventional bulk unit cell.

Detailed information is available in the following literature article. Schultz, Peter A. “First-principles calculations of metal surfaces. I. Slab-consistent bulk reference for convergent surface properties.” Physical Review B 103.19 (2021): 195426.

How this functionality works:

1. the difference between the x,y and z coordinates of an atom and its corresponding periodic image will be a linear combination of the cell vectors a,b and c.

2. This idea is extended here in case of a slab model where the linear combination check is considered only w.r.t to the a and b vectors

Args:

  • slab: ASE Atoms object representing the slab structure

  • cutoff_distance: Cutoff distance for identifying neighboring atoms (default is 10.0).

This will be used to

perform further checks on atoms which satisfy the condition of linear combination.

Returns:

  • ASE Atoms object representing the new bulk structure

carmm.build.slab_consistent_bulk_generator.is_close_to_integer(arr, tolerance=1e-05)

Check if elements in the array are close to integers within a given tolerance. This will be used to check if the x- and y-coordinate of an atom is a linear combination the ‘a’ and ‘b’ cell vector of the slab model

Args: - arr: Numpy array - tolerance: Tolerance for closeness to integers (default is 1e-5)

Returns: - Boolean array indicating if elements are close to integers

carmm.build.unwrap module

Module contents