| Trees | Indices | Help |
|
|---|
|
|
1 # $Id: Recap.py 1742 2011-06-10 11:13:41Z glandrum $
2 #
3 # Copyright (c) 2007, Novartis Institutes for BioMedical Research Inc.
4 # All rights reserved.
5 #
6 # Redistribution and use in source and binary forms, with or without
7 # modification, are permitted provided that the following conditions are
8 # met:
9 #
10 # * Redistributions of source code must retain the above copyright
11 # notice, this list of conditions and the following disclaimer.
12 # * Redistributions in binary form must reproduce the above
13 # copyright notice, this list of conditions and the following
14 # disclaimer in the documentation and/or other materials provided
15 # with the distribution.
16 # * Neither the name of Novartis Institutes for BioMedical Research Inc.
17 # nor the names of its contributors may be used to endorse or promote
18 # products derived from this software without specific prior written permission.
19 #
20 # THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 # "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 # LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR
23 # A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT
24 # OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
25 # SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT
26 # LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE,
27 # DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY
28 # THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
29 # (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE
30 # OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
31 #
32 """ Implementation of the RECAP algorithm from Lewell et al. JCICS *38* 511-522 (1998)
33
34 The published algorithm is implemented more or less without
35 modification. The results are returned as a hierarchy of nodes instead
36 of just as a set of fragments. The hope is that this will allow a bit
37 more flexibility in working with the results.
38
39 For example:
40 >>> m = Chem.MolFromSmiles('C1CC1Oc1ccccc1-c1ncc(OC)cc1')
41 >>> res = Recap.RecapDecompose(m)
42 >>> res
43 <Chem.Recap.RecapHierarchyNode object at 0x00CDB5D0>
44 >>> res.children.keys()
45 ['[*]C1CC1', '[*]c1ccccc1-c1ncc(OC)cc1', '[*]c1ccc(OC)cn1', '[*]c1ccccc1OC1CC1']
46 >>> res.GetAllChildren().keys()
47 ['[*]c1ccccc1[*]', '[*]C1CC1', '[*]c1ccccc1-c1ncc(OC)cc1', '[*]c1ccc(OC)cn1', '[*]c1ccccc1OC1CC1']
48
49
50 To get the standard set of RECAP results, use GetLeaves():
51 >>> leaves=res.GetLeaves()
52 >>> leaves.keys()
53 ['[*]c1ccccc1[*]', '[*]c1ccc(OC)cn1', '[*]C1CC1']
54 >>> leaf = leaves['[*]C1CC1']
55 >>> leaf.mol
56 <Chem.rdchem.Mol object at 0x00CBE0F0>
57
58
59 """
60 import weakref
61 from rdkit import Chem
62 from rdkit.Chem import rdChemReactions as Reactions
63 import sys
64
65 # These are the definitions that will be applied to fragment molecules:
66 reactionDefs = (
67 "[#7;+0;D2,D3:1]!@C(!@=O)!@[#7;+0;D2,D3:2]>>[*][#7:1].[#7:2][*]", # urea
68
69 "[C;!$(C([#7])[#7]):1](=!@[O:2])!@[#7;+0;!D1:3]>>[*][C:1]=[O:2].[*][#7:3]", # amide
70
71 "[C:1](=!@[O:2])!@[O;+0:3]>>[*][C:1]=[O:2].[O:3][*]", # ester
72
73 "[N;!D1;+0;!$(N-C=[#7,#8,#15,#16])](-!@[*:1])-!@[*:2]>>[*][*:1].[*:2][*]", # amines
74 #"[N;!D1](!@[*:1])!@[*:2]>>[*][*:1].[*:2][*]", # amines
75
76 # again: what about aromatics?
77 "[#7;R;D3;+0:1]-!@[*:2]>>[*][#7:1].[*:2][*]", # cyclic amines
78
79 "[#6:1]-!@[O;+0]-!@[#6:2]>>[#6:1][*].[*][#6:2]", # ether
80
81 "[C:1]=!@[C:2]>>[C:1][*].[*][C:2]", # olefin
82
83 "[n;+0:1]-!@[C:2]>>[n:1][*].[C:2][*]", # aromatic nitrogen - aliphatic carbon
84
85 "[O:3]=[C:4]-@[N;+0:1]-!@[C:2]>>[O:3]=[C:4]-[N:1][*].[C:2][*]", # lactam nitrogen - aliphatic carbon
86
87 "[c:1]-!@[c:2]>>[c:1][*].[*][c:2]", # aromatic carbon - aromatic carbon
88
89 "[n;+0:1]-!@[c:2]>>[n:1][*].[*][c:2]", # aromatic nitrogen - aromatic carbon *NOTE* this is not part of the standard recap set.
90
91 "[#7;+0;D2,D3:1]-!@[S:2](=[O:3])=[O:4]>>[#7:1][*].[*][S:2](=[O:3])=[O:4]", # sulphonamide
92 )
93
94 reactions = tuple([Reactions.ReactionFromSmarts(x) for x in reactionDefs])
95
96
98 """ This class is used to hold the Recap hiearchy
99 """
100 mol=None
101 children=None
102 parents=None
103 smiles = None
108
110 " returns a dictionary, keyed by SMILES, of children "
111 res = {}
112 for smi,child in self.children.iteritems():
113 res[smi] = child
114 child._gacRecurse(res,terminalOnly=False)
115 return res
116
118 " returns a dictionary, keyed by SMILES, of leaf (terminal) nodes "
119 res = {}
120 for smi,child in self.children.iteritems():
121 if not len(child.children):
122 res[smi] = child
123 else:
124 child._gacRecurse(res,terminalOnly=True)
125 return res
126
128 """ returns all the nodes in the hierarchy tree that contain this
129 node as a child
130 """
131 if not self.parents:
132 res = [self]
133 else:
134 res = []
135 for p in self.parents.values():
136 for uP in p.getUltimateParents():
137 if uP not in res:
138 res.append(uP)
139 return res
140
142 for smi,child in self.children.iteritems():
143 if not terminalOnly or not len(child.children):
144 res[smi] = child
145 child._gacRecurse(res,terminalOnly=terminalOnly)
146
151
152
154 """ returns the recap decomposition for a molecule """
155 mSmi = Chem.MolToSmiles(mol,1)
156
157 if allNodes is None:
158 allNodes={}
159 if allNodes.has_key(mSmi):
160 return allNodes[mSmi]
161
162 res = RecapHierarchyNode(mol)
163 res.smiles =mSmi
164 activePool={mSmi:res}
165 allNodes[mSmi]=res
166 while activePool:
167 nSmi = activePool.keys()[0]
168 node = activePool.pop(nSmi)
169 if not node.mol: continue
170 for rxnIdx,reaction in enumerate(reactions):
171 if onlyUseReactions and rxnIdx not in onlyUseReactions:
172 continue
173 #print ' .',nSmi
174 #print ' !!!!',rxnIdx,nSmi,reactionDefs[rxnIdx]
175 ps = reaction.RunReactants((node.mol,))
176 #print ' ',len(ps)
177 if ps:
178 for prodSeq in ps:
179 seqOk=True
180 # we want to disqualify small fragments, so sort the product sequence by size
181 # and then look for "forbidden" fragments
182 prodSeq = [(prod.GetNumAtoms(onlyHeavy=True),prod) for prod in prodSeq]
183 prodSeq.sort()
184 for nats,prod in prodSeq:
185 try:
186 Chem.SanitizeMol(prod)
187 except:
188 continue
189 pSmi = Chem.MolToSmiles(prod,1)
190 if minFragmentSize>0:
191 nDummies = pSmi.count('*')
192 if nats-nDummies<minFragmentSize:
193 seqOk=False
194 break
195 # don't forget after replacing dummy atoms to remove any empty
196 # branches:
197 elif pSmi.replace('[*]','').replace('()','') in ('','C','CC','CCC'):
198 seqOk=False
199 break
200 prod.pSmi = pSmi
201 if seqOk:
202 for nats,prod in prodSeq:
203 pSmi = prod.pSmi
204 #print '\t',nats,pSmi
205 if not allNodes.has_key(pSmi):
206 pNode = RecapHierarchyNode(prod)
207 pNode.smiles=pSmi
208 pNode.parents[nSmi]=weakref.proxy(node)
209 node.children[pSmi]=pNode
210 activePool[pSmi] = pNode
211 allNodes[pSmi]=pNode
212 else:
213 pNode=allNodes[pSmi]
214 pNode.parents[nSmi]=weakref.proxy(node)
215 node.children[pSmi]=pNode
216 #print ' >>an:',allNodes.keys()
217 return res
218
219
220 # ------- ------- ------- ------- ------- ------- ------- -------
221 # Begin testing code
222 if __name__=='__main__':
223 import unittest
226 m = Chem.MolFromSmiles('C1CC1Oc1ccccc1-c1ncc(OC)cc1')
227 res = RecapDecompose(m)
228 self.failUnless(res)
229 self.failUnless(len(res.children.keys())==4)
230 self.failUnless(len(res.GetAllChildren().keys())==5)
231 self.failUnless(len(res.GetLeaves().keys())==3)
233 m = Chem.MolFromSmiles('CCCOCCC')
234 res = RecapDecompose(m)
235 self.failUnless(res)
236 self.failUnless(res.children=={})
238 allNodes={}
239 m = Chem.MolFromSmiles('c1ccccc1-c1ncccc1')
240 res = RecapDecompose(m,allNodes=allNodes)
241 self.failUnless(res)
242 self.failUnless(len(res.children.keys())==2)
243 self.failUnless(len(allNodes.keys())==3)
244
245 m = Chem.MolFromSmiles('COc1ccccc1-c1ncccc1')
246 res = RecapDecompose(m,allNodes=allNodes)
247 self.failUnless(res)
248 self.failUnless(len(res.children.keys())==2)
249 # we get two more nodes from that:
250 self.failUnless(len(allNodes.keys())==5)
251 self.failUnless(allNodes.has_key('[*]c1ccccc1OC'))
252 self.failUnless(allNodes.has_key('[*]c1ccccc1'))
253
254 m = Chem.MolFromSmiles('C1CC1Oc1ccccc1-c1ncccc1')
255 res = RecapDecompose(m,allNodes=allNodes)
256 self.failUnless(res)
257 self.failUnless(len(res.children.keys())==4)
258 self.failUnless(len(allNodes.keys())==10)
259
261 m = Chem.MolFromSmiles('c1ccccc1OC(Oc1ccccc1)Oc1ccccc1')
262 res = RecapDecompose(m)
263 self.failUnless(res)
264 self.failUnless(len(res.GetLeaves())==2)
265 ks = res.GetLeaves().keys()
266 self.failIf('[*]C([*])[*]' in ks)
267 self.failUnless('[*]c1ccccc1' in ks)
268 self.failUnless('[*]C([*])Oc1ccccc1' in ks)
269
271 m = Chem.MolFromSmiles('C1CCCCN1CCCC')
272 res = RecapDecompose(m)
273 self.failUnless(res)
274 self.failUnless(len(res.GetLeaves())==2)
275 ks = res.GetLeaves().keys()
276 self.failUnless('[*]N1CCCCC1' in ks)
277 self.failUnless('[*]CCCC' in ks)
278
280 m = Chem.MolFromSmiles('CCCOCCC')
281 res = RecapDecompose(m)
282 self.failUnless(res)
283 self.failUnless(res.children=={})
284 res = RecapDecompose(m,minFragmentSize=3)
285 self.failUnless(res)
286 self.failUnless(len(res.GetLeaves())==1)
287 ks = res.GetLeaves().keys()
288 self.failUnless('[*]CCC' in ks)
289
290 m = Chem.MolFromSmiles('CCCOCC')
291 res = RecapDecompose(m,minFragmentSize=3)
292 self.failUnless(res)
293 self.failUnless(res.children=={})
294
295 m = Chem.MolFromSmiles('CCCOCCOC')
296 res = RecapDecompose(m,minFragmentSize=2)
297 self.failUnless(res)
298 self.failUnless(len(res.GetLeaves())==2)
299 ks = res.GetLeaves().keys()
300 self.failUnless('[*]CCC' in ks)
301 ks = res.GetLeaves().keys()
302 self.failUnless('[*]CCOC' in ks)
303
305 m = Chem.MolFromSmiles('C1CC1C(=O)NC1OC1')
306 res = RecapDecompose(m,onlyUseReactions=[1])
307 self.failUnless(res)
308 self.failUnless(len(res.GetLeaves())==2)
309 ks = res.GetLeaves().keys()
310 self.failUnless('[*]C(=O)C1CC1' in ks)
311 self.failUnless('[*]NC1CO1' in ks)
312
313 m = Chem.MolFromSmiles('C1CC1C(=O)N(C)C1OC1')
314 res = RecapDecompose(m,onlyUseReactions=[1])
315 self.failUnless(res)
316 self.failUnless(len(res.GetLeaves())==2)
317 ks = res.GetLeaves().keys()
318 self.failUnless('[*]C(=O)C1CC1' in ks)
319 self.failUnless('[*]N(C)C1CO1' in ks)
320
321 m = Chem.MolFromSmiles('C1CC1C(=O)n1cccc1')
322 res = RecapDecompose(m,onlyUseReactions=[1])
323 self.failUnless(res)
324 self.failUnless(len(res.GetLeaves())==2)
325 ks = res.GetLeaves().keys()
326 self.failUnless('[*]C(=O)C1CC1' in ks)
327 self.failUnless('[*]n1cccc1' in ks)
328
329 m = Chem.MolFromSmiles('C1CC1C(=O)CC1OC1')
330 res = RecapDecompose(m,onlyUseReactions=[1])
331 self.failUnless(res)
332 self.failUnless(len(res.GetLeaves())==0)
333
334 m = Chem.MolFromSmiles('C1CCC(=O)NC1')
335 res = RecapDecompose(m,onlyUseReactions=[1])
336 self.failUnless(res)
337 self.failUnless(len(res.GetLeaves())==0)
338
339 m = Chem.MolFromSmiles('CC(=O)NC')
340 res = RecapDecompose(m,onlyUseReactions=[1])
341 self.failUnless(res)
342 self.failUnless(len(res.GetLeaves())==2)
343 ks = res.GetLeaves().keys()
344
345 m = Chem.MolFromSmiles('CC(=O)N')
346 res = RecapDecompose(m,onlyUseReactions=[1])
347 self.failUnless(res)
348 self.failUnless(len(res.GetLeaves())==0)
349
350 m = Chem.MolFromSmiles('C(=O)NCCNC(=O)CC')
351 res = RecapDecompose(m,onlyUseReactions=[1])
352 self.failUnless(res)
353 self.failUnless(len(res.children)==4)
354 self.failUnless(len(res.GetLeaves())==3)
355
357 m = Chem.MolFromSmiles('C1CC1C(=O)OC1OC1')
358 res = RecapDecompose(m,onlyUseReactions=[2])
359 self.failUnless(res)
360 self.failUnless(len(res.GetLeaves())==2)
361 ks = res.GetLeaves().keys()
362 self.failUnless('[*]C(=O)C1CC1' in ks)
363 self.failUnless('[*]OC1CO1' in ks)
364
365 m = Chem.MolFromSmiles('C1CC1C(=O)CC1OC1')
366 res = RecapDecompose(m,onlyUseReactions=[2])
367 self.failUnless(res)
368 self.failUnless(len(res.GetLeaves())==0)
369
370 m = Chem.MolFromSmiles('C1CCC(=O)OC1')
371 res = RecapDecompose(m,onlyUseReactions=[2])
372 self.failUnless(res)
373 self.failUnless(len(res.GetLeaves())==0)
374
376 m = Chem.MolFromSmiles('C1CC1NC(=O)NC1OC1')
377 res = RecapDecompose(m,onlyUseReactions=[0])
378 self.failUnless(res)
379 self.failUnless(len(res.GetLeaves())==2)
380 ks = res.GetLeaves().keys()
381 self.failUnless('[*]NC1CC1' in ks)
382 self.failUnless('[*]NC1CO1' in ks)
383
384 m = Chem.MolFromSmiles('C1CC1NC(=O)N(C)C1OC1')
385 res = RecapDecompose(m,onlyUseReactions=[0])
386 self.failUnless(res)
387 self.failUnless(len(res.GetLeaves())==2)
388 ks = res.GetLeaves().keys()
389 self.failUnless('[*]NC1CC1' in ks)
390 self.failUnless('[*]N(C)C1CO1' in ks)
391
392 m = Chem.MolFromSmiles('C1CCNC(=O)NC1C')
393 res = RecapDecompose(m,onlyUseReactions=[0])
394 self.failUnless(res)
395 self.failUnless(len(res.GetLeaves())==0)
396
397 m = Chem.MolFromSmiles('c1cccn1C(=O)NC1OC1')
398 res = RecapDecompose(m,onlyUseReactions=[0])
399 self.failUnless(res)
400 self.failUnless(len(res.GetLeaves())==2)
401 ks = res.GetLeaves().keys()
402 self.failUnless('[*]n1cccc1' in ks)
403 self.failUnless('[*]NC1CO1' in ks)
404
405 m = Chem.MolFromSmiles('c1cccn1C(=O)n1c(C)ccc1')
406 res = RecapDecompose(m,onlyUseReactions=[0])
407 self.failUnless(res)
408 self.failUnless(len(res.GetLeaves())==2)
409 ks = res.GetLeaves().keys()
410 self.failUnless('[*]n1c(C)ccc1' in ks)
411
413 m = Chem.MolFromSmiles('C1CC1N(C1NC1)C1OC1')
414 res = RecapDecompose(m)
415 self.failUnless(res)
416 self.failUnless(len(res.GetLeaves())==3)
417 ks = res.GetLeaves().keys()
418 self.failUnless('[*]C1CC1' in ks)
419 self.failUnless('[*]C1CO1' in ks)
420 self.failUnless('[*]C1CN1' in ks)
421
422 m = Chem.MolFromSmiles('c1ccccc1N(C1NC1)C1OC1')
423 res = RecapDecompose(m)
424 self.failUnless(res)
425 self.failUnless(len(res.GetLeaves())==3)
426 ks = res.GetLeaves().keys()
427 self.failUnless('[*]c1ccccc1' in ks)
428 self.failUnless('[*]C1CO1' in ks)
429 self.failUnless('[*]C1CN1' in ks)
430
431 m = Chem.MolFromSmiles('c1ccccc1N(c1ncccc1)C1OC1')
432 res = RecapDecompose(m)
433 self.failUnless(res)
434 self.failUnless(len(res.GetLeaves())==3)
435 ks = res.GetLeaves().keys()
436 self.failUnless('[*]c1ccccc1' in ks)
437 self.failUnless('[*]c1ncccc1' in ks)
438 self.failUnless('[*]C1CO1' in ks)
439
440 m = Chem.MolFromSmiles('c1ccccc1N(c1ncccc1)c1ccco1')
441 res = RecapDecompose(m)
442 self.failUnless(res)
443 self.failUnless(len(res.GetLeaves())==3)
444 ks = res.GetLeaves().keys()
445 self.failUnless('[*]c1ccccc1' in ks)
446 self.failUnless('[*]c1ncccc1' in ks)
447 self.failUnless('[*]c1occc1' in ks)
448
449 m = Chem.MolFromSmiles('C1CCCCN1C1CC1')
450 res = RecapDecompose(m)
451 self.failUnless(res)
452 self.failUnless(len(res.GetLeaves())==2)
453 ks = res.GetLeaves().keys()
454 self.failUnless('[*]N1CCCCC1' in ks)
455 self.failUnless('[*]C1CC1' in ks)
456
457 m = Chem.MolFromSmiles('C1CCC2N1CC2')
458 res = RecapDecompose(m)
459 self.failUnless(res)
460 self.failUnless(len(res.GetLeaves())==0)
461
463 m = Chem.MolFromSmiles('C1CC1OC1OC1')
464 res = RecapDecompose(m)
465 self.failUnless(res)
466 self.failUnless(len(res.GetLeaves())==2)
467 ks = res.GetLeaves().keys()
468 self.failUnless('[*]C1CC1' in ks)
469 self.failUnless('[*]C1CO1' in ks)
470
471 m = Chem.MolFromSmiles('C1CCCCO1')
472 res = RecapDecompose(m)
473 self.failUnless(res)
474 self.failUnless(len(res.GetLeaves())==0)
475
476 m = Chem.MolFromSmiles('c1ccccc1OC1OC1')
477 res = RecapDecompose(m)
478 self.failUnless(res)
479 self.failUnless(len(res.GetLeaves())==2)
480 ks = res.GetLeaves().keys()
481 self.failUnless('[*]c1ccccc1' in ks)
482 self.failUnless('[*]C1CO1' in ks)
483
484 m = Chem.MolFromSmiles('c1ccccc1Oc1ncccc1')
485 res = RecapDecompose(m)
486 self.failUnless(res)
487 self.failUnless(len(res.GetLeaves())==2)
488 ks = res.GetLeaves().keys()
489 self.failUnless('[*]c1ccccc1' in ks)
490 self.failUnless('[*]c1ncccc1' in ks)
491
493 m = Chem.MolFromSmiles('ClC=CBr')
494 res = RecapDecompose(m)
495 self.failUnless(res)
496 self.failUnless(len(res.GetLeaves())==2)
497 ks = res.GetLeaves().keys()
498 self.failUnless('[*]CCl' in ks)
499 self.failUnless('[*]CBr' in ks)
500
501 m = Chem.MolFromSmiles('C1CC=CC1')
502 res = RecapDecompose(m)
503 self.failUnless(res)
504 self.failUnless(len(res.GetLeaves())==0)
505
507 m = Chem.MolFromSmiles('c1cccn1CCCC')
508 res = RecapDecompose(m)
509 self.failUnless(res)
510 self.failUnless(len(res.GetLeaves())==2)
511 ks = res.GetLeaves().keys()
512 self.failUnless('[*]n1cccc1' in ks)
513 self.failUnless('[*]CCCC' in ks)
514
515 m = Chem.MolFromSmiles('c1ccc2n1CCCC2')
516 res = RecapDecompose(m)
517 self.failUnless(res)
518 self.failUnless(len(res.GetLeaves())==0)
519
521 m = Chem.MolFromSmiles('C1CC(=O)N1CCCC')
522 res = RecapDecompose(m,onlyUseReactions=[8])
523 self.failUnless(res)
524 self.failUnless(len(res.GetLeaves())==2)
525 ks = res.GetLeaves().keys()
526 self.failUnless('[*]N1C(=O)CC1' in ks)
527 self.failUnless('[*]CCCC' in ks)
528
529 m = Chem.MolFromSmiles('O=C1CC2N1CCCC2')
530 res = RecapDecompose(m)
531 self.failUnless(res)
532 self.failUnless(len(res.GetLeaves())==0)
533
535 m = Chem.MolFromSmiles('c1ccccc1c1ncccc1')
536 res = RecapDecompose(m)
537 self.failUnless(res)
538 self.failUnless(len(res.GetLeaves())==2)
539 ks = res.GetLeaves().keys()
540 self.failUnless('[*]c1ccccc1' in ks)
541 self.failUnless('[*]c1ncccc1' in ks)
542
543 m = Chem.MolFromSmiles('c1ccccc1C1CC1')
544 res = RecapDecompose(m)
545 self.failUnless(res)
546 self.failUnless(len(res.GetLeaves())==0)
547
549 m = Chem.MolFromSmiles('c1cccn1c1ccccc1')
550 res = RecapDecompose(m)
551 self.failUnless(res)
552 self.failUnless(len(res.GetLeaves())==2)
553 ks = res.GetLeaves().keys()
554 self.failUnless('[*]n1cccc1' in ks)
555 self.failUnless('[*]c1ccccc1' in ks)
556
558 m = Chem.MolFromSmiles('CCCNS(=O)(=O)CC')
559 res = RecapDecompose(m)
560 self.failUnless(res)
561 self.failUnless(len(res.GetLeaves())==2)
562 ks = res.GetLeaves().keys()
563 self.failUnless('[*]NCCC' in ks)
564 self.failUnless('[*]S(=O)(=O)CC' in ks)
565
566 m = Chem.MolFromSmiles('c1cccn1S(=O)(=O)CC')
567 res = RecapDecompose(m)
568 self.failUnless(res)
569 self.failUnless(len(res.GetLeaves())==2)
570 ks = res.GetLeaves().keys()
571 self.failUnless('[*]n1cccc1' in ks)
572 self.failUnless('[*]S(=O)(=O)CC' in ks)
573
574 m = Chem.MolFromSmiles('C1CNS(=O)(=O)CC1')
575 res = RecapDecompose(m)
576 self.failUnless(res)
577 self.failUnless(len(res.GetLeaves())==0)
578
580 m = Chem.MolFromSmiles('c1ccccc1n1cccc1')
581 res = RecapDecompose(m)
582 self.failUnless(res)
583 self.failUnless(len(res.GetLeaves())==2)
584 m = Chem.MolFromSmiles('c1ccccc1[n+]1ccccc1')
585 res = RecapDecompose(m)
586 self.failUnless(res)
587 self.failUnless(len(res.GetLeaves())==0)
588
589 m = Chem.MolFromSmiles('C1CC1NC(=O)CC')
590 res = RecapDecompose(m)
591 self.failUnless(res)
592 self.failUnless(len(res.GetLeaves())==2)
593 m = Chem.MolFromSmiles('C1CC1[NH+]C(=O)CC')
594 res = RecapDecompose(m)
595 self.failUnless(res)
596 self.failUnless(len(res.GetLeaves())==0)
597
598 m = Chem.MolFromSmiles('C1CC1NC(=O)NC1CCC1')
599 res = RecapDecompose(m)
600 self.failUnless(res)
601 self.failUnless(len(res.GetLeaves())==2)
602 m = Chem.MolFromSmiles('C1CC1[NH+]C(=O)[NH+]C1CCC1')
603 res = RecapDecompose(m)
604 self.failUnless(res)
605 self.failUnless(len(res.GetLeaves())==0)
606
607 unittest.main()
608
| Trees | Indices | Help |
|
|---|
| Generated by Epydoc 3.0.1 on Sat Jan 7 08:30:01 2012 | http://epydoc.sourceforge.net |