API Reference

MOF Deconstructions

mofstructure.mofdeconstructor.all_ferrocene_metals(ase_atom, graph)[source]

A function to find metals corresponding to ferrocene. These metals should not be considered during mof-constructions

parameters:

ase_atom: ASE atom graph : dictionary containing neigbour lists

returns:

list_of_metals: list of indices of ferrocene metals

mofstructure.mofdeconstructor.angle_tolerance_to_rad(angle, tolerance=5)[source]

Convert angle from degrees to radians with tolerance.

parameters:
angle: float

Angle in degrees.

tolerance: float

Tolerance in degrees, default is 5.

returns:

bond_angle_rad: float Adjusted angle in radians.

mofstructure.mofdeconstructor.ase_2_pybel(atoms)[source]

As simple script to convert from ase atom object to pybel. There are many functionalities like atom typing, bond typing and autmatic addition of hydrogen that can be performed on a pybel molecule object, which can not be directly performed on an ase atom object. E.g

add_hygrogen = pybel.addh()

remove_hydrogen = pybel.removeh()

https://openbabel.org/docs/dev/UseTheLibrary/Python_Pybel.html

parameters:

atoms : ASE atoms object

returns:

pybel: pybel molecular object.

mofstructure.mofdeconstructor.ase_2_xyz(atoms)[source]

Create an xyz string from an ase atom object to be compatible with pybel in order to perfom some cheminformatics.

parameters:

atoms : ASE atoms object

returns:

a_str : string block of the atom object in xyz format

mofstructure.mofdeconstructor.check_planarity(p1, p2, p3, p4)[source]

A simple procedure to check whether a point is planar to three other points. Important to distinguish porphyrin type metals

parameters:

p1, p2, p3, p4 : ndarray containing x,y,z values

returns:

Boolean True: planar False: noneplanar

mofstructure.mofdeconstructor.compute_ase_neighbour(ase_atom)[source]

Create a connectivity graph using ASE neigbour list.

parameters:

ASE atoms

returns:
  1. atom_neighbors: A python dictionary, wherein each atom index is key and the value are the indices of it neigbours. e.g. atom_neighbors ={0:[1,2,3,4], 1:[3,4,5]…}

  2. matrix: An adjacency matrix that wherein each row correspond to to an atom index and the colums correspond to the interaction between that atom to the other atoms. The entries in the matrix are 1 or 0. 1 implies bonded and 0 implies not bonded.

mofstructure.mofdeconstructor.compute_cheminformatic_from_rdkit(ase_atom)[source]

A function that converts and ase atom into rkdit mol to extract some cheminformatic information.The script begins by taking an ase atom an creating an xyz string of the molecule, which is stored to memory. The xyz string is then read using the MolFromXYZBlock function in rdkit, which reads an xyz string of the molecule. However, this does not include any of the bonding information. To do this the following commands are parsed to the rdkit mol,bond_moll = Chem.Mol(rdmol) rdDetermineBonds.DetermineBonds(bond_moll,charge=0) https://github.com/rdkit/UGM_2022/blob/main/Notebooks/Landrum_WhatsNew.ipynb

parameters:

ASE atoms

returns:

smile string, inchi, inchikey

mofstructure.mofdeconstructor.compute_inchis(obmol)[source]

A function to compute the cheminformatic inchi and inchikey using openbabel. The inchi and inchikey are IUPAC cheminformatic identifiers of molecules. They are important to easily search for molecules in databases. For MOFs, this is quite important to search for different secondary building units (sbu) and ligands found in the MOF. More about inchi can be found in the following link

https://iupac.org/who-we-are/divisions/division-details/inchi

parameters:

obmol: openbabel molecule object

returns:

inChi, inChiKey: iupac inchi and inchikeys for thesame molecule.

mofstructure.mofdeconstructor.compute_openbabel_cheminformatic(ase_atom)[source]

A function the returns all smiles, inchi and inchikeys from an ase_atom object. The function starts by wrapping the system, which is important for systems in which minimum image convertion has made atoms appear to be uncoordinated at random positions. The wrap functions coordinates all the atoms together.

parameters:

atoms : ASE atoms object

returns:

