Package rdkit :: Package Chem :: Package Draw :: Module SimilarityMaps
[hide private]
[frames] | no frames]

Module SimilarityMaps

source code

Functions [hide private]
 
GetAtomicWeightsForFingerprint(refMol, probeMol, fpFunction, metric=<Boost.Python.function object at 0x2717270>)
Calculates the atomic weights for the probe molecule based on a fingerprint function and a metric.
source code
 
GetAtomicWeightsForModel(probeMol, fpFunction, predictionFunction)
Calculates the atomic weights for the probe molecule based on a fingerprint function and the prediction function of a ML model.
source code
 
GetStandardizedWeights(weights)
Normalizes the weights, such that the absolute maximum weight equals 1.0.
source code
 
GetSimilarityMapFromWeights(mol, weights, colorMap=cm.PiYG, scale=-1, size=(250, 250), sigma=None, coordScale=1.5, step=0.01, colors='k', contourLines=10, alpha=0.5, **kwargs)
Generates the similarity map for a molecule given the atomic weights.
source code
 
GetSimilarityMapForFingerprint(refMol, probeMol, fpFunction, metric=<Boost.Python.function object at 0x2717270>, **kwargs)
Generates the similarity map for a given reference and probe molecule, fingerprint function and similarity metric.
source code
 
GetSimilarityMapForModel(probeMol, fpFunction, predictionFunction, **kwargs)
Generates the similarity map for a given ML model and probe molecule, and fingerprint function.
source code
 
GetAPFingerprint(mol, atomId=-1, fpType='normal', nBits=2048, minLength=1, maxLength=30, nBitsPerEntry=4, **kwargs)
Calculates the atom pairs fingerprint with the torsions of atomId removed.
source code
 
GetTTFingerprint(mol, atomId=-1, fpType='normal', nBits=2048, targetSize=4, nBitsPerEntry=4, **kwargs)
Calculates the topological torsion fingerprint with the pairs of atomId removed.
source code
 
GetMorganFingerprint(mol, atomId=-1, radius=2, fpType='bv', nBits=2048, useFeatures=False, **kwargs)
Calculates the Morgan fingerprint with the environments of atomId removed.
source code
 
GetRDKFingerprint(mol, atomId=-1, fpType='bv', nBits=2048, minPath=1, maxPath=5, nBitsPerHash=2, **kwargs)
Calculates the RDKit fingerprint with the paths of atomId removed.
source code
Variables [hide private]
  apDict = {}
  ttDict = {}
  __package__ = 'rdkit.Chem.Draw'

Imports: copy, math, cm, numpy, Chem, DataStructs, Draw, rdMD, iteritems


Function Details [hide private]

GetAtomicWeightsForFingerprint(refMol, probeMol, fpFunction, metric=<Boost.Python.function object at 0x2717270>)

source code 

Calculates the atomic weights for the probe molecule
based on a fingerprint function and a metric.

Parameters:
  refMol -- the reference molecule
  probeMol -- the probe molecule
  fpFunction -- the fingerprint function
  metric -- the similarity metric

Note:
  If fpFunction needs additional parameters, use a lambda construct

GetAtomicWeightsForModel(probeMol, fpFunction, predictionFunction)

source code 

Calculates the atomic weights for the probe molecule based on
a fingerprint function and the prediction function of a ML model.

Parameters:
  probeMol -- the probe molecule
  fpFunction -- the fingerprint function
  predictionFunction -- the prediction function of the ML model

GetStandardizedWeights(weights)

source code 

Normalizes the weights,
such that the absolute maximum weight equals 1.0.

Parameters:
  weights -- the list with the atomic weights

GetSimilarityMapFromWeights(mol, weights, colorMap=cm.PiYG, scale=-1, size=(250, 250), sigma=None, coordScale=1.5, step=0.01, colors='k', contourLines=10, alpha=0.5, **kwargs)

source code 

Generates the similarity map for a molecule given the atomic weights.

