Package rdkit :: Package Chem :: Package AtomPairs :: Module Sheridan
[hide private]
[frames] | no frames]

Module Sheridan

source code

Contains an implementation of Physicochemical property fingerprints, as
described in:
Kearsley, S. K. et al.
"Chemical Similarity Using Physiochemical Property Descriptors."
J. Chem.Inf. Model. 36, 118-127 (1996)

Functions [hide private]
 
_readPattyDefs(fname=os.path.join(RDConfig.RDDataDir,'SmartsLib','patty_rules.txt')) source code
 
AssignPattyTypes(mol, defns=None) source code
 
GetBPFingerprint(mol, fpfn=GetAtomPairFingerprint)
>>> from rdkit import Chem...
source code
 
GetBTFingerprint(mol, fpfn=GetTopologicalTorsionFingerprint)
>>> from rdkit import Chem...
source code
 
_test() source code
Variables [hide private]
  numPathBits = rdMolDescriptors.AtomPairsParameters.numPathBits
  _maxPathLen = 1 << numPathBits-1
  numFpBits = numPathBits+ 2* rdMolDescriptors.AtomPairsParamete...
  fpLen = 1 << numFpBits
  _pattyDefs = None
hash(x)
  typMap = dict(CAT= 1, ANI= 2, POL= 3, DON= 4, ACC= 5, HYD= 6, ...

Imports: IntSparseIntVect, Chem, rdMolDescriptors, DataStructs, GetAtomPairFingerprint, GetTopologicalTorsionFingerprint, os, re, RDConfig


Function Details [hide private]

AssignPattyTypes(mol, defns=None)

source code 


>>> from rdkit import Chem
>>> AssignPattyTypes(Chem.MolFromSmiles('OCC(=O)O'))
['POL', 'HYD', 'OTH', 'ANI', 'ANI']

GetBPFingerprint(mol, fpfn=GetAtomPairFingerprint)

source code 

>>> from rdkit import Chem
>>> fp = GetBPFingerprint(Chem.MolFromSmiles('OCC(=O)O'))
>>> fp.GetTotalVal()
10
>>> nze=fp.GetNonzeroElements()
>>> sorted([(k,v) for k,v in nze.items()])
[(32834, 1), (49219, 2), (98370, 2), (98401, 1), (114753, 2), (114786, 1), (114881, 1)]

GetBTFingerprint(mol, fpfn=GetTopologicalTorsionFingerprint)

source code 

>>> from rdkit import Chem
>>> mol = Chem.MolFromSmiles('OCC(N)O')
>>> AssignPattyTypes(mol)
['POL', 'HYD', 'HYD', 'CAT', 'POL']
>>> fp = GetBTFingerprint(mol)
>>> fp.GetTotalVal()
2
>>> nze=fp.GetNonzeroElements()
>>> sorted([(k,v) for k,v in nze.items()])
[(538446850L, 1), (538446852L, 1)]


Variables Details [hide private]

numFpBits

Value:
numPathBits+ 2* rdMolDescriptors.AtomPairsParameters.codeSize

typMap

Value:
dict(CAT= 1, ANI= 2, POL= 3, DON= 4, ACC= 5, HYD= 6, OTH= 7)