Surfgen Module¶
Submodules¶
surfgen.core module¶
-
surfgen.core.align_adsorbate_multiple_atoms(slab, first_layer, molecule, match_dict)[source]¶ This function can align adsorbate with multiple adsorption site on the surface with correct geometry.
Parameter:
- slab:
- Atoms Object. The slab that we want to study. It has to be the pure-surface
- first_layer:
- List. Indexes of atoms on the first layer.
- molecule:
- Atoms Object. The adsorbate
- match_dict:
Dictionary which contains the matching information between the molecule and the surface
Data Structure of match_dict: match_dict = {‘0’: {‘location’: a, ‘norm_vector’: b}}
Return:
The Atoms object of the molecule that aligned with the norm_vector
-
surfgen.core.align_adsorbate_one_atom(slab, first_layer, site, norm_vector, label, molecule, mol_index, hybridization)[source]¶ This function can align the molecule on particular site through norm_vector.
Parameter:
- slab:
- Atoms Object. The slab that we want to study. It has to be the pure-surface.
- first_layer:
- List of integers. the index of atoms in first layer of slab.
- site:
- Numpy Array (3x1). The location of adsorption site.
- norm_vector:
- Numpy Array (3x1). The norm_vector of the adsorption site (‘ontop’, ‘bridge’, ‘hollow’, ‘four_fold’)
- label:
- List of Strings. The index of atoms that construct the adsorption site
- molecule:
- Atoms Object. The adsorbate
- mol_index:
- Integer. The index of adsorbed atoms
- hybridization:
- String. ‘sp’, ‘sp2’, ‘sp3’. Different hybridzation means different bond length.
Return:
The Atoms object of the molecule that aligned with the norm_vector
-
surfgen.core.connection_matrix(slab, surface_atoms, bond_length)[source]¶ This function can help us determine the connection matrix in the slab, and it will return two things:
- The connection matrix (which contains the index information)
- All the coordinates (which contains all the information of coordinates)
Parameters:
- slab:
- The slab that we need
- surface_atoms:
- all the indexes of atoms on the surface
- Return:
Two dictionaries and one list.
Two dictionaries are:
- connector. It defines all the nearest neighbour of top surface layer atoms
- conn_coordinates. It defines all the coordinates of the atoms above
One list is:
- outer_atom_index. It defines all the index of atoms on the outside of our slab.
-
surfgen.core.find_all_ads_sites(slab, connector, conn_coordinates, surf_atoms, bond_length)[source]¶ This function can help us generate all possible adsorption sites on a surface.
Parameters:
- slab:
- slab that we want to find adsorption sites
- connector:
- the connection matrix of all atoms on top layer
- conn_coordiantes:
- coordinates of all atoms on top layer and its related atoms
- surf_atoms:
- index of all atoms on top layer
- bond_length:
- for determining the bridge site
- Return:
a dictionary which contains informations for ‘ontop’, ‘bridge’, ‘hollow’ and ‘four_fold’ sites
There are three pieces of information:
- location: contains the location of adsorption sites
- norm_vector: contains the norm_vector for each adsorption sites
- label: contains the index of atoms that construct the adsorption sites
-
surfgen.core.get_coordination_number(connector)[source]¶ This function can help us determine the coordination number of the top layer of the surface.
Parameters:
connector: the connector variable created by connection_matrix
- Return:
- A dictionary. {atom_index: coordination_number}
-
surfgen.core.slab_generator(bulk, miller_index, slab_height, vacuum, supercell)[source]¶ This function can generate slab from the ASE bulk object
Parameter:
- bulk:
- ASE Object. It is the bulk that we want to cut surface from
- miller_index:
- (3x1) tuple. It defines the miller index of the surface
- slab_height:
Real (in Angstrom). It defines how large your surface slab should be on z-direction
Note
You need to try slab_height first to generate the desired slab that you want to have. A role of thumb is that usually 10.0 is large enough for the calculations.
- vacuum:
- Real (in Angstrom). It defines the length of the vacuum layer
- supercell:
(3x1) tuple or (3x3) np.array. It defines the supercell of our surface slab
- (3x3) np.array is quite useful in generating surface like sqrt(3)*sqrt(3)
- Return:
- Atoms object, which is the slab that we want
-
surfgen.core.surf_atom_finder(slab, upper_lower='upper')[source]¶ This function is used to get the surface atom
Parameters:
- slab:
- the surface slab that you want, it must be clean surface (without adsorption)
- upper_lower:
string, to state whether you want the upper layer or lower layer or both.
e.g. ‘upper’, ‘lower’, ‘both’. In default it is ‘upper’
- ‘upper’: return indexes of atoms on upper surface layer
- ‘lower’: return indexes of atoms on lower surface layer
- ‘both’: return both of the surface layers
surfgen.math module¶
-
surfgen.math.matrix_multiple(matrix_list)[source]¶ This function can help to return the multiplication of a series of matrixes.
Parameters:
- matrix_list:
- The list of matrix that we want to multiply
- Return:
- Numpy array. The result that we want
-
surfgen.math.rotation_along_two_points(p, p1, p2, theta)[source]¶ This function can rotate a point (p) at an angle (theta) along a line generated by two point (p1, p2)
Parameters:
- p:
- Numpy Array (3x1). The point that we want to rotate
- p1:
- Numpy Array (3x1). First point of the line
- p2:
- Numpy Array (3x1). Second point of the line
- theta:
- Real. Angle that we want to rotate (0 - 360 degree)
- Return:
- Numpy Array (3x1). The point from rotation.
Reference:
The algorithm of this function is copied from two sources:
- https://robotics.stackexchange.com/questions/12782/how-rotate-a-point-around-an-arbitrary-line-in-3d
- http://paulbourke.net/geometry/rotate/
Many thanks to the authors
-
surfgen.math.rotation_matrix_from_vectors(vec1, vec2)[source]¶ Find the rotation matrix that aligns vec1 to vec2
This function is taken from Peter’s answer on stackoverflow: Question Link
Parameters:
- vec1:
- Numpy array. A 3d “source” vector
- vec2:
- Numpy array. A 3d “destination” vector
- Return:
- A transform matrix (3x3) which when applied to vec1, aligns it with vec2.
surfgen.util module¶
-
surfgen.util.average(v)[source]¶ This function can calculate the average of a list of vectors.
Parameter:
- v:
- A list of (3x1) np.array. It contains several vectors
- Return:
- (3x1) tuple. The averaged vector.
-
surfgen.util.check_4fold(v1, v2, v3, v4)[source]¶ This function is used to check whether the four atoms can form a square size or not
Parameters:
- v1:
- (3x1) np.array. Coordinates of the first point
- v2:
- (3x1) np.array. Coordinates of the second point
- v3:
- (3x1) np.array. Coordinates of the third point
- v4:
- (3x1) np.array. Coordinates of the fourth point
- Return:
- Boolean. True or False
-
surfgen.util.check_above(site, molecule)[source]¶ This function can check whether all atoms in molecule are above the adsorption site
Parameters:
- site:
- Numpy Array. The adsorption site that we want
- molecule:
- Atoms Object. The adsorbates.
- Return:
- Boolean. Shows whether the molecule is above the surface or not.
-
surfgen.util.check_distance(molecule, match_dict)[source]¶ This function can calculate the distance between the molecule and the adsorption sites.
Parameters:
- molecule:
- The molecule that we are interested in
- match_dict:
- The dictionary of the adsorbed information.
- Return:
- One real number. Shows the distance between molecule and adsorption sites.
-
surfgen.util.check_hollow(v1, v2, v3)[source]¶ This function can determine whether three points can form a 3-fold hollow site or not
Parameters:
- v1:
- (3x1) np.array. Coordinates of the first point
- v2:
- (3x1) np.array. Coordinates of the second point
- v3:
- (3x1) np.array. Coordinates of the third point
- Return:
- Boolean. True or False
-
surfgen.util.check_hybridization(slab, first_layer, norm_vector, molecule, mol_index, hybridization)[source]¶ This function can check whether the hybridization is confirmed or not
Parameters:
- slab:
- Atoms Object. The slab that we want to study. It has to be the pure-surface.
- first_layer:
- List of Integers. The indexes of first layer of slab.
- molecule:
- Atoms Object. The molecule that we are interested in
- mol_index:
- Integer. Index of adsorbed atom
- hybridization:
- String. ‘sp’, ‘sp2’, ‘sp3’
- Return:
- Atoms Object. With correct direction of hybridization
-
surfgen.util.check_match(molecule, match_dict, threshold=2)[source]¶ This function can check whether the molecule and the adsorption site are matched
Parameters:
- molecule:
- The molecule that
- match_dict:
- Dictionary that contains the match information of molecule and the adsorption sites.
- threshold:
- Real number. It shows how large the mismatch can be. The default number is 0.8
- Return:
- Boolean. Check whether the molecule and the adsorption sites are matched.
-
surfgen.util.dist_molecule_surface(slab, first_layer, molecule)[source]¶ This function can calculate the distance between the molecule and the slab
Parameters:
- slab:
- Atoms Object. The surface slab that we are interested in
- first_layer:
- List of integers. The indexes of the first layer of slab
- molecule:
- Atoms Object. The adsorbates that we want
- Return:
- shortest distance between the slab and the molecule
-
surfgen.util.find_center_molecule_adsSites(molecule, match_dict)[source]¶ This function can help us find the center of molecule and the center of the adsorption sites
Parameters:
- molecule:
- Atoms Object. The molecules we are interested in.
- match_dict:
- The matching dictionary of the molecule and the adsorption site
- Return:
Four Numpy Arraies (3x1). Contains:
- mol_center_location: the center of adsorbed atom in molecule
- mol_center_norm_vector_normalized: the averaged normalized norm vector of all adsorbed atom in molecule
- adsSite_center_location: the center of adsorption sites
- adsSites_center_norm_vector_normalized: the averaged normalized norm vector of all adsorption sites
-
surfgen.util.find_connection_atoms(molecule, mol_index, distance_range=[0.5, 1.8])[source]¶ This function can find all the closest atoms respect to the adsorbed atom. We only need to use this function when we only have one adsorbed atom. If there are more than one adsorbed atoms, then we only need to find the furthest atom.
Parameters:
- molecule:
- ASE Atoms Object. The molecule that we want to adsorb on the surface
- mol_index:
- Integer. The index of the adsorbed atom of molecules
- distange_range:
- List of Reals. Define the shortest and longest bond length in the molecule. Default is 0.5->1.8. [a,b]: a < b needed.
- Return:
- The function will retun a list of indexes in the molecule.