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 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537
|
RDKit Cookbook
%%%%%%%%%%%%%%
What is this?
*************
This document provides examples of how to carry out particular tasks
using the RDKit functionality from Python. The contents have been
contributed by the RDKit community.
If you find mistakes, or have suggestions for improvements, please
either fix them yourselves in the source document (the .rst file) or
send them to the mailing list: rdkit-discuss@lists.sourceforge.net
(you will need to subscribe first)
Miscellaneous Topics
********************
Using a different aromaticity model
-----------------------------------
By default, the RDKit applies its own model of aromaticity (explained
in the RDKit Theory Book) when it reads in molecules. It is, however,
fairly easy to override this and use your own aromaticity model.
The easiest way to do this is it provide the molecules as SMILES with
the aromaticity set as you would prefer to have it. For example,
consider indole:
.. image:: images/indole1.png
By default the RDKit considers both rings to be aromatic:
>>> from rdkit import Chem
>>> m = Chem.MolFromSmiles('N1C=Cc2ccccc12')
>>> m.GetSubstructMatches(Chem.MolFromSmarts('c'))
((1,), (2,), (3,), (4,), (5,), (6,), (7,), (8,))
If you'd prefer to treat the five-membered ring as aliphatic, which is
how the input SMILES is written, you just need to do a partial
sanitization that skips the kekulization and aromaticity perception
steps:
>>> m2 = Chem.MolFromSmiles('N1C=Cc2ccccc12',sanitize=False)
>>> Chem.SanitizeMol(m2,sanitizeOps=Chem.SanitizeFlags.SANITIZE_ALL^Chem.SanitizeFlags.SANITIZE_KEKULIZE^Chem.SanitizeFlags.SANITIZE_SETAROMATICITY)
rdkit.Chem.rdmolops.SanitizeFlags.SANITIZE_NONE
>>> m2.GetSubstructMatches(Chem.MolFromSmarts('c'))
((3,), (4,), (5,), (6,), (7,), (8,))
It is, of course, also possible to write your own aromaticity
perception function, but that is beyond the scope of this document.
Manipulating Molecules
**********************
Cleaning up heterocycles
------------------------
Mailing list discussions:
* http://www.mail-archive.com/rdkit-discuss@lists.sourceforge.net/msg01185.html
* http://www.mail-archive.com/rdkit-discuss@lists.sourceforge.net/msg01162.html
* http://www.mail-archive.com/rdkit-discuss@lists.sourceforge.net/msg01900.html
* http://www.mail-archive.com/rdkit-discuss@lists.sourceforge.net/msg01901.html
The code:
.. testcode::
""" sanifix4.py
Contribution from James Davidson
"""
from rdkit import Chem
from rdkit.Chem import AllChem
def _FragIndicesToMol(oMol,indices):
em = Chem.EditableMol(Chem.Mol())
newIndices={}
for i,idx in enumerate(indices):
em.AddAtom(oMol.GetAtomWithIdx(idx))
newIndices[idx]=i
for i,idx in enumerate(indices):
at = oMol.GetAtomWithIdx(idx)
for bond in at.GetBonds():
if bond.GetBeginAtomIdx()==idx:
oidx = bond.GetEndAtomIdx()
else:
oidx = bond.GetBeginAtomIdx()
# make sure every bond only gets added once:
if oidx<idx:
continue
em.AddBond(newIndices[idx],newIndices[oidx],bond.GetBondType())
res = em.GetMol()
res.ClearComputedProps()
Chem.GetSymmSSSR(res)
res.UpdatePropertyCache(False)
res._idxMap=newIndices
return res
def _recursivelyModifyNs(mol,matches,indices=None):
if indices is None:
indices=[]
res=None
while len(matches) and res is None:
tIndices=indices[:]
nextIdx = matches.pop(0)
tIndices.append(nextIdx)
nm = Chem.Mol(mol)
nm.GetAtomWithIdx(nextIdx).SetNoImplicit(True)
nm.GetAtomWithIdx(nextIdx).SetNumExplicitHs(1)
cp = Chem.Mol(nm)
try:
Chem.SanitizeMol(cp)
except ValueError:
res,indices = _recursivelyModifyNs(nm,matches,indices=tIndices)
else:
indices=tIndices
res=cp
return res,indices
def AdjustAromaticNs(m,nitrogenPattern='[n&D2&H0;r5,r6]'):
"""
default nitrogen pattern matches Ns in 5 rings and 6 rings in order to be able
to fix: O=c1ccncc1
"""
Chem.GetSymmSSSR(m)
m.UpdatePropertyCache(False)
# break non-ring bonds linking rings:
em = Chem.EditableMol(m)
linkers = m.GetSubstructMatches(Chem.MolFromSmarts('[r]!@[r]'))
plsFix=set()
for a,b in linkers:
em.RemoveBond(a,b)
plsFix.add(a)
plsFix.add(b)
nm = em.GetMol()
for at in plsFix:
at=nm.GetAtomWithIdx(at)
if at.GetIsAromatic() and at.GetAtomicNum()==7:
at.SetNumExplicitHs(1)
at.SetNoImplicit(True)
# build molecules from the fragments:
fragLists = Chem.GetMolFrags(nm)
frags = [_FragIndicesToMol(nm,x) for x in fragLists]
# loop through the fragments in turn and try to aromatize them:
ok=True
for i,frag in enumerate(frags):
cp = Chem.Mol(frag)
try:
Chem.SanitizeMol(cp)
except ValueError:
matches = [x[0] for x in frag.GetSubstructMatches(Chem.MolFromSmarts(nitrogenPattern))]
lres,indices=_recursivelyModifyNs(frag,matches)
if not lres:
#print 'frag %d failed (%s)'%(i,str(fragLists[i]))
ok=False
break
else:
revMap={}
for k,v in frag._idxMap.iteritems():
revMap[v]=k
for idx in indices:
oatom = m.GetAtomWithIdx(revMap[idx])
oatom.SetNoImplicit(True)
oatom.SetNumExplicitHs(1)
if not ok:
return None
return m
Examples of using it:
.. testcode::
smis= ('O=c1ccc2ccccc2n1',
'Cc1nnnn1C',
'CCc1ccc2nc(=O)c(cc2c1)Cc1nnnn1C1CCCCC1',
'c1cnc2cc3ccnc3cc12',
'c1cc2cc3ccnc3cc2n1',
'O=c1ccnc(c1)-c1cnc2cc3ccnc3cc12',
'O=c1ccnc(c1)-c1cc1',
)
for smi in smis:
m = Chem.MolFromSmiles(smi,False)
try:
m.UpdatePropertyCache(False)
cp = Chem.Mol(m)
Chem.SanitizeMol(cp)
m = cp
print 'fine:',Chem.MolToSmiles(m)
except ValueError:
nm=AdjustAromaticNs(m)
if nm is not None:
Chem.SanitizeMol(nm)
print 'fixed:',Chem.MolToSmiles(nm)
else:
print 'still broken:',smi
This produces:
.. testoutput::
fixed: O=c1ccc2ccccc2[nH]1
fine: Cc1nnnn1C
fixed: CCc1ccc2[nH]c(=O)c(Cc3nnnn3C3CCCCC3)cc2c1
fine: C1=Cc2cc3c(cc2=N1)C=CN=3
fine: C1=Cc2cc3c(cc2=N1)N=CC=3
fixed: O=c1cc[nH]c(C2=CN=c3cc4c(cc32)=NC=C4)c1
still broken: O=c1ccnc(c1)-c1cc1
Parallel conformation generation
--------------------------------
Mailing list discussion:
http://www.mail-archive.com/rdkit-discuss@lists.sourceforge.net/msg02648.html
The code::
""" contribution from Andrew Dalke """
import sys
from rdkit import Chem
from rdkit.Chem import AllChem
# Download this from http://pypi.python.org/pypi/futures
from concurrent import futures
# Download this from http://pypi.python.org/pypi/progressbar
import progressbar
## On my machine, it takes 39 seconds with 1 worker and 10 seconds with 4.
## 29.055u 0.102s 0:28.68 101.6% 0+0k 0+3io 0pf+0w
#max_workers=1
## With 4 threads it takes 11 seconds.
## 34.933u 0.188s 0:10.89 322.4% 0+0k 125+1io 0pf+0w
max_workers=4
# (The "u"ser time includes time spend in the children processes.
# The wall-clock time is 28.68 and 10.89 seconds, respectively.)
# This function is called in the subprocess.
# The parameters (molecule and number of conformers) are passed via a Python
def generateconformations(m, n):
m = Chem.AddHs(m)
ids=AllChem.EmbedMultipleConfs(m, numConfs=n)
for id in ids:
AllChem.UFFOptimizeMolecule(m, confId=id)
# EmbedMultipleConfs returns a Boost-wrapped type which
# cannot be pickled. Convert it to a Python list, which can.
return m, list(ids)
smi_input_file, sdf_output_file = sys.argv[1:3]
n = int(sys.argv[3])
writer = Chem.SDWriter(sdf_output_file)
suppl = Chem.SmilesMolSupplier(smi_input_file, titleLine=False)
with futures.ProcessPoolExecutor(max_workers=max_workers) as executor:
# Submit a set of asynchronous jobs
jobs = []
for mol in suppl:
if mol:
job = executor.submit(generateconformations, mol, n)
jobs.append(job)
widgets = ["Generating conformations; ", progressbar.Percentage(), " ",
progressbar.ETA(), " ", progressbar.Bar()]
pbar = progressbar.ProgressBar(widgets=widgets, maxval=len(jobs))
for job in pbar(futures.as_completed(jobs)):
mol,ids=job.result()
for id in ids:
writer.write(mol, confId=id)
writer.close()
Neutralizing Charged Molecules
------------------------------
Mailing list discussion:
http://www.mail-archive.com/rdkit-discuss@lists.sourceforge.net/msg02648.html
Wiki page:
http://code.google.com/p/rdkit/wiki/NeutralisingCompounds
The code:
.. testcode::
""" contribution from Hans de Winter """
from rdkit import Chem
from rdkit.Chem import AllChem
def _InitialiseNeutralisationReactions():
patts= (
# Imidazoles
('[n+;H]','n'),
# Amines
('[N+;!H0]','N'),
# Carboxylic acids and alcohols
('[$([O-]);!$([O-][#7])]','O'),
# Thiols
('[S-;X1]','S'),
# Sulfonamides
('[$([N-;X2]S(=O)=O)]','N'),
# Enamines
('[$([N-;X2][C,N]=C)]','N'),
# Tetrazoles
('[n-]','[nH]'),
# Sulfoxides
('[$([S-]=O)]','S'),
# Amides
('[$([N-]C=O)]','N'),
)
return [(Chem.MolFromSmarts(x),Chem.MolFromSmiles(y,False)) for x,y in patts]
_reactions=None
def NeutraliseCharges(smiles, reactions=None):
global _reactions
if reactions is None:
if _reactions is None:
_reactions=_InitialiseNeutralisationReactions()
reactions=_reactions
mol = Chem.MolFromSmiles(smiles)
replaced = False
for i,(reactant, product) in enumerate(reactions):
while mol.HasSubstructMatch(reactant):
replaced = True
rms = AllChem.ReplaceSubstructs(mol, reactant, product)
mol = rms[0]
if replaced:
return (Chem.MolToSmiles(mol,True), True)
else:
return (smiles, False)
Examples of using it:
.. testcode::
smis=("c1cccc[nH+]1",
"C[N+](C)(C)C","c1ccccc1[NH3+]",
"CC(=O)[O-]","c1ccccc1[O-]",
"CCS",
"C[N-]S(=O)(=O)C",
"C[N-]C=C","C[N-]N=C",
"c1ccc[n-]1",
"CC[N-]C(=O)CC")
for smi in smis:
(molSmiles, neutralised) = NeutraliseCharges(smi)
print smi,"->",molSmiles
This produces:
.. testoutput::
c1cccc[nH+]1 -> c1ccncc1
C[N+](C)(C)C -> C[N+](C)(C)C
c1ccccc1[NH3+] -> Nc1ccccc1
CC(=O)[O-] -> CC(=O)O
c1ccccc1[O-] -> Oc1ccccc1
CCS -> CCS
C[N-]S(=O)(=O)C -> CNS(C)(=O)=O
C[N-]C=C -> C=CNC
C[N-]N=C -> C=NNC
c1ccc[n-]1 -> c1cc[nH]c1
CC[N-]C(=O)CC -> CCNC(=O)CC
Using scikit-learn with RDKit
-----------------------------
scikit-learn is a machine-learning library for Python
containing a variety of supervised and unsupervised methods.
The documention can be found here:
http://scikit-learn.org/stable/user_guide.html
RDKit fingerprints can be used to train machine-learning
models from scikit-learn. Here is an example for random forest:
The code::
from rdkit import Chem, DataStructs
from rdkit.Chem import AllChem
from sklearn.ensemble import RandomForestClassifier
import numpy
# generate four molecules
m1 = Chem.MolFromSmiles('c1ccccc1')
m2 = Chem.MolFromSmiles('c1ccccc1CC')
m3 = Chem.MolFromSmiles('c1ccncc1')
m4 = Chem.MolFromSmiles('c1ccncc1CC')
mols = [m1, m2, m3, m4]
# generate fingeprints: Morgan fingerprint with radius 2
fps = [AllChem.GetMorganFingerprintAsBitVect(m, 2) for m in mols]
# convert the RDKit explicit vectors into numpy arrays
np_fps = []
for fp in fps:
arr = numpy.zeros((1,))
DataStructs.ConvertToNumpyArray(fp, arr)
np_fps.append(arr)
# get a random forest classifiert with 100 trees
rf = RandomForestClassifier(n_estimators=100, random_state=1123)
# train the random forest
# with the first two molecules being actives (class 1) and
# the last two being inactives (class 0)
ys_fit = [1, 1, 0, 0]
rf.fit(np_fps, ys_fit)
# use the random forest to predict a new molecule
m5 = Chem.MolFromSmiles('c1ccccc1O')
fp = numpy.zeros((1,))
DataStructs.ConvertToNumpyArray(AllChem.GetMorganFingerprintAsBitVect(m5, 2), fp)
print rf.predict(fp)
print rf.predict_proba(fp)
The output with scikit-learn version 0.13 is:
[1]
[[ 0.14 0.86]]
Generating a similarity map for this model.
The code::
from rdkit.Chem.Draw import SimilarityMaps
# helper function
def getProba(fp, predictionFunction):
return predictionFunction(fp)[0][1]
m5 = Chem.MolFromSmiles('c1ccccc1O')
fig, maxweight = SimilarityMaps.GetSimilarityMapForModel(m5, SimilarityMaps.GetMorganFingerprint, lambda x: getProba(x, rf.predict_proba))
This produces:
.. image:: images/similarity_map_rf.png
Using custom MCS atom types
---------------------------
Mailing list discussion:
http://www.mail-archive.com/rdkit-discuss@lists.sourceforge.net/msg03676.html
IPython notebook:
http://nbviewer.ipython.org/gist/greglandrum/8351725
https://gist.github.com/greglandrum/8351725
The goal is to be able to use custom atom types in the MCS code, yet
still be able to get a readable SMILES for the MCS. We will use the
MCS code's option to use isotope information in the matching and then
set bogus isotope values that contain our isotope information.
The code:
.. testcode::
# our test molecules:
smis=["COc1ccc(C(Nc2nc3c(ncn3COCC=O)c(=O)[nH]2)(c2ccccc2)c2ccccc2)cc1",
"COc1ccc(C(Nc2nc3c(ncn3COC(CO)(CO)CO)c(=O)[nH]2)(c2ccccc2)c2ccccc2)cc1"]
ms = [Chem.MolFromSmiles(x) for x in smis]
from rdkit import Chem
from rdkit.Chem import MCS
def label(a):
" a simple hash combining atom number and hybridization "
return 100*int(a.GetHybridization())+a.GetAtomicNum()
# copy the molecules, since we will be changing them
nms = [Chem.Mol(x) for x in ms]
for nm in nms:
for at in nm.GetAtoms():
at.SetIsotope(label(at))
mcs=MCS.FindMCS(nms,atomCompare='isotopes')
print mcs.smarts
This generates the following output:
.. testoutput::
[406*]-[308*]-[306*]:1:[306*]:[306*]:[306*](:[306*]:[306*]:1)-[406*](-[306*]:1:[306*]:[306*]:[306*]:[306*]:[306*]:1)(-[306*]:1:[306*]:[306*]:[306*]:[306*]:[306*]:1)-[307*]-[306*]:1:[307*]:[306*]:2:[306*](:[306*](:[307*]:1)=[308*]):[307*]:[306*]:[307*]:2-[406*]-[408*]-[406*]
That's what we asked for, but it's not exactly readable. We can get to a more readable form in a two step process:
1) Do a substructure match of the MCS onto a copied molecule
2) Generate SMILES for the original molecule, using only the atoms that matched in the copy.
This works because we know that the atom indices in the copies and the original molecules are the same.
.. testcode::
def getMCSSmiles(mol,labelledMol,mcs):
mcsp = Chem.MolFromSmarts(mcs.smarts)
match = labelledMol.GetSubstructMatch(mcsp)
return Chem.MolFragmentToSmiles(ms[0],atomsToUse=match,
isomericSmiles=True,
canonical=False)
print getMCSSmiles(ms[0],nms[0],mcs)
.. testoutput::
COc1ccc(C(Nc2nc3c(ncn3COC)c(=O)[nH]2)(c2ccccc2)c2ccccc2)cc1
That's what we were looking for.
License
*******
This document is copyright (C) 2012-2014 by Greg Landrum
This work is licensed under the Creative Commons Attribution-ShareAlike 3.0 License.
To view a copy of this license, visit http://creativecommons.org/licenses/by-sa/3.0/ or send a letter to Creative Commons, 543 Howard Street, 5th Floor, San Francisco, California, 94105, USA.
The intent of this license is similar to that of the RDKit itself.
In simple words: “Do whatever you want with it, but please give us some credit.”
|