File: __init__.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 (179 lines) | stat: -rwxr-xr-x 5,579 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
#
#  Copyright (C) 2000-2017  greg Landrum and Rational Discovery LLC
#
#   @@ All Rights Reserved @@
#  This file is part of the RDKit.
#  The contents are covered by the terms of the BSD license
#  which is included in the file license.txt, found at the root
#  of the RDKit source tree.
#
""" A module for molecules and stuff

 see Chem/index.html in the doc tree for documentation

"""
from rdkit import rdBase
from rdkit import RDConfig

from rdkit import DataStructs
from rdkit.Geometry import rdGeometry
from rdkit.Chem import rdchem
_HasSubstructMatchStr = rdchem._HasSubstructMatchStr
from rdkit.Chem.rdchem import *
from rdkit.Chem.rdmolfiles import *
from rdkit.Chem.rdmolops import *
from rdkit.Chem.rdCIPLabeler import *
from rdkit.Chem.inchi import *
try:
  # This is an optional component of the build
  from rdkit.Chem.rdMolInterchange import *
except ImportError:
  pass

# Coordgen needs to know where its template file is.
# The default install puts it in RDDataDir
try:
  from rdkit.Chem import rdCoordGen
except ImportError:
  pass
else:
  templDir = RDConfig.RDDataDir
  if templDir[-1] != '/':
    templDir += '/'
  rdCoordGen.SetDefaultTemplateFileDir(templDir)


def QuickSmartsMatch(smi, sma, unique=True, display=False):
  m = MolFromSmiles(smi)
  p = MolFromSmarts(sma)
  res = m.GetSubstructMatches(p, unique)
  if display:
    pass
  return res


def CanonSmiles(smi, useChiral=1):
  m = MolFromSmiles(smi)
  return MolToSmiles(m, useChiral)


def SupplierFromFilename(fileN, delim='', **kwargs):
  ext = fileN.split('.')[-1].lower()
  if ext == 'sdf':
    suppl = SDMolSupplier(fileN, **kwargs)
  elif ext == 'csv':
    if not delim:
      delim = ','
    suppl = SmilesMolSupplier(fileN, delimiter=delim, **kwargs)
  elif ext == 'txt':
    if not delim:
      delim = '\t'
    suppl = SmilesMolSupplier(fileN, delimiter=delim, **kwargs)
  elif ext == 'tdt':
    suppl = TDTMolSupplier(fileN, delimiter=delim, **kwargs)
  else:
    raise ValueError("unrecognized extension: %s" % ext)

  return suppl


def FindMolChiralCenters(mol, force=True, includeUnassigned=False, includeCIP=True,
                         useLegacyImplementation=True):
  """
    >>> from rdkit import Chem
    >>> mol = Chem.MolFromSmiles('[C@H](Cl)(F)Br')
    >>> FindMolChiralCenters(mol)
    [(0, 'R')]
    >>> mol = Chem.MolFromSmiles('[C@@H](Cl)(F)Br')
    >>> FindMolChiralCenters(mol)
    [(0, 'S')]

    >>> FindMolChiralCenters(Chem.MolFromSmiles('CCC'))
    []

    By default unassigned stereo centers are not reported:

    >>> mol = Chem.MolFromSmiles('C[C@H](F)C(F)(Cl)Br')
    >>> FindMolChiralCenters(mol,force=True)
    [(1, 'S')]

    but this can be changed:

    >>> FindMolChiralCenters(mol,force=True,includeUnassigned=True)
    [(1, 'S'), (3, '?')]

    The handling of unassigned stereocenters for dependent stereochemistry is not correct 
    using the legacy implementation:

    >>> Chem.FindMolChiralCenters(Chem.MolFromSmiles('C1CC(C)C(C)C(C)C1'),includeUnassigned=True)
    [(2, '?'), (6, '?')]
    >>> Chem.FindMolChiralCenters(Chem.MolFromSmiles('C1C[C@H](C)C(C)[C@H](C)C1'),includeUnassigned=True)
    [(2, 'S'), (4, '?'), (6, 'R')]

    But works with the new implementation:

    >>> Chem.FindMolChiralCenters(Chem.MolFromSmiles('C1CC(C)C(C)C(C)C1'),includeUnassigned=True, useLegacyImplementation=False)
    [(2, '?'), (4, '?'), (6, '?')]

    Note that the new implementation also gets the correct descriptors for para-stereochemistry:

    >>> Chem.FindMolChiralCenters(Chem.MolFromSmiles('C1C[C@H](C)[C@H](C)[C@H](C)C1'),useLegacyImplementation=False)
    [(2, 'S'), (4, 's'), (6, 'R')]

    With the new implementation, if you don't care about the CIP labels of stereocenters, you can save
    some time by disabling those:

    >>> Chem.FindMolChiralCenters(Chem.MolFromSmiles('C1C[C@H](C)[C@H](C)[C@H](C)C1'), includeCIP=False, useLegacyImplementation=False)
    [(2, 'Tet_CCW'), (4, 'Tet_CCW'), (6, 'Tet_CCW')]

  """
  if useLegacyImplementation:
    AssignStereochemistry(mol, force=force, flagPossibleStereoCenters=includeUnassigned)
    centers = []
    for atom in mol.GetAtoms():
      if atom.HasProp('_CIPCode'):
        centers.append((atom.GetIdx(), atom.GetProp('_CIPCode')))
      elif includeUnassigned and atom.HasProp('_ChiralityPossible'):
        centers.append((atom.GetIdx(), '?'))
  else:
    centers = []
    itms = FindPotentialStereo(mol)
    if includeCIP:
      atomsToLabel = []
      bondsToLabel = []
      for si in itms:
        if si.type == StereoType.Atom_Tetrahedral:
          atomsToLabel.append(si.centeredOn)
        elif si.type == StereoType.Bond_Double:
          bondsToLabel.append(si.centeredOn)
      AssignCIPLabels(mol, atomsToLabel=atomsToLabel, bondsToLabel=bondsToLabel)
    for si in itms:
      if si.type == StereoType.Atom_Tetrahedral and (includeUnassigned
                                                     or si.specified == StereoSpecified.Specified):
        idx = si.centeredOn
        atm = mol.GetAtomWithIdx(idx)
        if includeCIP and atm.HasProp("_CIPCode"):
          code = atm.GetProp("_CIPCode")
        else:
          if si.specified:
            code = str(si.descriptor)
          else:
            code = '?'
            atm.SetIntProp('_ChiralityPossible',1)
        centers.append((idx, code))
  return centers


#------------------------------------
#
#  doctest boilerplate
#
def _test():
  import doctest, sys
  return doctest.testmod(sys.modules["__main__"])


if __name__ == '__main__':
  import sys
  failed, tried = _test()
  sys.exit(failed)