rdkit.Chem.Descriptors module

rdkit.Chem.Descriptors.CalcMolDescriptors(mol, missingVal=None, silent=True)

calculate the full set of descriptors for a molecule

Parameters:
  • mol (RDKit molecule) –

  • missingVal (float, optional) – This will be used if a particular descriptor cannot be calculated

  • silent (bool, optional) – if True then exception messages from descriptors will be displayed

Returns:

A dictionary with decriptor names as keys and the descriptor values as values

Return type:

dict

rdkit.Chem.Descriptors.ExactMolWt(*x, **y)

The exact molecular weight of the molecule

>>> ExactMolWt(Chem.MolFromSmiles('CC'))
30.04...
>>> ExactMolWt(Chem.MolFromSmiles('[13CH3]C'))
31.05...
rdkit.Chem.Descriptors.FpDensityMorgan1(x)
rdkit.Chem.Descriptors.FpDensityMorgan2(x)
rdkit.Chem.Descriptors.FpDensityMorgan3(x)
rdkit.Chem.Descriptors.HeavyAtomMolWt(x)

The average molecular weight of the molecule ignoring hydrogens

>>> HeavyAtomMolWt(Chem.MolFromSmiles('CC'))
24.02...
>>> HeavyAtomMolWt(Chem.MolFromSmiles('[NH4+].[Cl-]'))
49.46
rdkit.Chem.Descriptors.MaxAbsPartialCharge(mol, force=False)
rdkit.Chem.Descriptors.MaxPartialCharge(mol, force=False)
rdkit.Chem.Descriptors.MinAbsPartialCharge(mol, force=False)
rdkit.Chem.Descriptors.MinPartialCharge(mol, force=False)
rdkit.Chem.Descriptors.MolWt(*x, **y)

The average molecular weight of the molecule

>>> MolWt(Chem.MolFromSmiles('CC'))
30.07
>>> MolWt(Chem.MolFromSmiles('[NH4+].[Cl-]'))
53.49...
rdkit.Chem.Descriptors.NumRadicalElectrons(mol)
The number of radical electrons the molecule has

(says nothing about spin state)

>>> NumRadicalElectrons(Chem.MolFromSmiles('CC'))
0
>>> NumRadicalElectrons(Chem.MolFromSmiles('C[CH3]'))
0
>>> NumRadicalElectrons(Chem.MolFromSmiles('C[CH2]'))
1
>>> NumRadicalElectrons(Chem.MolFromSmiles('C[CH]'))
2
>>> NumRadicalElectrons(Chem.MolFromSmiles('C[C]'))
3
rdkit.Chem.Descriptors.NumValenceElectrons(mol)

The number of valence electrons the molecule has

>>> NumValenceElectrons(Chem.MolFromSmiles('CC'))
14
>>> NumValenceElectrons(Chem.MolFromSmiles('C(=O)O'))
18
>>> NumValenceElectrons(Chem.MolFromSmiles('C(=O)[O-]'))
18
>>> NumValenceElectrons(Chem.MolFromSmiles('C(=O)'))
12
class rdkit.Chem.Descriptors.PropertyFunctor((object)self, (object)callback, (str)name, (str)version)

Bases: PythonPropertyFunctor

Creates a python based property function that can be added to the global property list. To use, subclass this class and override the __call__ method. Then create an instance and add it to the registry. The __call__ method should return a numeric value.

Example

class NumAtoms(Descriptors.PropertyFunctor):
def __init__(self):

Descriptors.PropertyFunctor.__init__(self, “NumAtoms”, “1.0.0”)

def __call__(self, mol):

return mol.GetNumAtoms()

numAtoms = NumAtoms() rdMolDescriptors.Properties.RegisterProperty(numAtoms)

C++ signature :

void __init__(_object*,_object*,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >,std::__cxx11::basic_string<char, std::char_traits<char>, std::allocator<char> >)

rdkit.Chem.Descriptors.setupAUTOCorrDescriptors()

Adds AUTOCORR descriptors to the default descriptor lists