-
Notifications
You must be signed in to change notification settings - Fork 874
/
functional_groups.py
364 lines (285 loc) · 13.7 KB
/
functional_groups.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
"""Determine functional groups present in a Molecule."""
from __future__ import annotations
import copy
from pymatgen.analysis.graphs import MoleculeGraph
from pymatgen.analysis.local_env import OpenBabelNN
from pymatgen.core.structure import Molecule
from pymatgen.io.babel import BabelMolAdaptor
try:
import networkx as nx
import networkx.algorithms.isomorphism as iso
except ImportError:
raise ImportError("pymatgen.analysis.functional_groups requires the NetworkX graph library to be installed.")
__author__ = "Evan Spotte-Smith"
__version__ = "0.1"
__maintainer__ = "Evan Spotte-Smith"
__email__ = "[email protected]"
__status__ = "Beta"
__date__ = "July 2018"
__credit__ = "Peiyuan Yu"
class FunctionalGroupExtractor:
"""
This class is used to algorithmically parse a molecule (represented by an
instance of pymatgen.analysis.graphs.MoleculeGraph) and determine arbitrary
functional groups.
"""
def __init__(self, molecule, optimize=False):
"""
Instantiation method for FunctionalGroupExtractor.
Args:
molecule: Either a filename, a pymatgen.core.structure.Molecule
object, or a pymatgen.analysis.graphs.MoleculeGraph object.
optimize: Default False. If True, then the input molecule will be
modified, adding Hydrogens, performing a simple conformer search, etc.
"""
self.molgraph = None
if isinstance(molecule, str):
try:
if optimize:
ob_mol = BabelMolAdaptor.from_file(molecule, file_format="mol")
# OBMolecule does not contain pymatgen Molecule information
# So, we need to wrap the ob_mol in a BabelMolAdapter
ob_mol.add_hydrogen()
ob_mol.make3d()
ob_mol.localopt()
self.molecule = ob_mol.pymatgen_mol
else:
self.molecule = Molecule.from_file(molecule)
except OSError:
raise ValueError("Input must be a valid molecule file, a Molecule object, or a MoleculeGraph object.")
elif isinstance(molecule, Molecule):
if optimize:
ob_mol = BabelMolAdaptor(molecule)
ob_mol.add_hydrogen()
ob_mol.make3d()
ob_mol.localopt()
self.molecule = ob_mol.pymatgen_mol
else:
self.molecule = molecule
elif isinstance(molecule, MoleculeGraph):
if optimize:
ob_mol = BabelMolAdaptor(molecule.molecule)
ob_mol.add_hydrogen()
ob_mol.make3d()
ob_mol.localopt()
self.molecule = ob_mol.pymatgen_mol
else:
self.molecule = molecule.molecule
self.molgraph = molecule
else:
raise TypeError("Input to FunctionalGroupExtractor must be str, Molecule, or MoleculeGraph.")
if self.molgraph is None:
self.molgraph = MoleculeGraph.from_local_env_strategy(self.molecule, OpenBabelNN())
# Assign a specie and coordinates to each node in the graph,
# corresponding to the Site in the Molecule object
self.molgraph.set_node_attributes()
self.species = nx.get_node_attributes(self.molgraph.graph, "specie")
def get_heteroatoms(self, elements=None):
"""
Identify non-H, non-C atoms in the MoleculeGraph, returning a list of
their node indices.
Args:
elements: List of elements to identify (if only certain
functional groups are of interest).
Returns:
set of ints representing node indices
"""
hetero_atoms = set()
for node in self.molgraph.graph.nodes():
if elements is not None:
if str(self.species[node]) in elements:
hetero_atoms.add(node)
elif str(self.species[node]) not in ["C", "H"]:
hetero_atoms.add(node)
return hetero_atoms
def get_special_carbon(self, elements=None):
"""
Identify Carbon atoms in the MoleculeGraph that fit the characteristics
defined Ertl (2017), returning a list of their node indices.
The conditions for marking carbon atoms are (quoted from Ertl):
"- atoms connected by non-aromatic double or triple bond to any
heteroatom
- atoms in nonaromatic carbon-carbon double or triple bonds
- acetal carbons, i.e. sp3 carbons connected to two or more oxygens,
nitrogens or sulfurs; these O, N or S atoms must have only single bonds
- all atoms in oxirane, aziridine and thiirane rings"
Args:
elements: List of elements that will qualify a carbon as special
(if only certain functional groups are of interest).
Default None.
Returns:
set of ints representing node indices
"""
specials = set()
# For this function, only carbons are considered
carbons = [n for n in self.molgraph.graph.nodes if str(self.species[n]) == "C"]
# Condition one: double/triple bonds to heteroatoms
for node in carbons:
neighbors = self.molgraph.graph[node]
for neighbor, attributes in neighbors.items():
if elements is not None:
if str(self.species[neighbor]) in elements and int(attributes[0]["weight"]) in [2, 3]:
specials.add(node)
elif str(self.species[neighbor]) not in ["C", "H"] and int(attributes[0]["weight"]) in [2, 3]:
specials.add(node)
# Condition two: carbon-carbon double & triple bonds
for node in carbons:
neighbors = self.molgraph.graph[node]
for neighbor, attributes in neighbors.items():
if str(self.species[neighbor]) == "C" and int(attributes[0]["weight"]) in [2, 3]:
specials.add(node)
specials.add(neighbor)
# Condition three: Acetal carbons
for node in carbons:
neighbors = self.molgraph.graph[node]
neighbor_spec = [str(self.species[n]) for n in neighbors]
ons = len([n for n in neighbor_spec if n in ["O", "N", "S"]])
if len(neighbors) == 4 and ons >= 2:
specials.add(node)
# Condition four: oxirane/aziridine/thiirane rings
rings = self.molgraph.find_rings()
rings_indices = [set(sum(ring, ())) for ring in rings]
for ring in rings_indices:
ring_spec = sorted(str(self.species[node]) for node in ring)
# All rings of interest are three-member rings
if len(ring) == 3 and ring_spec in [
["C", "C", "O"],
["C", "C", "N"],
["C", "C", "S"],
]:
for node in ring:
if node in carbons:
specials.add(node)
return specials
def link_marked_atoms(self, atoms):
"""
Take a list of marked "interesting" atoms (heteroatoms, special carbons)
and attempt to connect them, returning a list of disjoint groups of
special atoms (and their connected hydrogens).
Args:
atoms: set of marked "interesting" atoms, presumably identified
using other functions in this class.
Returns:
list of sets of ints, representing groups of connected atoms
"""
# We will add hydrogens to functional groups
hydrogens = {n for n in self.molgraph.graph.nodes if str(self.species[n]) == "H"}
# Graph representation of only marked atoms
subgraph = self.molgraph.graph.subgraph(list(atoms)).to_undirected()
func_groups = []
for func_grp in nx.connected_components(subgraph):
grp_hs = set()
for node in func_grp:
neighbors = self.molgraph.graph[node]
for neighbor in neighbors:
# Add all associated hydrogens into the functional group
if neighbor in hydrogens:
grp_hs.add(neighbor)
func_grp = func_grp | grp_hs
func_groups.append(func_grp)
return func_groups
def get_basic_functional_groups(self, func_groups=None):
"""
Identify functional groups that cannot be identified by the Ertl method
of get_special_carbon and get_heteroatoms, such as benzene rings, methyl
groups, and ethyl groups.
TODO: Think of other functional groups that are important enough to be
added (ex: do we need ethyl, butyl, propyl?)
Args:
func_groups: List of strs representing the functional groups of
interest. Default to None, meaning that all of the functional groups
defined in this function will be sought.
Returns:
list of sets of ints, representing groups of connected atoms
"""
strat = OpenBabelNN()
hydrogens = {n for n in self.molgraph.graph.nodes if str(self.species[n]) == "H"}
carbons = [n for n in self.molgraph.graph.nodes if str(self.species[n]) == "C"]
if func_groups is None:
func_groups = ["methyl", "phenyl"]
results = []
if "methyl" in func_groups:
for node in carbons:
neighbors = strat.get_nn_info(self.molecule, node)
hs = {n["site_index"] for n in neighbors if n["site_index"] in hydrogens}
# Methyl group is CH3, but this will also catch methane
if len(hs) >= 3:
hs.add(node)
results.append(hs)
if "phenyl" in func_groups:
rings_indices = [set(sum(ring, ())) for ring in self.molgraph.find_rings()]
possible_phenyl = [r for r in rings_indices if len(r) == 6]
for ring in possible_phenyl:
# Phenyl group should have only one (0 for benzene) member whose
# neighbors are not two carbons and one hydrogen
num_deviants = 0
for node in ring:
neighbors = strat.get_nn_info(self.molecule, node)
neighbor_spec = sorted(str(self.species[n["site_index"]]) for n in neighbors)
if neighbor_spec != ["C", "C", "H"]:
num_deviants += 1
if num_deviants <= 1:
ring_group = []
for node in ring:
ring_group = copy.deepcopy(ring)
neighbors = self.molgraph.graph[node]
# Add hydrogens to the functional group
for neighbor in neighbors:
if neighbor in hydrogens:
ring_group.add(neighbor)
results.append(ring_group)
return results
def get_all_functional_groups(self, elements=None, func_groups=None, catch_basic=True):
"""
Identify all functional groups (or all within a certain subset) in the
molecule, combining the methods described above.
Args:
elements: List of elements that will qualify a carbon as special
(if only certain functional groups are of interest).
Default None.
func_groups: List of strs representing the functional groups of
interest. Default to None, meaning that all of the functional groups
defined in this function will be sought.
catch_basic: bool. If True, use get_basic_functional_groups and
other methods
Returns:
list of sets of ints, representing groups of connected atoms
"""
heteroatoms = self.get_heteroatoms(elements=elements)
special_cs = self.get_special_carbon(elements=elements)
groups = self.link_marked_atoms(heteroatoms | special_cs)
if catch_basic:
groups += self.get_basic_functional_groups(func_groups=func_groups)
return groups
def categorize_functional_groups(self, groups):
"""
Determine classes of functional groups present in a set.
Args:
groups: Set of functional groups.
Returns:
dict containing representations of the groups, the indices of
where the group occurs in the MoleculeGraph, and how many of each
type of group there is.
"""
categories = {}
em = iso.numerical_edge_match("weight", 1)
nm = iso.categorical_node_match("specie", "C")
for group in groups:
atoms = [self.molecule[a] for a in group]
species = [a.specie for a in atoms]
coords = [a.coords for a in atoms]
adaptor = BabelMolAdaptor(Molecule(species, coords))
# Use Canonical SMILES to ensure uniqueness
smiles = adaptor.pybel_mol.write("can").strip()
if smiles in categories:
this_subgraph = self.molgraph.graph.subgraph(list(group)).to_undirected()
for other in categories[smiles]["groups"]:
other_subgraph = self.molgraph.graph.subgraph(list(other)).to_undirected()
if not nx.is_isomorphic(this_subgraph, other_subgraph, edge_match=em, node_match=nm):
break
if group not in categories[smiles]["groups"]:
categories[smiles]["groups"].append(group)
categories[smiles]["count"] += 1
else:
categories[smiles] = {"groups": [group], "count": 1}
return categories