smi, inChi, inChiKey

mofstructure.mofdeconstructor.compute_smi(obmol)[source]

A function to compute the SMILES (Simplified molecular-input line-entry system) notation of a molecule using openbabel. The SMILES is a line notation method that uses ASCII strings to represent molecules. More about SMILES can be found in the following link https://en.wikipedia.org/wiki/Simplified_molecular-input_line-entry_system

parameters:

obmol : openbabel molecule object

returns:

smi: SMILES notatation of the molecule.

mofstructure.mofdeconstructor.connected_components(graph)[source]

Find the connected fragments in a graph. Should work for any graph defined as a dictionary

parameters:

A graph in the form of dictionary e.g.:: graph = {1:[0,1,3], 2:[2,4,5]}

returns:

returns a python list of list of connected components These correspond to individual molecular fragments. list_of_connected_components = [[1,2],[1,3,4]]

mofstructure.mofdeconstructor.covalent_radius(element)[source]

A function that returns the covalent radius from the chemical symbol of an element.

parameters:

element: chemical symbol of an atom type.string

returns:

covalent_radii: covalent radius of the elements. type.float

mofstructure.mofdeconstructor.dfsutil_graph_method(graph, temp, node, visited)[source]

Depth-first search graph algorithm for traversing graph data structures. I starts at the root ‘node’ and explores as far as possible along each branch before backtracking. It is used here a a util for searching connected components in the MOF graph

parameters:

graph: any python dictionary temp: a python list to hold nodes that have been visited node: a key in the python dictionary (graph), which is used as the starting or root node visited: python list containing nodes that have been traversed.

returns:

python dictionary

mofstructure.mofdeconstructor.find_COS(ase_atom, graph)[source]

A simple aglorimth to search for COS.

O
|
-C
|
S
parameters:

ase_atom: ASE atom

returns:

dictionary of key = carbon index and values = sulphur index

mofstructure.mofdeconstructor.find_carbonyl_sulphate(ase_atom, graph)[source]

A simple aglorimth to search for Carbonyl sulphate found in the system.

   O
   |
-C-S
   |
   O
parameters:

ase_atom: ASE atom

returns:

dictionary of key = carbon index and values = oxygen index

mofstructure.mofdeconstructor.find_carboxylates(ase_atom, graph)[source]

A simple aglorimth to search for carboxylates found in the system.

parameters:
  • ase_atom: ASE atom

returns:
  • dictionary of key = carbon index and values = oxygen index

mofstructure.mofdeconstructor.find_key_or_value(key_or_value, dictionary)[source]

Search for a key or value in a dictionary. If the key or value is found, returns the corresponding value or key.

parameters:

key_or_value : int or str dictionary : dict

returns:

The corresponding value or key if found, otherwise None.

mofstructure.mofdeconstructor.find_phosphate(ase_atom, graph)[source]

A simple algorithm to search for Carbonyl sulphate found in the system.

 O
 |
-P-o
 |
 O
parameters:

ase_atom: ASE atom

returns

dictionary of key = carbon index and values = oxygen index

mofstructure.mofdeconstructor.find_phosphite(ase_atom, graph)[source]

A simple aglorimth to search for sulfides.

 P
 |
-C
 |
 P
parameters:

ase_atom: ASE atom

returns:

dictionary of key = carbon index and values = phosphorous index

mofstructure.mofdeconstructor.find_sulfides(ase_atom, graph)[source]

A simple aglorimth to search for sulfides.

 S
 |
-C
 |
 S
parameters:

ase_atom: ASE atom

returns:

dictionary of key = carbon index and values = sulphur index

mofstructure.mofdeconstructor.find_unique_building_units(list_of_connected_components, atom_pairs_at_breaking_point, ase_atom, porphyrin_checker, all_regions, wrap_system=True, cheminfo=False, add_dummy=False)[source]

A function that identifies and processes unique building units within a MOF. This function deconstructs a MOF into its constituent building unit and identifies unique building units. It can also 1. wrap the system into the unit cell, 2. add dummy atoms to neutralize the building units, 3. compute cheminformatics data using Open Babel.

parameters:

