File: EnumerateStereoisomers.py

package info (click to toggle)
rdkit 202009.4-1
  • links: PTS, VCS
  • area: main
  • in suites: bullseye
  • size: 129,624 kB
  • sloc: cpp: 288,030; python: 75,571; java: 6,999; ansic: 5,481; sql: 1,968; yacc: 1,842; lex: 1,254; makefile: 572; javascript: 461; xml: 229; fortran: 183; sh: 134; cs: 93
file content (353 lines) | stat: -rw-r--r-- 13,527 bytes parent folder | download | duplicates (2)
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
import random
from rdkit import Chem
from rdkit.Chem.rdDistGeom import EmbedMolecule


class StereoEnumerationOptions(object):
  """
          - tryEmbedding: if set the process attempts to generate a standard RDKit distance geometry
            conformation for the stereisomer. If this fails, we assume that the stereoisomer is
            non-physical and don't return it. NOTE that this is computationally expensive and is
            just a heuristic that could result in stereoisomers being lost.

          - onlyUnassigned: if set (the default), stereocenters which have specified stereochemistry
            will not be perturbed unless they are part of a relative stereo
            group.

          - maxIsomers: the maximum number of isomers to yield, if the
            number of possible isomers is greater than maxIsomers, a
            random subset will be yielded. If 0, all isomers are
            yielded. Since every additional stereo center doubles the
            number of results (and execution time) it's important to
            keep an eye on this.

          - onlyStereoGroups: Only find stereoisomers that differ at the
            StereoGroups associated with the molecule.
    """
  __slots__ = ('tryEmbedding', 'onlyUnassigned', 'onlyStereoGroups', 'maxIsomers', 'rand', 'unique')

  def __init__(self, tryEmbedding=False, onlyUnassigned=True, maxIsomers=1024, rand=None,
               unique=True, onlyStereoGroups=False):
    self.tryEmbedding = tryEmbedding
    self.onlyUnassigned = onlyUnassigned
    self.onlyStereoGroups = onlyStereoGroups
    self.maxIsomers = maxIsomers
    self.rand = rand
    self.unique = unique


class _BondFlipper(object):

  def __init__(self, bond):
    self.bond = bond

  def flip(self, flag):
    if flag:
      self.bond.SetStereo(Chem.BondStereo.STEREOCIS)
    else:
      self.bond.SetStereo(Chem.BondStereo.STEREOTRANS)


class _AtomFlipper(object):

  def __init__(self, atom):
    self.atom = atom

  def flip(self, flag):
    if flag:
      self.atom.SetChiralTag(Chem.ChiralType.CHI_TETRAHEDRAL_CW)
    else:
      self.atom.SetChiralTag(Chem.ChiralType.CHI_TETRAHEDRAL_CCW)


class _StereoGroupFlipper(object):

  def __init__(self, group):
    self._original_parities = [(a, a.GetChiralTag()) for a in group.GetAtoms()]

  def flip(self, flag):
    if flag:
      for a, original_parity in self._original_parities:
        a.SetChiralTag(original_parity)
    else:
      for a, original_parity in self._original_parities:
        if original_parity == Chem.ChiralType.CHI_TETRAHEDRAL_CW:
          a.SetChiralTag(Chem.ChiralType.CHI_TETRAHEDRAL_CCW)
        elif original_parity == Chem.ChiralType.CHI_TETRAHEDRAL_CCW:
          a.SetChiralTag(Chem.ChiralType.CHI_TETRAHEDRAL_CW)


def _getFlippers(mol, options):
  Chem.FindPotentialStereoBonds(mol)
  flippers = []
  if not options.onlyStereoGroups:
    for atom in mol.GetAtoms():
      if atom.HasProp("_ChiralityPossible"):
        if (not options.onlyUnassigned or atom.GetChiralTag() == Chem.ChiralType.CHI_UNSPECIFIED):
          flippers.append(_AtomFlipper(atom))

    for bond in mol.GetBonds():
      bstereo = bond.GetStereo()
      if bstereo != Chem.BondStereo.STEREONONE:
        if (not options.onlyUnassigned or bstereo == Chem.BondStereo.STEREOANY):
          flippers.append(_BondFlipper(bond))

  if options.onlyUnassigned:
    # otherwise these will be counted twice
    for group in mol.GetStereoGroups():
      if group.GetGroupType() != Chem.StereoGroupType.STEREO_ABSOLUTE:
        flippers.append(_StereoGroupFlipper(group))

  return flippers


class _RangeBitsGenerator(object):

  def __init__(self, nCenters):
    self.nCenters = nCenters

  def __iter__(self):
    for val in range(2**self.nCenters):
      yield val


class _UniqueRandomBitsGenerator(object):

  def __init__(self, nCenters, maxIsomers, rand):
    self.nCenters = nCenters
    self.maxIsomers = maxIsomers
    self.rand = rand
    self.already_seen = set()

  def __iter__(self):
    # note: important that this is not 'while True' otherwise it
    # would be possible to have an infinite loop caused by all
    # isomers failing the embedding process
    while len(self.already_seen) < 2**self.nCenters:
      bits = self.rand.getrandbits(self.nCenters)
      if bits in self.already_seen:
        continue

      self.already_seen.add(bits)
      yield bits


