| Trees | Indices | Help |
|
|---|
|
|
1 # $Id: Recap.py 826 2008-09-08 11:44:49Z 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 import Chem
62 from 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:1;+0;D2,D3]!@C(!@=O)!@[#7:2;+0;D2,D3]>>[*][#7:1].[#7:2][*]", # urea
68
69 "[C:1;!$(C([#7])[#7])](=!@[O:2])!@[#7:3;+0;!D1]>>[*][C:1]=[O:2].[*][#7:3]", # amide
70
71 "[C:1](=!@[O:2])!@[O:3;+0]>>[*][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:1;R;D3;+0]-!@[*: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:1;+0]-!@[C:2]>>[n:1][*].[C:2][*]", # aromatic nitrogen - aliphatic carbon
84
85 "[O:3]=[C:4]-@[N:1;+0]-!@[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:1;+0]-!@[c:2]>>[n:1][*].[*][c:2]", # aromatic nitrogen - aromatic carbon *NOTE* this is not part of the standard recap set.
90
91 "[#7:1;+0;D2,D3]-!@[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 pSmi = Chem.MolToSmiles(prod,1)
186 if minFragmentSize>0:
187 nDummies = pSmi.count('*')
188 if nats-nDummies<minFragmentSize:
189 seqOk=False
190 break
191 # don't forget after replacing dummy atoms to remove any empty
192 # branches:
193 elif pSmi.replace('[*]','').replace('()','') in ('','C','CC','CCC'):
194 seqOk=False
195 break
196 prod.pSmi = pSmi
197 if seqOk:
198 for nats,prod in prodSeq:
199 pSmi = prod.pSmi
200 #print '\t',nats,pSmi
201 if not allNodes.has_key(pSmi):
202 pNode = RecapHierarchyNode(prod)
203 pNode.smiles=pSmi
204 pNode.parents[nSmi]=weakref.proxy(node)
205 node.children[pSmi]=pNode
206 activePool[pSmi] = pNode
207 allNodes[pSmi]=pNode
208 else:
209 pNode=allNodes[pSmi]
210 pNode.parents[nSmi]=weakref.proxy(node)
211 node.children[pSmi]=pNode
212 #print ' >>an:',allNodes.keys()
213 return res
214
215
216 # ------- ------- ------- ------- ------- ------- ------- -------
217 # Begin testing code
218 if __name__=='__main__':
219 import unittest
222 m = Chem.MolFromSmiles('C1CC1Oc1ccccc1-c1ncc(OC)cc1')
223 res = RecapDecompose(m)
224 self.failUnless(res)
225 self.failUnless(len(res.children.keys())==4)
226 self.failUnless(len(res.GetAllChildren().keys())==5)
227 self.failUnless(len(res.GetLeaves().keys())==3)
229 m = Chem.MolFromSmiles('CCCOCCC')
230 res = RecapDecompose(m)
231 self.failUnless(res)
232 self.failUnless(res.children=={})
234 allNodes={}
235 m = Chem.MolFromSmiles('c1ccccc1-c1ncccc1')
236 res = RecapDecompose(m,allNodes=allNodes)
237 self.failUnless(res)
238 self.failUnless(len(res.children.keys())==2)
239 self.failUnless(len(allNodes.keys())==3)
240
241 m = Chem.MolFromSmiles('COc1ccccc1-c1ncccc1')
242 res = RecapDecompose(m,allNodes=allNodes)
243 self.failUnless(res)
244 self.failUnless(len(res.children.keys())==2)
245 # we get two more nodes from that:
246 self.failUnless(len(allNodes.keys())==5)
247 self.failUnless(allNodes.has_key('[*]c1ccccc1OC'))
248 self.failUnless(allNodes.has_key('[*]c1ccccc1'))
249
250 m = Chem.MolFromSmiles('C1CC1Oc1ccccc1-c1ncccc1')
251 res = RecapDecompose(m,allNodes=allNodes)
252 self.failUnless(res)
253 self.failUnless(len(res.children.keys())==4)
254 self.failUnless(len(allNodes.keys())==10)
255
257 m = Chem.MolFromSmiles('c1ccccc1OC(Oc1ccccc1)Oc1ccccc1')
258 res = RecapDecompose(m)
259 self.failUnless(res)
260 self.failUnless(len(res.GetLeaves())==2)
261 ks = res.GetLeaves().keys()
262 self.failIf('[*]C([*])[*]' in ks)
263 self.failUnless('[*]c1ccccc1' in ks)
264 self.failUnless('[*]C([*])Oc1ccccc1' in ks)
265
267 m = Chem.MolFromSmiles('C1CCCCN1CCCC')
268 res = RecapDecompose(m)
269 self.failUnless(res)
270 self.failUnless(len(res.GetLeaves())==2)
271 ks = res.GetLeaves().keys()
272 self.failUnless('[*]N1CCCCC1' in ks)
273 self.failUnless('[*]CCCC' in ks)
274
276 m = Chem.MolFromSmiles('CCCOCCC')
277 res = RecapDecompose(m)
278 self.failUnless(res)
279 self.failUnless(res.children=={})
280 res = RecapDecompose(m,minFragmentSize=3)
281 self.failUnless(res)
282 self.failUnless(len(res.GetLeaves())==1)
283 ks = res.GetLeaves().keys()
284 self.failUnless('[*]CCC' in ks)
285
286 m = Chem.MolFromSmiles('CCCOCC')
287 res = RecapDecompose(m,minFragmentSize=3)
288 self.failUnless(res)
289 self.failUnless(res.children=={})
290
291 m = Chem.MolFromSmiles('CCCOCCOC')
292 res = RecapDecompose(m,minFragmentSize=2)
293 self.failUnless(res)
294 self.failUnless(len(res.GetLeaves())==2)
295 ks = res.GetLeaves().keys()
296 self.failUnless('[*]CCC' in ks)
297 ks = res.GetLeaves().keys()
298 self.failUnless('[*]CCOC' in ks)
299
301 m = Chem.MolFromSmiles('C1CC1C(=O)NC1OC1')
302 res = RecapDecompose(m,onlyUseReactions=[1])
303 self.failUnless(res)
304 self.failUnless(len(res.GetLeaves())==2)
305 ks = res.GetLeaves().keys()
306 self.failUnless('[*]C(=O)C1CC1' in ks)
307 self.failUnless('[*]NC1CO1' in ks)
308
309 m = Chem.MolFromSmiles('C1CC1C(=O)N(C)C1OC1')
310 res = RecapDecompose(m,onlyUseReactions=[1])
311 self.failUnless(res)
312 self.failUnless(len(res.GetLeaves())==2)
313 ks = res.GetLeaves().keys()
314 self.failUnless('[*]C(=O)C1CC1' in ks)
315 self.failUnless('[*]N(C)C1CO1' in ks)
316
317 m = Chem.MolFromSmiles('C1CC1C(=O)n1cccc1')
318 res = RecapDecompose(m,onlyUseReactions=[1])
319 self.failUnless(res)
320 self.failUnless(len(res.GetLeaves())==2)
321 ks = res.GetLeaves().keys()
322 self.failUnless('[*]C(=O)C1CC1' in ks)
323 self.failUnless('[*]n1cccc1' in ks)
324
325 m = Chem.MolFromSmiles('C1CC1C(=O)CC1OC1')
326 res = RecapDecompose(m,onlyUseReactions=[1])
327 self.failUnless(res)
328 self.failUnless(len(res.GetLeaves())==0)
329
330 m = Chem.MolFromSmiles('C1CCC(=O)NC1')
331 res = RecapDecompose(m,onlyUseReactions=[1])
332 self.failUnless(res)
333 self.failUnless(len(res.GetLeaves())==0)
334
335 m = Chem.MolFromSmiles('CC(=O)NC')
336 res = RecapDecompose(m,onlyUseReactions=[1])
337 self.failUnless(res)
338 self.failUnless(len(res.GetLeaves())==2)
339 ks = res.GetLeaves().keys()
340
341 m = Chem.MolFromSmiles('CC(=O)N')
342 res = RecapDecompose(m,onlyUseReactions=[1])
343 self.failUnless(res)
344 self.failUnless(len(res.GetLeaves())==0)
345
346 m = Chem.MolFromSmiles('C(=O)NCCNC(=O)CC')
347 res = RecapDecompose(m,onlyUseReactions=[1])
348 self.failUnless(res)
349 self.failUnless(len(res.children)==4)
350 self.failUnless(len(res.GetLeaves())==3)
351
353 m = Chem.MolFromSmiles('C1CC1C(=O)OC1OC1')
354 res = RecapDecompose(m,onlyUseReactions=[2])
355 self.failUnless(res)
356 self.failUnless(len(res.GetLeaves())==2)
357 ks = res.GetLeaves().keys()
358 self.failUnless('[*]C(=O)C1CC1' in ks)
359 self.failUnless('[*]OC1CO1' in ks)
360
361 m = Chem.MolFromSmiles('C1CC1C(=O)CC1OC1')
362 res = RecapDecompose(m,onlyUseReactions=[2])
363 self.failUnless(res)
364 self.failUnless(len(res.GetLeaves())==0)
365
366 m = Chem.MolFromSmiles('C1CCC(=O)OC1')
367 res = RecapDecompose(m,onlyUseReactions=[2])
368 self.failUnless(res)
369 self.failUnless(len(res.GetLeaves())==0)
370
372 m = Chem.MolFromSmiles('C1CC1NC(=O)NC1OC1')
373 res = RecapDecompose(m,onlyUseReactions=[0])
374 self.failUnless(res)
375 self.failUnless(len(res.GetLeaves())==2)
376 ks = res.GetLeaves().keys()
377 self.failUnless('[*]NC1CC1' in ks)
378 self.failUnless('[*]NC1CO1' in ks)
379
380 m = Chem.MolFromSmiles('C1CC1NC(=O)N(C)C1OC1')
381 res = RecapDecompose(m,onlyUseReactions=[0])
382 self.failUnless(res)
383 self.failUnless(len(res.GetLeaves())==2)
384 ks = res.GetLeaves().keys()
385 self.failUnless('[*]NC1CC1' in ks)
386 self.failUnless('[*]N(C)C1CO1' in ks)
387
388 m = Chem.MolFromSmiles('C1CCNC(=O)NC1C')
389 res = RecapDecompose(m,onlyUseReactions=[0])
390 self.failUnless(res)
391 self.failUnless(len(res.GetLeaves())==0)
392
393 m = Chem.MolFromSmiles('c1cccn1C(=O)NC1OC1')
394 res = RecapDecompose(m,onlyUseReactions=[0])
395 self.failUnless(res)
396 self.failUnless(len(res.GetLeaves())==2)
397 ks = res.GetLeaves().keys()
398 self.failUnless('[*]n1cccc1' in ks)
399 self.failUnless('[*]NC1CO1' in ks)
400
401 m = Chem.MolFromSmiles('c1cccn1C(=O)n1c(C)ccc1')
402 res = RecapDecompose(m,onlyUseReactions=[0])
403 self.failUnless(res)
404 self.failUnless(len(res.GetLeaves())==2)
405 ks = res.GetLeaves().keys()
406 self.failUnless('[*]n1c(C)ccc1' in ks)
407
409 m = Chem.MolFromSmiles('C1CC1N(C1NC1)C1OC1')
410 res = RecapDecompose(m)
411 self.failUnless(res)
412 self.failUnless(len(res.GetLeaves())==3)
413 ks = res.GetLeaves().keys()
414 self.failUnless('[*]C1CC1' in ks)
415 self.failUnless('[*]C1CO1' in ks)
416 self.failUnless('[*]C1CN1' in ks)
417
418 m = Chem.MolFromSmiles('c1ccccc1N(C1NC1)C1OC1')
419 res = RecapDecompose(m)
420 self.failUnless(res)
421 self.failUnless(len(res.GetLeaves())==3)
422 ks = res.GetLeaves().keys()
423 self.failUnless('[*]c1ccccc1' in ks)
424 self.failUnless('[*]C1CO1' in ks)
425 self.failUnless('[*]C1CN1' in ks)
426
427 m = Chem.MolFromSmiles('c1ccccc1N(c1ncccc1)C1OC1')
428 res = RecapDecompose(m)
429 self.failUnless(res)
430 self.failUnless(len(res.GetLeaves())==3)
431 ks = res.GetLeaves().keys()
432 self.failUnless('[*]c1ccccc1' in ks)
433 self.failUnless('[*]c1ncccc1' in ks)
434 self.failUnless('[*]C1CO1' in ks)
435
436 m = Chem.MolFromSmiles('c1ccccc1N(c1ncccc1)c1ccco1')
437 res = RecapDecompose(m)
438 self.failUnless(res)
439 self.failUnless(len(res.GetLeaves())==3)
440 ks = res.GetLeaves().keys()
441 self.failUnless('[*]c1ccccc1' in ks)
442 self.failUnless('[*]c1ncccc1' in ks)
443 self.failUnless('[*]c1occc1' in ks)
444
445 m = Chem.MolFromSmiles('C1CCCCN1C1CC1')
446 res = RecapDecompose(m)
447 self.failUnless(res)
448 self.failUnless(len(res.GetLeaves())==2)
449 ks = res.GetLeaves().keys()
450 self.failUnless('[*]N1CCCCC1' in ks)
451 self.failUnless('[*]C1CC1' in ks)
452
453 m = Chem.MolFromSmiles('C1CCC2N1CC2')
454 res = RecapDecompose(m)
455 self.failUnless(res)
456 self.failUnless(len(res.GetLeaves())==0)
457
459 m = Chem.MolFromSmiles('C1CC1OC1OC1')
460 res = RecapDecompose(m)
461 self.failUnless(res)
462 self.failUnless(len(res.GetLeaves())==2)
463 ks = res.GetLeaves().keys()
464 self.failUnless('[*]C1CC1' in ks)
465 self.failUnless('[*]C1CO1' in ks)
466
467 m = Chem.MolFromSmiles('C1CCCCO1')
468 res = RecapDecompose(m)
469 self.failUnless(res)
470 self.failUnless(len(res.GetLeaves())==0)
471
472 m = Chem.MolFromSmiles('c1ccccc1OC1OC1')
473 res = RecapDecompose(m)
474 self.failUnless(res)
475 self.failUnless(len(res.GetLeaves())==2)
476 ks = res.GetLeaves().keys()
477 self.failUnless('[*]c1ccccc1' in ks)
478 self.failUnless('[*]C1CO1' in ks)
479
480 m = Chem.MolFromSmiles('c1ccccc1Oc1ncccc1')
481 res = RecapDecompose(m)
482 self.failUnless(res)
483 self.failUnless(len(res.GetLeaves())==2)
484 ks = res.GetLeaves().keys()
485 self.failUnless('[*]c1ccccc1' in ks)
486 self.failUnless('[*]c1ncccc1' in ks)
487
489 m = Chem.MolFromSmiles('ClC=CBr')
490 res = RecapDecompose(m)
491 self.failUnless(res)
492 self.failUnless(len(res.GetLeaves())==2)
493 ks = res.GetLeaves().keys()
494 self.failUnless('[*]CCl' in ks)
495 self.failUnless('[*]CBr' in ks)
496
497 m = Chem.MolFromSmiles('C1CC=CC1')
498 res = RecapDecompose(m)
499 self.failUnless(res)
500 self.failUnless(len(res.GetLeaves())==0)
501
503 m = Chem.MolFromSmiles('c1cccn1CCCC')
504 res = RecapDecompose(m)
505 self.failUnless(res)
506 self.failUnless(len(res.GetLeaves())==2)
507 ks = res.GetLeaves().keys()
508 self.failUnless('[*]n1cccc1' in ks)
509 self.failUnless('[*]CCCC' in ks)
510
511 m = Chem.MolFromSmiles('c1ccc2n1CCCC2')
512 res = RecapDecompose(m)
513 self.failUnless(res)
514 self.failUnless(len(res.GetLeaves())==0)
515
517 m = Chem.MolFromSmiles('C1CC(=O)N1CCCC')
518 res = RecapDecompose(m,onlyUseReactions=[8])
519 self.failUnless(res)
520 self.failUnless(len(res.GetLeaves())==2)
521 ks = res.GetLeaves().keys()
522 self.failUnless('[*]N1C(=O)CC1' in ks)
523 self.failUnless('[*]CCCC' in ks)
524
525 m = Chem.MolFromSmiles('O=C1CC2N1CCCC2')
526 res = RecapDecompose(m)
527 self.failUnless(res)
528 self.failUnless(len(res.GetLeaves())==0)
529
531 m = Chem.MolFromSmiles('c1ccccc1c1ncccc1')
532 res = RecapDecompose(m)
533 self.failUnless(res)
534 self.failUnless(len(res.GetLeaves())==2)
535 ks = res.GetLeaves().keys()
536 self.failUnless('[*]c1ccccc1' in ks)
537 self.failUnless('[*]c1ncccc1' in ks)
538
539 m = Chem.MolFromSmiles('c1ccccc1C1CC1')
540 res = RecapDecompose(m)
541 self.failUnless(res)
542 self.failUnless(len(res.GetLeaves())==0)
543
545 m = Chem.MolFromSmiles('c1cccn1c1ccccc1')
546 res = RecapDecompose(m)
547 self.failUnless(res)
548 self.failUnless(len(res.GetLeaves())==2)
549 ks = res.GetLeaves().keys()
550 self.failUnless('[*]n1cccc1' in ks)
551 self.failUnless('[*]c1ccccc1' in ks)
552
554 m = Chem.MolFromSmiles('CCCNS(=O)(=O)CC')
555 res = RecapDecompose(m)
556 self.failUnless(res)
557 self.failUnless(len(res.GetLeaves())==2)
558 ks = res.GetLeaves().keys()
559 self.failUnless('[*]NCCC' in ks)
560 self.failUnless('[*]S(=O)(=O)CC' in ks)
561
562 m = Chem.MolFromSmiles('c1cccn1S(=O)(=O)CC')
563 res = RecapDecompose(m)
564 self.failUnless(res)
565 self.failUnless(len(res.GetLeaves())==2)
566 ks = res.GetLeaves().keys()
567 self.failUnless('[*]n1cccc1' in ks)
568 self.failUnless('[*]S(=O)(=O)CC' in ks)
569
570 m = Chem.MolFromSmiles('C1CNS(=O)(=O)CC1')
571 res = RecapDecompose(m)
572 self.failUnless(res)
573 self.failUnless(len(res.GetLeaves())==0)
574
576 m = Chem.MolFromSmiles('c1ccccc1n1cccc1')
577 res = RecapDecompose(m)
578 self.failUnless(res)
579 self.failUnless(len(res.GetLeaves())==2)
580 m = Chem.MolFromSmiles('c1ccccc1[n+]1ccccc1')
581 res = RecapDecompose(m)
582 self.failUnless(res)
583 self.failUnless(len(res.GetLeaves())==0)
584
585 m = Chem.MolFromSmiles('C1CC1NC(=O)CC')
586 res = RecapDecompose(m)
587 self.failUnless(res)
588 self.failUnless(len(res.GetLeaves())==2)
589 m = Chem.MolFromSmiles('C1CC1[NH+]C(=O)CC')
590 res = RecapDecompose(m)
591 self.failUnless(res)
592 self.failUnless(len(res.GetLeaves())==0)
593
594 m = Chem.MolFromSmiles('C1CC1NC(=O)NC1CCC1')
595 res = RecapDecompose(m)
596 self.failUnless(res)
597 self.failUnless(len(res.GetLeaves())==2)
598 m = Chem.MolFromSmiles('C1CC1[NH+]C(=O)[NH+]C1CCC1')
599 res = RecapDecompose(m)
600 self.failUnless(res)
601 self.failUnless(len(res.GetLeaves())==0)
602
603 unittest.main()
604
| Trees | Indices | Help |
|
|---|
| Generated by Epydoc 3.0beta1 on Tue Oct 7 06:26:53 2008 | http://epydoc.sourceforge.net |