list_of_connected_components : list of connected components atom_pairs_at_breaking_point : dictionary of atom pairs at breaking point ase_atom : ASE atoms object porphyrin_checker : list of indices of porphyrins all_regions : dictionary of regions wrap_system : boolean, default True cheminfo : boolean, default False add_dummy: boolean, default False add_dummy keyword enables the addition of dumnmy atoms which can then be replaced by hydrogen to neutralize the building blocks.

returns:

mof_metal: list of metal indices mof_linker: list of linker indices concentration: dictionary of concentration of building units

mofstructure.mofdeconstructor.inter_atomic_distance_check(ase_atom)[source]

A function that checks whether two atoms are within a distance 1.0 Amstrong unless it is an R-H bond

parameters:

ase_atom : ASE atoms object

returns

boolean :

mofstructure.mofdeconstructor.is_ferrocene(metal_sbu, graph)[source]

A simple script to check whether a metal_sbu is ferrocene

parameter:

metal_sbu : ase_atom graph : dictionary containing neigbour lists

returns:

bool : True if the metal_sbu is ferrocene, False otherwise

mofstructure.mofdeconstructor.is_irmof(ase_atom, graph)[source]

returns True if the atom is part of a IRMOF motif

parameter:

ase_atom : ase_atom graph : dictionary containing neigbour lists

returns:

bool : True if the metal_sbu is part of a IRMOF motif, False otherwise

mofstructure.mofdeconstructor.is_mof32(ase_atom, graph)[source]

returns True if the atom is part of a MOF32 motif

parameter:

ase_atom : ase_atom graph : dictionary containing neigbour lists

returns:

bool : True if the metal_sbu is part of a MOF32 motif, False otherwise

mofstructure.mofdeconstructor.is_paddlewheel(metal_sbu, graph)[source]

returns True if the atom is part of a paddlewheel motif

parameter:

metal_sbu : ase_atom graph : dictionary containing neigbour lists

returns:

bool : True if the metal_sbu is part of a paddlewheel motif, False otherwise

mofstructure.mofdeconstructor.is_paddlewheel_with_water(ase_atom, graph)[source]

returns True if the atom is part of a paddle wheel with water motif

parameter:

ase_atom : ase_atom graph : dictionary containing neigbour lists

returns:

bool : True if the metal_sbu is part of a paddle wheel with water motif, False otherwise

mofstructure.mofdeconstructor.is_rodlike(metal_sbu)[source]

Simple test to check whether a metal sbu is a rodlike MOF

parameter:

metal_sbu : ase_atom

returns:

bool : True if the metal sbu is a rodlike MOF, False otherwise

mofstructure.mofdeconstructor.is_uio66(ase_atom, graph)[source]

returns True if the atom is part of a UIO66 motif

parameter:

ase_atom : ase_atom graph : dictionary containing neigbour lists

returns:

bool : True if the metal_sbu is part of a UIO66 motif, False otherwise

mofstructure.mofdeconstructor.ligands_and_metal_clusters(ase_atom)[source]

Start by checking whether there are more than 2 layers if yes, select one Here we select the largest connected component

parameters:

ase_atom: ASE atom

returns:

list_of_connected_components: list of connected components, in which each list contains atom indices atom_pairs_at_breaking_point: Dictionary containing point of disconnection Porpyrin_checker: Boolean showing whether the metal is in the centre of a porpherin Regions: Dictionary of regions.

mofstructure.mofdeconstructor.longest_list(lst)[source]

return longest list in list of list

mofstructure.mofdeconstructor.matrix2dict(bond_matrix)[source]

A simple procedure to convert an adjacency matrix into a python dictionary.

parameters:

bond matrix : adjacency matrix, type: nxn ndarray

returns:

graph: python dictionary

mofstructure.mofdeconstructor.max_index(lists)[source]

Extract index of list with the maximum element

mofstructure.mofdeconstructor.metal_coordination_enviroment(ase_atom)[source]

Find the enviroment of a metal

paramters:

ase_atom : ASE atoms object

returns:

metal_elt : list of metal elements metal_enviroment : dict where key is metal element and value is list of neighboring elements

mofstructure.mofdeconstructor.metal_coordination_number(ase_atom)[source]

