Package Chem :: Module TemplateAlign
[hide private]
[frames] | no frames]

Source Code for Module Chem.TemplateAlign

 1  # $Id: TemplateAlign.py 2 2006-05-06 22:54:39Z glandrum $
 
 2  #
 
 3  #  Copyright (C) 2005-2006 Rational Discovery LLC
 
 4  #
 
 5  #   @@ All Rights Reserved  @@
 
 6  #
 
 7  import Chem 
 8  from Chem import rdDepictor 
 9  import Geometry 
10  
 
11 -def AlignMolToTemplate2D(mol,template,match=None,clearConfs=False, 12 templateConfId=-1,):
13 """ 14 Arguments: 15 16 - mol: the molecule to be aligned 17 - template: the template to align to 18 - match: If provided, this should be a sequence of integers 19 containing the indices of the atoms in mol that match 20 those in template. This is the result of calling: 21 mol.GetSubstructMatch(template) 22 - clearConfs: toggles removing any existing conformers on mol 23 24 Returns the confId of the conformer containing the depiction 25 26 >>> patt = Chem.MolFromSmiles('C1CC1') 27 >>> rdDepictor.Compute2DCoords(patt) 28 0 29 >>> mol = Chem.MolFromSmiles('OC1CC1CC1CCC1') 30 >>> rdDepictor.Compute2DCoords(mol) 31 0 32 >>> pc = patt.GetConformer(0) 33 >>> mc = mol.GetConformer(0) 34 35 We start out with the molecules not aligned: 36 >>> vs = [abs(pc.GetAtomPosition(i).x-mc.GetAtomPosition(i+1).x) for i in range(pc.GetNumAtoms())] 37 >>> [x<1e-4 for x in vs] 38 [False, False, False] 39 40 But then we can replace the conformer of mol: 41 >>> AlignMolToTemplate2D(mol,patt,clearConfs=True) 42 0 43 >>> mol.GetNumConformers() 44 1 45 >>> pc = patt.GetConformer(0) 46 >>> mc = mol.GetConformer(0) 47 >>> vs = [abs(pc.GetAtomPosition(i).x-mc.GetAtomPosition(i+1).x) for i in range(pc.GetNumAtoms())] 48 >>> [x<1e-4 for x in vs] 49 [True, True, True] 50 51 If we like, we can specify the atom map explicitly in order to align to the second 52 matching ring in the probe molecule: 53 >>> match = (5,6,7) 54 >>> AlignMolToTemplate2D(mol,patt,clearConfs=True,match=match) 55 0 56 >>> mol.GetNumConformers() 57 1 58 >>> pc = patt.GetConformer(0) 59 >>> mc = mol.GetConformer(0) 60 >>> vs = [abs(pc.GetAtomPosition(i).x-mc.GetAtomPosition(i+1).x) for i in range(pc.GetNumAtoms())] 61 >>> [x<1e-4 for x in vs] 62 [False, False, False] 63 >>> vs = [abs(pc.GetAtomPosition(i).x-mc.GetAtomPosition(i+5).x) for i in range(pc.GetNumAtoms())] 64 >>> [x<1e-4 for x in vs] 65 [True, True, True] 66 67 68 69 """ 70 if not match: 71 match = mol.GetSubstructMatch(template) 72 if not match: 73 raise ValueError,'no match between mol and template' 74 75 atomMap = {} 76 templateConf = template.GetConformer(templateConfId) 77 for i,idx in enumerate(match): 78 p = templateConf.GetAtomPosition(i) 79 atomMap[idx] = Geometry.Point2D(p.x,p.y) 80 molConfId = rdDepictor.Compute2DCoords(mol,clearConfs=clearConfs,coordMap=atomMap) 81 return molConfId
82 83 84 #------------------------------------ 85 # 86 # doctest boilerplate 87 #
88 -def _test():
89 import doctest,sys 90 return doctest.testmod(sys.modules["__main__"])
91 92 93 if __name__ == '__main__': 94 import sys 95 failed,tried = _test() 96 sys.exit(failed) 97