Parameters:
  mol -- the molecule of interest
  colorMap -- the matplotlib color map scheme
  scale -- the scaling: scale < 0 -> the absolute maximum weight is used as maximum scale
                        scale = double -> this is the maximum scale
  size -- the size of the figure
  sigma -- the sigma for the Gaussians
  coordScale -- scaling factor for the coordinates
  step -- the step for calcAtomGaussian
  colors -- color of the contour lines
  contourLines -- if integer number N: N contour lines are drawn
                  if list(numbers): contour lines at these numbers are drawn
  alpha -- the alpha blending value for the contour lines
  kwargs -- additional arguments for drawing

GetSimilarityMapForFingerprint(refMol, probeMol, fpFunction, metric=<Boost.Python.function object at 0x2717270>, **kwargs)

source code 

Generates the similarity map for a given reference and probe molecule,
fingerprint function and similarity metric.

Parameters:
  refMol -- the reference molecule
  probeMol -- the probe molecule
  fpFunction -- the fingerprint function
  metric -- the similarity metric.
  kwargs -- additional arguments for drawing

GetSimilarityMapForModel(probeMol, fpFunction, predictionFunction, **kwargs)

source code 

Generates the similarity map for a given ML model and probe molecule,
and fingerprint function.

Parameters:
  probeMol -- the probe molecule
  fpFunction -- the fingerprint function
  predictionFunction -- the prediction function of the ML model
  kwargs -- additional arguments for drawing

GetAPFingerprint(mol, atomId=-1, fpType='normal', nBits=2048, minLength=1, maxLength=30, nBitsPerEntry=4, **kwargs)

source code 

Calculates the atom pairs fingerprint with the torsions of atomId removed.

Parameters:
  mol -- the molecule of interest
  atomId -- the atom to remove the pairs for (if -1, no pair is removed)
  fpType -- the type of AP fingerprint ('normal', 'hashed', 'bv')
  nBits -- the size of the bit vector (only for fpType='bv')
  minLength -- the minimum path length for an atom pair
  maxLength -- the maxmimum path length for an atom pair
  nBitsPerEntry -- the number of bits available for each pair

GetTTFingerprint(mol, atomId=-1, fpType='normal', nBits=2048, targetSize=4, nBitsPerEntry=4, **kwargs)

source code 

Calculates the topological torsion fingerprint with the pairs of atomId removed.

Parameters:
  mol -- the molecule of interest
  atomId -- the atom to remove the torsions for (if -1, no torsion is removed)
  fpType -- the type of TT fingerprint ('normal', 'hashed', 'bv')
  nBits -- the size of the bit vector (only for fpType='bv')
  minLength -- the minimum path length for an atom pair
  maxLength -- the maxmimum path length for an atom pair
  nBitsPerEntry -- the number of bits available for each torsion

any additional keyword arguments will be passed to the fingerprinting function.

GetMorganFingerprint(mol, atomId=-1, radius=2, fpType='bv', nBits=2048, useFeatures=False, **kwargs)

source code 

Calculates the Morgan fingerprint with the environments of atomId removed.

Parameters:
  mol -- the molecule of interest
  radius -- the maximum radius
  fpType -- the type of Morgan fingerprint: 'count' or 'bv'
  atomId -- the atom to remove the environments for (if -1, no environments is removed)
  nBits -- the size of the bit vector (only for fpType = 'bv')
  useFeatures -- if false: ConnectivityMorgan, if true: FeatureMorgan

any additional keyword arguments will be passed to the fingerprinting function.

GetRDKFingerprint(mol, atomId=-1, fpType='bv', nBits=2048, minPath=1, maxPath=5, nBitsPerHash=2, **kwargs)

source code 

Calculates the RDKit fingerprint with the paths of atomId removed.

Parameters:
  mol -- the molecule of interest
  atomId -- the atom to remove the paths for (if -1, no path is removed)
  fpType -- the type of RDKit fingerprint: 'bv'
  nBits -- the size of the bit vector
  minPath -- minimum path length
  maxPath -- maximum path length
  nBitsPerHash -- number of to set per path