Package rdkit :: Package Chem :: Package Pharm2D :: Module Generate
[hide private]
[frames] | no frames]

Source Code for Module rdkit.Chem.Pharm2D.Generate

  1  # $Id: Generate.py 1528 2010-09-26 17:04:37Z glandrum $ 
  2  # 
  3  # Copyright (C) 2002-2008 greg Landrum and Rational Discovery LLC 
  4  # 
  5  #   @@ All Rights Reserved @@ 
  6  #  This file is part of the RDKit. 
  7  #  The contents are covered by the terms of the BSD license 
  8  #  which is included in the file license.txt, found at the root 
  9  #  of the RDKit source tree. 
 10  # 
 11  """ generation of 2D pharmacophores 
 12   
 13  **Notes** 
 14   
 15    - The terminology for this gets a bit rocky, so here's a glossary of what 
 16      terms used here mean: 
 17   
 18        1) *N-point pharmacophore* a combination of N features along with 
 19           distances betwen them. 
 20   
 21        2) *N-point proto-pharmacophore*: a combination of N feature 
 22           definitions without distances.  Each N-point 
 23           proto-pharmacophore defines a manifold of potential N-point 
 24           pharmacophores. 
 25   
 26        3) *N-point scaffold*: a collection of the distances defining 
 27           an N-point pharmacophore without feature identities. 
 28   
 29    See Docs/Chem/Pharm2D.triangles.jpg for an illustration of the way 
 30    pharmacophores are broken into triangles and labelled. 
 31   
 32    See Docs/Chem/Pharm2D.signatures.jpg for an illustration of bit 
 33    numbering 
 34   
 35  """ 
 36  from rdkit.Chem.Pharm2D import Utils,SigFactory 
 37  from rdkit.RDLogger import logger 
 38  logger = logger() 
 39   
 40  _verbose = 0 
 41   
42 -def _ShortestPathsMatch(match,featureSet,sig,dMat,sigFactory):
43 """ Internal use only 44 45 """ 46 if _verbose: 47 print 'match:',match 48 nPts = len(match) 49 distsToCheck = Utils.nPointDistDict[nPts] 50 nDists = len(distsToCheck) 51 dist = [0]*nDists 52 bins = sigFactory.GetBins() 53 minD,maxD = bins[0][0],bins[-1][1] 54 55 for i in range(nDists): 56 pt0,pt1 = distsToCheck[i] 57 minSeen=maxD 58 for idx1 in match[pt0]: 59 for idx2 in match[pt1]: 60 minSeen=min(minSeen, dMat[idx1,idx2]) 61 if minSeen==0 or minSeen<minD: return 62 # FIX: this won't be an int if we're using the bond order. 63 d = int(minSeen) 64 # do a quick distance filter 65 if d == 0 or d < minD or d >= maxD: 66 return 67 dist[i] = d 68 69 idx = sigFactory.GetBitIdx(featureSet,dist,sortIndices=False) 70 if _verbose: 71 print '\t',dist,minD,maxD,idx 72 73 if sigFactory.useCounts: 74 sig[idx] = sig[idx]+1 75 else: 76 sig.SetBit(idx)
77 78
79 -def Gen2DFingerprint(mol,sigFactory,perms=None,dMat=None):
80 """ generates a 2D fingerprint for a molecule using the 81 parameters in _sig_ 82 83 **Arguments** 84 85 - mol: the molecule for which the signature should be generated 86 87 - sigFactory : the SigFactory object with signature parameters 88 NOTE: no preprocessing is carried out for _sigFactory_. 89 It *must* be pre-initialized. 90 91 - perms: (optional) a sequence of permutation indices limiting which 92 pharmacophore combinations are allowed 93 94 - dMat: (optional) the distance matrix to be used 95 96 """ 97 if not isinstance(sigFactory,SigFactory.SigFactory): 98 raise ValueError,'bad factory' 99 featFamilies=sigFactory.GetFeatFamilies() 100 if _verbose: 101 print '* feat famillies:',featFamilies 102 nFeats = len(featFamilies) 103 minCount = sigFactory.minPointCount 104 maxCount = sigFactory.maxPointCount 105 if maxCount>3: 106 logger.warning(' Pharmacophores with more than 3 points are not currently supported.\nSetting maxCount to 3.') 107 maxCount=3 108 109 # generate the molecule's distance matrix, if required 110 if dMat is None: 111 from rdkit import Chem 112 useBO = sigFactory.includeBondOrder 113 dMat = Chem.GetDistanceMatrix(mol,useBO) 114 115 # generate the permutations, if required 116 if perms is None: 117 perms = [] 118 for count in range(minCount,maxCount+1): 119 perms += Utils.GetIndexCombinations(nFeats,count) 120 121 # generate the matches: 122 featMatches = sigFactory.GetMolFeats(mol) 123 if _verbose: 124 print ' featMatches:',featMatches 125 126 sig = sigFactory.GetSignature() 127 for perm in perms: 128 # the permutation is a combination of feature indices 129 # defining the feature set for a proto-pharmacophore 130 featClasses=[0] 131 for i in range(1,len(perm)): 132 if perm[i]==perm[i-1]: 133 featClasses.append(featClasses[-1]) 134 else: 135 featClasses.append(featClasses[-1]+1) 136 137 # Get a set of matches at each index of 138 # the proto-pharmacophore. 139 matchPerms = [featMatches[x] for x in perm] 140 if _verbose: 141 print '\n->Perm: %s'%(str(perm)) 142 print ' matchPerms: %s'%(str(matchPerms)) 143 144 # Get all unique combinations of those possible matches: 145 matchesToMap=Utils.GetUniqueCombinations(matchPerms,featClasses) 146 for i,entry in enumerate(matchesToMap): 147 entry = [x[1] for x in entry] 148 matchesToMap[i]=entry 149 if _verbose: 150 print ' mtM:',matchesToMap 151 152 for match in matchesToMap: 153 if sigFactory.shortestPathsOnly: 154 _ShortestPathsMatch(match,perm,sig,dMat,sigFactory) 155 return sig
156