| Trees | Indices | Help |
|
|---|
|
|
1 # $Id: Recap.py 594 2008-04-06 10:41:14Z 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']
46 >>> res.GetAllChildren().keys()
47 ['[*]c1ccccc1[*]', '[*]c1ccc(OC)cn1', '[*]C1CC1', '[*]c1ccccc1-c1ncc(OC)cc1']
48
49 To get the standard set of RECAP results, use GetLeaves():
50 >>> leaves=res.GetLeaves()
51 >>> leaves.keys()
52 ['[*]c1ccccc1[*]', '[*]c1ccc(OC)cn1', '[*]C1CC1']
53 >>> leaf = leaves['[*]C1CC1']
54 >>> leaf.mol
55 <Chem.rdchem.Mol object at 0x00CBE0F0>
56
57
58 """
59 import weakref
60 import Chem
61 from Chem import rdChemReactions as Reactions
62
63 # These are the definitions that will be applied to fragment molecules:
64 reactionDefs = (
65 "[#7:1;+0;D2,D3]!@C(!@=O)!@[#7:2;+0;D2,D3]>>[*][#7:1].[#7:2][*]", # urea
66
67 "[C:1](=!@[O:2])!@[#7:3;+0;!D1]>>[*][C:1]=[O:2].[*][#7:3]", # amide
68
69 "[C:1](=!@[O:2])!@[O:3;+0]>>[*][C:1]=[O:2].[O:3][*]", # ester
70
71 "[N;!D1;+0](-!@[*:1])-!@[*:2]>>[*][*:1].[*:2][*]", # amines
72 #"[N;!D1](!@[*:1])!@[*:2]>>[*][*:1].[*:2][*]", # amines
73
74 # again: what about aromatics?
75 "[#7:1;R;D3;+0]-!@[*:2]>>[*][#7:1].[*:2][*]", # cyclic amines
76
77 "[#6:1]-!@[O;+0]-!@[#6:2]>>[#6:1][*].[*][#6:2]", # ether
78
79 "[C:1]=!@[C:2]>>[C:1][*].[*][C:2]", # olefin
80
81 "[n:1;+0]-!@[C:2]>>[n:1][*].[C:2][*]", # aromatic nitrogen - aliphatic carbon
82
83 "O=C-@[N:1;+0]-!@[C:2]>>[N:1][*].[C:2][*]", # lactam nitrogen - aliphatic carbon
84
85 "[c:1]-!@[c:2]>>[c:1][*].[*][c:2]", # aromatic carbon - aromatic carbon
86
87 "[n:1;+0]-!@[c:2]>>[n:1][*].[*][c:2]", # aromatic nitrogen - aromatic carbon *NOTE* this is not part of the standard recap set.
88
89 "[#7:1;+0;D2,D3]-!@[S:2](=[O:3])=[O:4]>>[#7:1][*].[*][S:2](=[O:3])=[O:4]", # sulphonamide
90 )
91
92 reactions = tuple([Reactions.ReactionFromSmarts(x) for x in reactionDefs])
93
94
96 """ This class is used to hold the Recap hiearchy
97 """
98 mol=None
99 children=None
100 parents=None
101 smiles = None
106
108 " returns a dictionary, keyed by SMILES, of children "
109 res = {}
110 for smi,child in self.children.iteritems():
111 res[smi] = child
112 child._gacRecurse(res,terminalOnly=False)
113 return res
114
116 " returns a dictionary, keyed by SMILES, of leaf (terminal) nodes "
117 res = {}
118 for smi,child in self.children.iteritems():
119 if not len(child.children):
120 res[smi] = child
121 else:
122 child._gacRecurse(res,terminalOnly=True)
123 return res
124
126 """ returns all the nodes in the hierarchy tree that contain this
127 node as a child
128 """
129 if not self.parents:
130 res = [self]
131 else:
132 res = []
133 for p in self.parents.values():
134 for uP in p.getUltimateParents():
135 if uP not in res:
136 res.append(uP)
137 return res
138
140 for smi,child in self.children.iteritems():
141 if not terminalOnly or not len(child.children):
142 res[smi] = child
143 child._gacRecurse(res,terminalOnly=terminalOnly)
144
149
150
152 """ returns the recap decomposition for a molecule """
153 mSmi = Chem.MolToSmiles(mol,1)
154
155 if allNodes is None:
156 allNodes={}
157 if allNodes.has_key(mSmi):
158 return allNodes[mSmi]
159
160 res = RecapHierarchyNode(mol)
161 res.smiles =mSmi
162 activePool={mSmi:res}
163 allNodes[mSmi]=res
164 for rxnIdx,reaction in enumerate(reactions):
165 localRes = {}
166 #print '>',rxnIdx,len(activePool.keys())
167 missed = {}
168 while activePool:
169 nSmi = activePool.keys()[0]
170 node = activePool.pop(nSmi)
171 ps = reaction.RunReactants((node.mol,))
172 if not ps:
173 #print ' !',nSmi
174 localRes[nSmi]=node
175 else:
176 rxnApplied=False
177 #print ' .',nSmi
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 rxnApplied=True
199 for nats,prod in prodSeq:
200 pSmi = prod.pSmi
201 pNode = allNodes.get(pSmi,RecapHierarchyNode(prod))
202 pNode.smiles=pSmi
203 pNode.parents[nSmi]=weakref.proxy(node)
204 node.children[pSmi]=pNode
205 activePool[pSmi] = pNode
206 allNodes[pSmi]=pNode
207 if not rxnApplied:
208 localRes[nSmi]=node
209 activePool=localRes
210 return res
211
212
213 # ------- ------- ------- ------- ------- ------- ------- -------
214 # Begin testing code
215 if __name__=='__main__':
216 import unittest
219 m = Chem.MolFromSmiles('C1CC1Oc1ccccc1-c1ncc(OC)cc1')
220 res = RecapDecompose(m)
221 self.failUnless(res)
222 self.failUnless(len(res.children.keys())==2)
223 self.failUnless(len(res.GetAllChildren().keys())==4)
224 self.failUnless(len(res.GetLeaves().keys())==3)
226 m = Chem.MolFromSmiles('CCCOCCC')
227 res = RecapDecompose(m)
228 self.failUnless(res)
229 self.failUnless(res.children=={})
231 allNodes={}
232 m = Chem.MolFromSmiles('c1ccccc1-c1ncccc1')
233 res = RecapDecompose(m,allNodes=allNodes)
234 self.failUnless(res)
235 self.failUnless(len(res.children.keys())==2)
236 self.failUnless(len(allNodes.keys())==3)
237
238 m = Chem.MolFromSmiles('COc1ccccc1-c1ncccc1')
239 res = RecapDecompose(m,allNodes=allNodes)
240 self.failUnless(res)
241 self.failUnless(len(res.children.keys())==2)
242 # we get two more nodes from that:
243 self.failUnless(len(allNodes.keys())==5)
244 self.failUnless(allNodes.has_key('[*]c1ccccc1OC'))
245 self.failUnless(allNodes.has_key('[*]c1ccccc1'))
246
247 m = Chem.MolFromSmiles('C1CC1Oc1ccccc1-c1ncccc1')
248 res = RecapDecompose(m,allNodes=allNodes)
249 self.failUnless(res)
250 self.failUnless(len(res.children.keys())==2)
251 self.failUnless(len(allNodes.keys())==9)
252
254 m = Chem.MolFromSmiles('c1ccccc1OC(Oc1ccccc1)Oc1ccccc1')
255 res = RecapDecompose(m)
256 self.failUnless(res)
257 self.failUnless(len(res.GetLeaves())==2)
258 ks = res.GetLeaves().keys()
259 self.failIf('[*]C([*])[*]' in ks)
260 self.failUnless('[*]c1ccccc1' in ks)
261 self.failUnless('[*]C([*])Oc1ccccc1' in ks)
262
264 m = Chem.MolFromSmiles('C1CCCCN1CCCC')
265 res = RecapDecompose(m)
266 self.failUnless(res)
267 self.failUnless(len(res.GetLeaves())==2)
268 ks = res.GetLeaves().keys()
269 self.failUnless('[*]N1CCCCC1' in ks)
270 self.failUnless('[*]CCCC' in ks)
271
273 m = Chem.MolFromSmiles('CCCOCCC')
274 res = RecapDecompose(m)
275 self.failUnless(res)
276 self.failUnless(res.children=={})
277 res = RecapDecompose(m,minFragmentSize=3)
278 self.failUnless(res)
279 self.failUnless(len(res.GetLeaves())==1)
280 ks = res.GetLeaves().keys()
281 self.failUnless('[*]CCC' in ks)
282
283 m = Chem.MolFromSmiles('CCCOCC')
284 res = RecapDecompose(m,minFragmentSize=3)
285 self.failUnless(res)
286 self.failUnless(res.children=={})
287
288 m = Chem.MolFromSmiles('CCCOCCOC')
289 res = RecapDecompose(m,minFragmentSize=2)
290 self.failUnless(res)
291 self.failUnless(len(res.GetLeaves())==2)
292 ks = res.GetLeaves().keys()
293 self.failUnless('[*]CCC' in ks)
294 ks = res.GetLeaves().keys()
295 self.failUnless('[*]CCOC' in ks)
296
298 m = Chem.MolFromSmiles('C1CC1C(=O)NC1OC1')
299 res = RecapDecompose(m)
300 self.failUnless(res)
301 self.failUnless(len(res.GetLeaves())==2)
302 ks = res.GetLeaves().keys()
303 self.failUnless('[*]C(=O)C1CC1' in ks)
304 self.failUnless('[*]NC1CO1' in ks)
305
306 m = Chem.MolFromSmiles('C1CC1C(=O)N(C)C1OC1')
307 res = RecapDecompose(m)
308 self.failUnless(res)
309 self.failUnless(len(res.GetLeaves())==2)
310 ks = res.GetLeaves().keys()
311 self.failUnless('[*]C(=O)C1CC1' in ks)
312 self.failUnless('[*]N(C)C1CO1' in ks)
313
314 m = Chem.MolFromSmiles('C1CC1C(=O)n1cccc1')
315 res = RecapDecompose(m)
316 self.failUnless(res)
317 self.failUnless(len(res.GetLeaves())==2)
318 ks = res.GetLeaves().keys()
319 self.failUnless('[*]C(=O)C1CC1' in ks)
320 self.failUnless('[*]n1cccc1' in ks)
321
322 m = Chem.MolFromSmiles('C1CC1C(=O)CC1OC1')
323 res = RecapDecompose(m)
324 self.failUnless(res)
325 self.failUnless(len(res.GetLeaves())==0)
326
327 m = Chem.MolFromSmiles('C1CCC(=O)NC1')
328 res = RecapDecompose(m)
329 self.failUnless(res)
330 self.failUnless(len(res.GetLeaves())==0)
331
332 m = Chem.MolFromSmiles('CC(=O)NC')
333 res = RecapDecompose(m)
334 self.failUnless(res)
335 self.failUnless(len(res.GetLeaves())==2)
336 ks = res.GetLeaves().keys()
337
338 m = Chem.MolFromSmiles('CC(=O)N')
339 res = RecapDecompose(m)
340 self.failUnless(res)
341 self.failUnless(len(res.GetLeaves())==0)
342
343 m = Chem.MolFromSmiles('C(=O)NCCNC(=O)CC')
344 res = RecapDecompose(m)
345 self.failUnless(res)
346 self.failUnless(len(res.children)==4)
347 self.failUnless(len(res.GetLeaves())==3)
348
350 m = Chem.MolFromSmiles('C1CC1C(=O)OC1OC1')
351 res = RecapDecompose(m)
352 self.failUnless(res)
353 self.failUnless(len(res.GetLeaves())==2)
354 ks = res.GetLeaves().keys()
355 self.failUnless('[*]C(=O)C1CC1' in ks)
356 self.failUnless('[*]OC1CO1' in ks)
357
358 m = Chem.MolFromSmiles('C1CC1C(=O)CC1OC1')
359 res = RecapDecompose(m)
360 self.failUnless(res)
361 self.failUnless(len(res.GetLeaves())==0)
362
363 m = Chem.MolFromSmiles('C1CCC(=O)OC1')
364 res = RecapDecompose(m)
365 self.failUnless(res)
366 self.failUnless(len(res.GetLeaves())==0)
367
369 m = Chem.MolFromSmiles('C1CC1NC(=O)NC1OC1')
370 res = RecapDecompose(m)
371 self.failUnless(res)
372 self.failUnless(len(res.GetLeaves())==2)
373 ks = res.GetLeaves().keys()
374 self.failUnless('[*]NC1CC1' in ks)
375 self.failUnless('[*]NC1CO1' in ks)
376
377 m = Chem.MolFromSmiles('C1CC1NC(=O)N(C)C1OC1')
378 res = RecapDecompose(m)
379 self.failUnless(res)
380 self.failUnless(len(res.GetLeaves())==2)
381 ks = res.GetLeaves().keys()
382 self.failUnless('[*]NC1CC1' in ks)
383 self.failUnless('[*]N(C)C1CO1' in ks)
384
385 m = Chem.MolFromSmiles('C1CCNC(=O)NC1C')
386 res = RecapDecompose(m)
387 self.failUnless(res)
388 self.failUnless(len(res.GetLeaves())==0)
389
390 m = Chem.MolFromSmiles('c1cccn1C(=O)NC1OC1')
391 res = RecapDecompose(m)
392 self.failUnless(res)
393 self.failUnless(len(res.GetLeaves())==2)
394 ks = res.GetLeaves().keys()
395 self.failUnless('[*]n1cccc1' in ks)
396 self.failUnless('[*]NC1CO1' in ks)
397
398 m = Chem.MolFromSmiles('c1cccn1C(=O)n1c(C)ccc1')
399 res = RecapDecompose(m)
400 self.failUnless(res)
401 self.failUnless(len(res.GetLeaves())==2)
402 ks = res.GetLeaves().keys()
403 self.failUnless('[*]n1c(C)ccc1' in ks)
404
406 m = Chem.MolFromSmiles('C1CC1N(C1NC1)C1OC1')
407 res = RecapDecompose(m)
408 self.failUnless(res)
409 self.failUnless(len(res.GetLeaves())==3)
410 ks = res.GetLeaves().keys()
411 self.failUnless('[*]C1CC1' in ks)
412 self.failUnless('[*]C1CO1' in ks)
413 self.failUnless('[*]C1CN1' in ks)
414
415 m = Chem.MolFromSmiles('c1ccccc1N(C1NC1)C1OC1')
416 res = RecapDecompose(m)
417 self.failUnless(res)
418 self.failUnless(len(res.GetLeaves())==3)
419 ks = res.GetLeaves().keys()
420 self.failUnless('[*]c1ccccc1' in ks)
421 self.failUnless('[*]C1CO1' in ks)
422 self.failUnless('[*]C1CN1' in ks)
423
424 m = Chem.MolFromSmiles('c1ccccc1N(c1ncccc1)C1OC1')
425 res = RecapDecompose(m)
426 self.failUnless(res)
427 self.failUnless(len(res.GetLeaves())==3)
428 ks = res.GetLeaves().keys()
429 self.failUnless('[*]c1ccccc1' in ks)
430 self.failUnless('[*]c1ncccc1' in ks)
431 self.failUnless('[*]C1CO1' in ks)
432
433 m = Chem.MolFromSmiles('c1ccccc1N(c1ncccc1)c1ccco1')
434 res = RecapDecompose(m)
435 self.failUnless(res)
436 self.failUnless(len(res.GetLeaves())==3)
437 ks = res.GetLeaves().keys()
438 self.failUnless('[*]c1ccccc1' in ks)
439 self.failUnless('[*]c1ncccc1' in ks)
440 self.failUnless('[*]c1occc1' in ks)
441
442 m = Chem.MolFromSmiles('C1CCCCN1C1CC1')
443 res = RecapDecompose(m)
444 self.failUnless(res)
445 self.failUnless(len(res.GetLeaves())==2)
446 ks = res.GetLeaves().keys()
447 self.failUnless('[*]N1CCCCC1' in ks)
448 self.failUnless('[*]C1CC1' in ks)
449
450 m = Chem.MolFromSmiles('C1CCC2N1CC2')
451 res = RecapDecompose(m)
452 self.failUnless(res)
453 self.failUnless(len(res.GetLeaves())==0)
454
456 m = Chem.MolFromSmiles('C1CC1OC1OC1')
457 res = RecapDecompose(m)
458 self.failUnless(res)
459 self.failUnless(len(res.GetLeaves())==2)
460 ks = res.GetLeaves().keys()
461 self.failUnless('[*]C1CC1' in ks)
462 self.failUnless('[*]C1CO1' in ks)
463
464 m = Chem.MolFromSmiles('C1CCCCO1')
465 res = RecapDecompose(m)
466 self.failUnless(res)
467 self.failUnless(len(res.GetLeaves())==0)
468
469 m = Chem.MolFromSmiles('c1ccccc1OC1OC1')
470 res = RecapDecompose(m)
471 self.failUnless(res)
472 self.failUnless(len(res.GetLeaves())==2)
473 ks = res.GetLeaves().keys()
474 self.failUnless('[*]c1ccccc1' in ks)
475 self.failUnless('[*]C1CO1' in ks)
476
477 m = Chem.MolFromSmiles('c1ccccc1Oc1ncccc1')
478 res = RecapDecompose(m)
479 self.failUnless(res)
480 self.failUnless(len(res.GetLeaves())==2)
481 ks = res.GetLeaves().keys()
482 self.failUnless('[*]c1ccccc1' in ks)
483 self.failUnless('[*]c1ncccc1' in ks)
484
486 m = Chem.MolFromSmiles('ClC=CBr')
487 res = RecapDecompose(m)
488 self.failUnless(res)
489 self.failUnless(len(res.GetLeaves())==2)
490 ks = res.GetLeaves().keys()
491 self.failUnless('[*]CCl' in ks)
492 self.failUnless('[*]CBr' in ks)
493
494 m = Chem.MolFromSmiles('C1CC=CC1')
495 res = RecapDecompose(m)
496 self.failUnless(res)
497 self.failUnless(len(res.GetLeaves())==0)
498
500 m = Chem.MolFromSmiles('c1cccn1CCCC')
501 res = RecapDecompose(m)
502 self.failUnless(res)
503 self.failUnless(len(res.GetLeaves())==2)
504 ks = res.GetLeaves().keys()
505 self.failUnless('[*]n1cccc1' in ks)
506 self.failUnless('[*]CCCC' in ks)
507
508 m = Chem.MolFromSmiles('c1ccc2n1CCCC2')
509 res = RecapDecompose(m)
510 self.failUnless(res)
511 self.failUnless(len(res.GetLeaves())==0)
512
514 m = Chem.MolFromSmiles('C1CC(=O)N1CCCC')
515 res = RecapDecompose(m)
516 self.failUnless(res)
517 self.failUnless(len(res.GetLeaves())==2)
518 ks = res.GetLeaves().keys()
519 self.failUnless('[*]N1C(=O)CC1' in ks)
520 self.failUnless('[*]CCCC' in ks)
521
522 m = Chem.MolFromSmiles('O=C1CC2N1CCCC2')
523 res = RecapDecompose(m)
524 self.failUnless(res)
525 self.failUnless(len(res.GetLeaves())==0)
526
528 m = Chem.MolFromSmiles('c1ccccc1c1ncccc1')
529 res = RecapDecompose(m)
530 self.failUnless(res)
531 self.failUnless(len(res.GetLeaves())==2)
532 ks = res.GetLeaves().keys()
533 self.failUnless('[*]c1ccccc1' in ks)
534 self.failUnless('[*]c1ncccc1' in ks)
535
536 m = Chem.MolFromSmiles('c1ccccc1C1CC1')
537 res = RecapDecompose(m)
538 self.failUnless(res)
539 self.failUnless(len(res.GetLeaves())==0)
540
542 m = Chem.MolFromSmiles('c1cccn1c1ccccc1')
543 res = RecapDecompose(m)
544 self.failUnless(res)
545 self.failUnless(len(res.GetLeaves())==2)
546 ks = res.GetLeaves().keys()
547 self.failUnless('[*]n1cccc1' in ks)
548 self.failUnless('[*]c1ccccc1' in ks)
549
551 m = Chem.MolFromSmiles('CCCNS(=O)(=O)CC')
552 res = RecapDecompose(m)
553 self.failUnless(res)
554 self.failUnless(len(res.GetLeaves())==2)
555 ks = res.GetLeaves().keys()
556 self.failUnless('[*]NCCC' in ks)
557 self.failUnless('[*]S(=O)(=O)CC' in ks)
558
559 m = Chem.MolFromSmiles('c1cccn1S(=O)(=O)CC')
560 res = RecapDecompose(m)
561 self.failUnless(res)
562 self.failUnless(len(res.GetLeaves())==2)
563 ks = res.GetLeaves().keys()
564 self.failUnless('[*]n1cccc1' in ks)
565 self.failUnless('[*]S(=O)(=O)CC' in ks)
566
567 m = Chem.MolFromSmiles('C1CNS(=O)(=O)CC1')
568 res = RecapDecompose(m)
569 self.failUnless(res)
570 self.failUnless(len(res.GetLeaves())==0)
571
573 m = Chem.MolFromSmiles('c1ccccc1n1cccc1')
574 res = RecapDecompose(m)
575 self.failUnless(res)
576 self.failUnless(len(res.GetLeaves())==2)
577 m = Chem.MolFromSmiles('c1ccccc1[n+]1ccccc1')
578 res = RecapDecompose(m)
579 self.failUnless(res)
580 self.failUnless(len(res.GetLeaves())==0)
581
582 m = Chem.MolFromSmiles('C1CC1NC(=O)CC')
583 res = RecapDecompose(m)
584 self.failUnless(res)
585 self.failUnless(len(res.GetLeaves())==2)
586 m = Chem.MolFromSmiles('C1CC1[NH+]C(=O)CC')
587 res = RecapDecompose(m)
588 self.failUnless(res)
589 self.failUnless(len(res.GetLeaves())==0)
590
591 m = Chem.MolFromSmiles('C1CC1NC(=O)NC1CCC1')
592 res = RecapDecompose(m)
593 self.failUnless(res)
594 self.failUnless(len(res.GetLeaves())==2)
595 m = Chem.MolFromSmiles('C1CC1[NH+]C(=O)[NH+]C1CCC1')
596 res = RecapDecompose(m)
597 self.failUnless(res)
598 self.failUnless(len(res.GetLeaves())==0)
599
600 unittest.main()
601
| Trees | Indices | Help |
|
|---|
| Generated by Epydoc 3.0beta1 on Sat May 24 08:37:18 2008 | http://epydoc.sourceforge.net |