rdkit.Chem.rdShapeHelpers module¶
Module containing functions to encode and compare the shapes of molecules
- rdkit.Chem.rdShapeHelpers.ComputeConfBox((Conformer)conf[, (AtomPairsParameters)trans=None[, (float)padding=2.0]]) tuple : ¶
Compute the lower and upper corners of a cuboid that will fit the conformer
- C++ signature :
boost::python::tuple ComputeConfBox(RDKit::Conformer [,boost::python::api::object=None [,double=2.0]])
- rdkit.Chem.rdShapeHelpers.ComputeConfDimsAndOffset((Conformer)conf[, (AtomPairsParameters)trans=None[, (float)padding=2.0]]) tuple : ¶
- Compute the size of the box that can fit the conformations, and offset
of the box from the origin
- C++ signature :
boost::python::tuple ComputeConfDimsAndOffset(RDKit::Conformer [,boost::python::api::object=None [,double=2.0]])
- rdkit.Chem.rdShapeHelpers.ComputeUnionBox((tuple)box1, (tuple)box2) tuple : ¶
- Compute the union of two boxes, so that all the points in both boxes are
contained in the new box
- C++ signature :
boost::python::tuple ComputeUnionBox(boost::python::tuple,boost::python::tuple)
- rdkit.Chem.rdShapeHelpers.EncodeShape((Mol)mol, (UniformGrid3D_)grid[, (int)confId=-1[, (AtomPairsParameters)trans=None[, (float)vdwScale=0.8[, (float)stepSize=0.25[, (int)maxLayers=-1[, (bool)ignoreHs=True]]]]]]) None : ¶
Encode the shape of a molecule (one of its conformer) onto a grid
ARGUMENTS:
mol : the molecule of interest
grid : grid onto which the encoding is written
confId : id of the conformation of interest on mol (defaults to the first one)
trans : any transformation that needs to be used to encode onto the grid (note the molecule remains unchanged)
- vdwScaleScaling factor for the radius of the atoms to determine the base radius
used in the encoding - grid points inside this sphere carry the maximum occupancy
- setpSizethickness of the layers outside the base radius, the occupancy value is decreased
from layer to layer from the maximum value
- maxLayersthe maximum number of layers - defaults to the number of bits
used per grid point - e.g. two bits per grid point will allow 3 layers
ignoreHs : when set, the contribution of Hs to the shape will be ignored
- C++ signature :
void EncodeShape(RDKit::ROMol,RDGeom::UniformGrid3D {lvalue} [,int=-1 [,boost::python::api::object=None [,double=0.8 [,double=0.25 [,int=-1 [,bool=True]]]]]])
- rdkit.Chem.rdShapeHelpers.ShapeProtrudeDist((Mol)mol1, (Mol)mol2[, (int)confId1=-1[, (int)confId2=-1[, (float)gridSpacing=0.5[, (DiscreteValueType)bitsPerPoint=rdkit.DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE[, (float)vdwScale=0.8[, (float)stepSize=0.25[, (int)maxLayers=-1[, (bool)ignoreHs=True[, (bool)allowReordering=True]]]]]]]]]) float : ¶
Compute the shape protrude distance between two molecule based on a predefined alignment
- ARGUMENTS:
mol1 : The first molecule of interest
mol2 : The second molecule of interest
confId1 : Conformer in the first molecule (defaults to first conformer)
confId2 : Conformer in the second molecule (defaults to first conformer)
gridSpacing : resolution of the grid used to encode the molecular shapes
- bitsPerPointnumber of bit used to encode the occupancy at each grid point
defaults to two bits per grid point
- vdwScaleScaling factor for the radius of the atoms to determine the base radius
used in the encoding - grid points inside this sphere carry the maximum occupancy
- stepSizethickness of the each layer outside the base radius, the occupancy value is decreased
from layer to layer from the maximum value
- maxLayersthe maximum number of layers - defaults to the number of bits
used per grid point - e.g. two bits per grid point will allow 3 layers
ignoreHs : when set, the contribution of Hs to the shape will be ignored
- allowReorderingwhen set, the order will be automatically updated so that the value calculated
is the protrusion of the smaller shape from the larger one.
- C++ signature :
double ShapeProtrudeDist(RDKit::ROMol,RDKit::ROMol [,int=-1 [,int=-1 [,double=0.5 [,RDKit::DiscreteValueVect::DiscreteValueType=rdkit.DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [,double=0.8 [,double=0.25 [,int=-1 [,bool=True [,bool=True]]]]]]]]])
- rdkit.Chem.rdShapeHelpers.ShapeTanimotoDist((Mol)mol1, (Mol)mol2[, (int)confId1=-1[, (int)confId2=-1[, (float)gridSpacing=0.5[, (DiscreteValueType)bitsPerPoint=rdkit.DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE[, (float)vdwScale=0.8[, (float)stepSize=0.25[, (int)maxLayers=-1[, (bool)ignoreHs=True]]]]]]]]) float : ¶
Compute the shape tanimoto distance between two molecule based on a predefined alignment
- ARGUMENTS:
mol1 : The first molecule of interest
mol2 : The second molecule of interest
confId1 : Conformer in the first molecule (defaults to first conformer)
confId2 : Conformer in the second molecule (defaults to first conformer)
gridSpacing : resolution of the grid used to encode the molecular shapes
- bitsPerPointnumber of bits used to encode the occupancy at each grid point
defaults to two bits per grid point
- vdwScaleScaling factor for the radius of the atoms to determine the base radius
used in the encoding - grid points inside this sphere carry the maximum occupancy
- stepSizethickness of the each layer outside the base radius, the occupancy value is decreased
from layer to layer from the maximum value
- maxLayersthe maximum number of layers - defaults to the number of bits
used per grid point - e.g. two bits per grid point will allow 3 layers
ignoreHs : when set, the contribution of Hs to the shape will be ignored
- C++ signature :
double ShapeTanimotoDist(RDKit::ROMol,RDKit::ROMol [,int=-1 [,int=-1 [,double=0.5 [,RDKit::DiscreteValueVect::DiscreteValueType=rdkit.DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [,double=0.8 [,double=0.25 [,int=-1 [,bool=True]]]]]]]])
- rdkit.Chem.rdShapeHelpers.ShapeTverskyIndex((Mol)mol1, (Mol)mol2, (float)alpha, (float)beta[, (int)confId1=-1[, (int)confId2=-1[, (float)gridSpacing=0.5[, (DiscreteValueType)bitsPerPoint=rdkit.DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE[, (float)vdwScale=0.8[, (float)stepSize=0.25[, (int)maxLayers=-1[, (bool)ignoreHs=True]]]]]]]]) float : ¶
Compute the shape tversky index between two molecule based on a predefined alignment
- ARGUMENTS:
mol1 : The first molecule of interest
mol2 : The second molecule of interest
alpha : first parameter of the Tversky index
beta : second parameter of the Tversky index
confId1 : Conformer in the first molecule (defaults to first conformer)
confId2 : Conformer in the second molecule (defaults to first conformer)
gridSpacing : resolution of the grid used to encode the molecular shapes
- bitsPerPointnumber of bits used to encode the occupancy at each grid point
defaults to two bits per grid point
- vdwScaleScaling factor for the radius of the atoms to determine the base radius
used in the encoding - grid points inside this sphere carry the maximum occupancy
- stepSizethickness of the each layer outside the base radius, the occupancy value is decreased
from layer to layer from the maximum value
- maxLayersthe maximum number of layers - defaults to the number of bits
used per grid point - e.g. two bits per grid point will allow 3 layers
ignoreHs : when set, the contribution of Hs to the shape will be ignored
- C++ signature :
double ShapeTverskyIndex(RDKit::ROMol,RDKit::ROMol,double,double [,int=-1 [,int=-1 [,double=0.5 [,RDKit::DiscreteValueVect::DiscreteValueType=rdkit.DataStructs.cDataStructs.DiscreteValueType.TWOBITVALUE [,double=0.8 [,double=0.25 [,int=-1 [,bool=True]]]]]]]])