Extract coordination number of central metal

paramters:

ase_atom : ASE atoms object

returns:

metal_elt : list of metal elements

mofstructure.mofdeconstructor.metal_in_porphyrin(ase_atom, graph)[source]

A funtion to check whether a metal is found at the centre of a phorphirin. The function specifically identifies metal atoms that are coordinated with four nitrogen atoms, as typically seen in porphyrin complexes. These atoms must also satisfy the planarity condition of the porphyrin structure.

https://en.wikipedia.org/wiki/Transition_metal_porphyrin_complexes

parameters:

ase_atom: ASE atom graph: python dictionary containing neigbours

returns:

list of indices consisting of index of metal atoms found in the ASE atom

mofstructure.mofdeconstructor.metal_in_porphyrin2(ase_atom, graph)[source]

https://en.wikipedia.org/wiki/Transition_metal_porphyrin_complexes Check whether a metal is found at the centre of a phorphirin

parameters:

ase_atom: ASE atom graph: python dictionary containing neigbours

returns:

list of indices consisting of index of metal atoms found in the ASE atom

mofstructure.mofdeconstructor.mof_regions(ase_atom, list_of_connected_components, atom_pairs_at_breaking_point)[source]

A function to map all atom indices to exact position in which the find themselves in the MOF. This function is used to partition a MOF into regions that correspond to unique unique building units.

parameters:

ase_atom: ASE atom list_of_connected_components : list of list, wherein each list correspond to atom indices of a specific building unit atom_pairs_at_breaking_point: dictionary containing pairs of atoms from which the bonds were broken

returns:

Move an index from any position in the list to the front The function is important to set the cell of a rodmof to point in the a-axis. Such that the system can be grow along this axis

mofstructure.mofdeconstructor.move2front(index_value, coords)[source]

Move an index from any position in the list to the front The function is important to set the cell of a rodmof to point in thea-axis. Such that the system can be grow along this axis

parameters:
  • index_value: index of item to move to the from

  • coords: list of coords

returns:

list of coords

mofstructure.mofdeconstructor.obmol_2_rdkit(obmol)[source]

A simple function to convert openbabel mol to RDkit mol A function that takes an openbabel molecule object and converts it to an RDkit molecule object. The importance of this function lies in the fact that there is no direct were to convert from an ase atom type to an rdkit molecule object. Consequently, the easiest approach it is firt convert the system to an openbabel molecule object using the function ase_2_pybel to convert to the pybel molecule object. Then using the following code obmol = pybel_mol.OBMol to obtain the openbabel molecule object that can then be used to convert to the rdkit molecule object.

parameters:

obmol : openbabel molecule object

returns:

rdmol : rdkit molecule object.

mofstructure.mofdeconstructor.remove_unbound_guest(ase_atom)[source]

A simple script to remove guest from a metal organic framework. 1) It begins by computing a connected graph component of all the fragments in the system using ASE neighbour list. 2) Secondly it selects indices of connected components which contain a metal 3) if the there are two or more components, we create a pytmagen graph for each components and filter out all components that are not polymeric 4) If there are two or more polymeric components, we check whether these systems there are identical or different and select all polymeric components

parameters:

ASE atoms

returns:

mof_indices: indices of the guest-free system. The guest-free ase_atom object can be obtained as follows; E.g. guest_free_system = ase_atom[mof_indices]

mofstructure.mofdeconstructor.remove_unbound_guest_and_return_unique(ase_atom)[source]

A simple script to remove guest from a metal organic framework. 1)It begins by computing a connected graph component of all the fragments in the system using ASE neighbour list. 2) Secondly it selects indicies of connected components which contain a metal 3)if the there are two or more components, we create a pytmagen graph for each components and filter out all components that are not polymeric 4) If there are two or more polymeric components, we check wether these systems there are identical or different and select only unique polymeric components

parameters:

ASE atoms

returns:

mof_indices : indinces of the guest free system. The guest free ase_atom object can be obtain as follows; E.g. guest_free_system = ase_atom[mof_indices]

mofstructure.mofdeconstructor.rod_manipulation(ase_atom, checker)[source]