def GetStereoisomerCount(m, options=StereoEnumerationOptions()):
  """ returns an estimate (upper bound) of the number of possible stereoisomers for a molecule

   Arguments:
      - m: the molecule to work with
      - options: parameters controlling the enumeration


    >>> from rdkit import Chem
    >>> from rdkit.Chem.EnumerateStereoisomers import EnumerateStereoisomers, StereoEnumerationOptions
    >>> m = Chem.MolFromSmiles('BrC(Cl)(F)CCC(O)C')
    >>> GetStereoisomerCount(m)
    4
    >>> m = Chem.MolFromSmiles('CC(Cl)(O)C')
    >>> GetStereoisomerCount(m)
    1

    double bond stereochemistry is also included:

    >>> m = Chem.MolFromSmiles('BrC(Cl)(F)C=CC(O)C')
    >>> GetStereoisomerCount(m)
    8

    """
  tm = Chem.Mol(m)
  flippers = _getFlippers(tm, options)
  return 2**len(flippers)


def EnumerateStereoisomers(m, options=StereoEnumerationOptions(), verbose=False):
  """ returns a generator that yields possible stereoisomers for a molecule

    Arguments:
      - m: the molecule to work with
      - options: parameters controlling the enumeration
      - verbose: toggles how verbose the output is

    If m has stereogroups, they will be expanded

    A small example with 3 chiral atoms and 1 chiral bond (16 theoretical stereoisomers):

    >>> from rdkit import Chem
    >>> from rdkit.Chem.EnumerateStereoisomers import EnumerateStereoisomers, StereoEnumerationOptions
    >>> m = Chem.MolFromSmiles('BrC=CC1OC(C2)(F)C2(Cl)C1')
    >>> isomers = tuple(EnumerateStereoisomers(m))
    >>> len(isomers)
    16
    >>> for smi in sorted(Chem.MolToSmiles(x, isomericSmiles=True) for x in isomers):
    ...     print(smi)
    ...
    F[C@@]12C[C@@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@@]12C[C@@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@@]12C[C@@]1(Cl)C[C@H](/C=C/Br)O2
    F[C@@]12C[C@@]1(Cl)C[C@H](/C=C\\Br)O2
    F[C@@]12C[C@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@@]12C[C@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@@]12C[C@]1(Cl)C[C@H](/C=C/Br)O2
    F[C@@]12C[C@]1(Cl)C[C@H](/C=C\\Br)O2
    F[C@]12C[C@@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@]12C[C@@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@]12C[C@@]1(Cl)C[C@H](/C=C/Br)O2
    F[C@]12C[C@@]1(Cl)C[C@H](/C=C\\Br)O2
    F[C@]12C[C@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@]12C[C@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@]12C[C@]1(Cl)C[C@H](/C=C/Br)O2
    F[C@]12C[C@]1(Cl)C[C@H](/C=C\\Br)O2

    Because the molecule is constrained, not all of those isomers can
    actually exist. We can check that:

    >>> opts = StereoEnumerationOptions(tryEmbedding=True)
    >>> isomers = tuple(EnumerateStereoisomers(m, options=opts))
    >>> len(isomers)
    8
    >>> for smi in sorted(Chem.MolToSmiles(x,isomericSmiles=True) for x in isomers):
    ...     print(smi)
    ...
    F[C@@]12C[C@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@@]12C[C@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@@]12C[C@]1(Cl)C[C@H](/C=C/Br)O2
    F[C@@]12C[C@]1(Cl)C[C@H](/C=C\\Br)O2
    F[C@]12C[C@@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@]12C[C@@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@]12C[C@@]1(Cl)C[C@H](/C=C/Br)O2
    F[C@]12C[C@@]1(Cl)C[C@H](/C=C\\Br)O2

    Or we can force the output to only give us unique isomers:

    >>> m = Chem.MolFromSmiles('FC(Cl)C=CC=CC(F)Cl')
    >>> opts = StereoEnumerationOptions(unique=True)
    >>> isomers = tuple(EnumerateStereoisomers(m, options=opts))
    >>> len(isomers)
    10
    >>> for smi in sorted(Chem.MolToSmiles(x,isomericSmiles=True) for x in isomers):
    ...     print(smi)
    ...
    F[C@@H](Cl)/C=C/C=C/[C@@H](F)Cl
    F[C@@H](Cl)/C=C/C=C/[C@H](F)Cl
    F[C@@H](Cl)/C=C/C=C\\[C@H](F)Cl
    F[C@@H](Cl)/C=C\\C=C/[C@@H](F)Cl
    F[C@@H](Cl)/C=C\\C=C/[C@H](F)Cl
    F[C@@H](Cl)/C=C\\C=C\\[C@@H](F)Cl
    F[C@@H](Cl)/C=C\\C=C\\[C@H](F)Cl
    F[C@H](Cl)/C=C/C=C/[C@H](F)Cl
    F[C@H](Cl)/C=C\\C=C/[C@H](F)Cl
    F[C@H](Cl)/C=C\\C=C\\[C@H](F)Cl

    By default the code only expands unspecified stereocenters:

    >>> m = Chem.MolFromSmiles('BrC=C[C@H]1OC(C2)(F)C2(Cl)C1')
    >>> isomers = tuple(EnumerateStereoisomers(m))
    >>> len(isomers)
    8
    >>> for smi in sorted(Chem.MolToSmiles(x,isomericSmiles=True) for x in isomers):
    ...     print(smi)
    ...
    F[C@@]12C[C@@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@@]12C[C@@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@@]12C[C@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@@]12C[C@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@]12C[C@@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@]12C[C@@]1(Cl)C[C@@H](/C=C\\Br)O2
    F[C@]12C[C@]1(Cl)C[C@@H](/C=C/Br)O2
    F[C@]12C[C@]1(Cl)C[C@@H](/C=C\\Br)O2

    But we can change that behavior:

    >>> opts = StereoEnumerationOptions(onlyUnassigned=False)
    >>> isomers = tuple(EnumerateStereoisomers(m, options=opts))
    >>> len(isomers)
    16

    Since the result is a generator, we can allow exploring at least parts of very
    large result sets:

    >>> m = Chem.MolFromSmiles('Br' + '[CH](Cl)' * 20 + 'F')
    >>> opts = StereoEnumerationOptions(maxIsomers=0)
    >>> isomers = EnumerateStereoisomers(m, options=opts)
    >>> for x in range(5):
    ...   print(Chem.MolToSmiles(next(isomers),isomericSmiles=True))
    F[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)Br
    F[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)Br
    F[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)Br
    F[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@H](Cl)Br
    F[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)Br

    Or randomly sample a small subset. Note that if we want that sampling to be consistent
    across python versions we need to provide a random number seed:

    >>> m = Chem.MolFromSmiles('Br' + '[CH](Cl)' * 20 + 'F')
    >>> opts = StereoEnumerationOptions(maxIsomers=3,rand=0xf00d)
    >>> isomers = EnumerateStereoisomers(m, options=opts)
    >>> for smi in isomers: #sorted(Chem.MolToSmiles(x, isomericSmiles=True) for x in isomers):
    ...     print(Chem.MolToSmiles(smi))
    F[C@@H](Cl)[C@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)Br
    F[C@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)[C@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)Br
    F[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)[C@@H](Cl)[C@H](Cl)[C@@H](Cl)Br

    """
  tm = Chem.Mol(m)
  for atom in tm.GetAtoms():
    atom.ClearProp("_CIPCode")
  flippers = _getFlippers(tm, options)
  nCenters = len(flippers)
  if not nCenters:
    yield tm
    return

  if (options.maxIsomers == 0 or 2**nCenters <= options.maxIsomers):
    bitsource = _RangeBitsGenerator(nCenters)
  else:
    if options.rand is None:
      # deterministic random seed invariant to input atom order
      seed = hash(tuple(sorted([(a.GetDegree(), a.GetAtomicNum()) for a in tm.GetAtoms()])))
      rand = random.Random(seed)
    elif isinstance(options.rand, random.Random):
      # other implementations of Python random number generators
      # can inherit from this class to pick up utility methods
      rand = options.rand
    else:
      rand = random.Random(options.rand)

    bitsource = _UniqueRandomBitsGenerator(nCenters, options.maxIsomers, rand)

  isomersSeen = set()
  numIsomers = 0
  for bitflag in bitsource:
    for i in range(nCenters):
      flag = bool(bitflag & (1 << i))
      flippers[i].flip(flag)
    isomer = Chem.Mol(tm)
    Chem.SetDoubleBondNeighborDirections(isomer)
    isomer.ClearComputedProps()

    Chem.AssignStereochemistry(isomer, cleanIt=True, force=True, flagPossibleStereoCenters=True)
    if options.unique:
      cansmi = Chem.MolToSmiles(isomer, isomericSmiles=True)
      if cansmi in isomersSeen:
        continue

      isomersSeen.add(cansmi)

    if options.tryEmbedding:
      ntm = Chem.AddHs(isomer)
      cid = EmbedMolecule(ntm, randomSeed=bitflag)
      if cid >= 0:
        conf = Chem.Conformer(isomer.GetNumAtoms())
        for aid in range(isomer.GetNumAtoms()):
          conf.SetAtomPosition(aid, ntm.GetConformer().GetAtomPosition(aid))
        isomer.AddConformer(conf)
    else:
      cid = 1
    if cid >= 0:
      yield isomer
      numIsomers += 1
      if options.maxIsomers != 0 and numIsomers >= options.maxIsomers:
        break
    elif verbose:
      print("%s    failed to embed" % (Chem.MolToSmiles(isomer, isomericSmiles=True)))