def pyScorePair(at1, at2, dist, atomCodes=None): """ Returns a score for an individual atom pair. >>> m = Chem.MolFromSmiles('CCCCC') >>> c1 = Utils.GetAtomCode(m.GetAtomWithIdx(0)) >>> c2 = Utils.GetAtomCode(m.GetAtomWithIdx(1)) >>> c3 = Utils.GetAtomCode(m.GetAtomWithIdx(2)) >>> t = 1 | min(c1,c2)<<numPathBits | max(c1,c2)<<(rdMolDescriptors.AtomPairsParameters.codeSize+numPathBits) >>> pyScorePair(m.GetAtomWithIdx(0),m.GetAtomWithIdx(1),1)==t 1 >>> pyScorePair(m.GetAtomWithIdx(1),m.GetAtomWithIdx(0),1)==t 1 >>> t = 2 | min(c1,c3)<<numPathBits | max(c1,c3)<<(rdMolDescriptors.AtomPairsParameters.codeSize+numPathBits) >>> pyScorePair(m.GetAtomWithIdx(0),m.GetAtomWithIdx(2),2)==t 1 >>> pyScorePair(m.GetAtomWithIdx(0),m.GetAtomWithIdx(2),2, ... atomCodes=(Utils.GetAtomCode(m.GetAtomWithIdx(0)),Utils.GetAtomCode(m.GetAtomWithIdx(2))))==t 1 """ if not atomCodes: code1 = Utils.GetAtomCode(at1) code2 = Utils.GetAtomCode(at2) else: code1, code2 = atomCodes accum = int(dist) % _maxPathLen accum |= min(code1, code2) << numPathBits accum |= max(code1, code2) << ( rdMolDescriptors.AtomPairsParameters.codeSize + numPathBits) return accum
def pyScorePath(mol, path, size, atomCodes=None): """ Returns a score for an individual path. >>> from rdkit import Chem >>> m = Chem.MolFromSmiles('CCCCC') >>> c1 = Utils.GetAtomCode(m.GetAtomWithIdx(0),1) >>> c2 = Utils.GetAtomCode(m.GetAtomWithIdx(1),2) >>> c3 = Utils.GetAtomCode(m.GetAtomWithIdx(2),2) >>> c4 = Utils.GetAtomCode(m.GetAtomWithIdx(3),1) >>> t = c1 | (c2 << rdMolDescriptors.AtomPairsParameters.codeSize) | (c3 << (rdMolDescriptors.AtomPairsParameters.codeSize*2)) | (c4 << (rdMolDescriptors.AtomPairsParameters.codeSize*3)) >>> pyScorePath(m,(0,1,2,3),4)==t 1 The scores are path direction independent: >>> pyScorePath(m,(3,2,1,0),4)==t 1 >>> m = Chem.MolFromSmiles('C=CC(=O)O') >>> c1 = Utils.GetAtomCode(m.GetAtomWithIdx(0),1) >>> c2 = Utils.GetAtomCode(m.GetAtomWithIdx(1),2) >>> c3 = Utils.GetAtomCode(m.GetAtomWithIdx(2),2) >>> c4 = Utils.GetAtomCode(m.GetAtomWithIdx(4),1) >>> t = c1 | (c2 << rdMolDescriptors.AtomPairsParameters.codeSize) | (c3 << (rdMolDescriptors.AtomPairsParameters.codeSize*2)) | (c4 << (rdMolDescriptors.AtomPairsParameters.codeSize*3)) >>> pyScorePath(m,(0,1,2,4),4)==t 1 """ codes = [None] * size for i in range(size): if i == 0 or i == (size - 1): sub = 1 else: sub = 2 if not atomCodes: codes[i] = Utils.GetAtomCode(mol.GetAtomWithIdx(path[i]), sub) else: base = atomCodes[path[i]] codes[i] = base - sub # "canonicalize" the code vector: beg = 0 end = len(codes) - 1 while (beg < end): if codes[beg] > codes[end]: codes.reverse() break elif codes[beg] == codes[end]: beg += 1 end -= 1 else: break accum = 0 for i in range(size): accum |= codes[i] << (rdMolDescriptors.AtomPairsParameters.codeSize * i) return accum