Script to adjust Rodlike sbus. 1) Its collects the axis responsible for expanding the rod 2) It shifts all coordinates to the axis 3) It rotates the rod to lie in the directions of expansion,

parameter:

ase_atom : ASE atoms object checker : list of indices of atoms responsible for rod expansion

returns:

ase_atom : ASE atoms object with adjusted coordinates

mofstructure.mofdeconstructor.secondary_building_units(ase_atom)[source]
  1. Search for all carboxylate that are connected to a metal. Cut at the position between the carboxylate carbon and the adjecent carbon.

  2. Find all Nitrogen connected to metal. Check whether the nitrogen is in the centre of a porphirin ring. If no, cut at nitrogen metal bond.

  3. Look for oxygen that is connected to metal and two carbon. cut at metal oxygen bond

parameters:

ase_atom: ASE atom

returns:

list_of_connected_components: list of connected components,in which each list contains atom indices atom_pairs_at_breaking_point : Dictionary containing point of disconnection Porphyrin_checker: Boolean showing whether the metal is in the centre of a porpherin Regions: Dictionary of regions.

mofstructure.mofdeconstructor.transition_metals()[source]

A function that returns the list of all chemical symbols of metallic elements present in the periodic table

mofstructure.mofdeconstructor.wrap_systems_in_unit_cell(ase_atom, max_iter=30, skin=0.3)[source]

A simple aglorithm to reconnnect all atoms wrapped in a periodic boundary condition such that all atoms outside the box will appear reconnected.

parameters:

ase_atom : ASE atoms object max_iter : int, optional (default=30) Maximum number of iterations for reconnection skin : float, optional (default=0.3) Skin distance for bond reconnection

returns:
  • ase_atom : ASE atoms object with reconnected atoms wrapped in a periodic boundary condition

Computing Porosity

mofstructure.porosity.ase_to_zeoobject(ase_atom)[source]

Converts an ase atom type to a zeo++ Cssr object In zeo++ the xyz coordinate system is rotated to a zyx format.

parameter:

ase_atom: ase atom object

returns:

cssr_object: string representing zeo++ cssr object

mofstructure.porosity.zeo_calculation(ase_atom, probe_radius=1.86, number_of_steps=10000, high_accuracy=True, rad_file=None)[source]

Main script to compute geometric structure of porous systems. The focus here is on MOF, but the script can run on any porous periodic system. The script computes the accesible surface area, accessible volume and the pore geometry. There are many more outputs which can be extracted from ,vol_str and sa_str. Moreover there are also other computation that can be done. Check out the test directory in dependencies/pyzeo/test. Else contact bafgreat@gmail.com. if you need more output and can’t figure it out.

parameter:

ase_atom: ase atom object probe_radius: The radius of the probe. Here 1.86 is used as default number_of_steps: Number of GCMC simulation cycles high_accuracy: key to determine where to perform high accuracy computation

return:

python dictionary containing 1) AV_Volume_fraction: Accessible volume void fraction 2) AV_A^3: Accessible volume in A^2 3) AV_cm^3/g: Accessible volume in cm^3/g. This values is often infinity because it literatly divides the given value by the avogadro’s number 4) ASA_A^2: Accessible surface area A^2 5) ASA_m^2/cm^3: Accessible surface area in m^2/cm^3 6) Number_of_channels: Number of channels present in the porous system, which correspond to the number of pores within the system 7) LCD_A: The largest cavity diameter is the largest sphere that can be inserted in a porous system without overlapping with any of the atoms in the system. 8) lfpd_A:The largest included sphere along free sphere path is largest sphere that can be inserted in the pore 9)PLD_A:The pore limiting diameter is the largest sphere that can freely diffuse through the porous network without overlapping with any of the atoms in the system

Computing stacking in 2D systems and inter-layer height

mofstructure.cof_stacking.compute_cof_stacking(ase_atom)[source]

A a simple function to compute the stacking pattern of COFs or layered materials like graphene

parameter:

ase_atom : ASE Atoms object

returns

layers : list of list wherei each list correspond to a layar lateral_offsets : list of list where each list contains the lateral offsets between two layers interlayer_height : list of list where each list contains the interlayer heights between two layers