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)

The fingerprints can be accessed through the following functions:
- GetBPFingerprint
- GetBTFingerprint

Functions [hide private]
 
_readPattyDefs(fname='/scratch/RDKit_git/Data/SmartsLib/patty_rules.txt') source code
 
AssignPattyTypes(mol, defns=None) source code
 
GetBPFingerprint(mol, fpfn=<Boost.Python.function object at 0x2ddd4a0>)
>>> from rdkit import Chem...
source code
 
GetBTFingerprint(mol, fpfn=<Boost.Python.function object at 0x2ddda40>)
>>> from rdkit import Chem...
source code
 
_runDoctests(verbose=None) source code
Variables [hide private]
  numPathBits = 5
  _maxPathLen = 31
  numFpBits = 23
  fpLen = 8388608
  _pattyDefs = None
hash(x)
  typMap = {'ACC': 5, 'ANI': 2, 'CAT': 1, 'DON': 4, 'HYD': 6, 'O...
  __package__ = 'rdkit.Chem.AtomPairs'

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


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=<Boost.Python.function object at 0x2ddd4a0>)

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=<Boost.Python.function object at 0x2ddda40>)

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()])
[(538446850..., 1), (538446852..., 1)]


Variables Details [hide private]

typMap

Value:
{'ACC': 5, 'ANI': 2, 'CAT': 1, 'DON': 4, 'HYD': 6, 'OTH': 7, 'POL': 3}