Classes | |
| class | CrippenParams |
| a class used to store Crippen parameters More... | |
| class | CrippenParamCollection |
| singleton class for retrieving Crippen parameters More... | |
Namespaces | |
| namespace | AtomPairs |
Functions | |
| void | getCrippenAtomContribs (const ROMol &mol, std::vector< double > &logpContribs, std::vector< double > &mrContribs, bool force=false) |
| void | CalcCrippenDescriptors (const ROMol &mol, double &logp, double &mr, bool includeHs=true, bool force=false) |
| generate Wildman-Crippen LogP and MR estimates for a molecule | |
| double | CalcAMW (const ROMol &mol, bool onlyHeavy=false) |
| double | getLabuteAtomContribs (const ROMol &mol, std::vector< double > &Vi, double &hContrib, bool includeHs=true, bool force=false) |
| calculates atomic contributions to Labute's Approximate Surface Area | |
| double | calcLabuteASA (const ROMol &mol, bool includeHs=true, bool force=false) |
| calculates Labute's Approximate Surface Area (ASA from MOE) | |
Variables | |
| const std::string | crippenVersion = "1.1.0" |
| const std::string | labuteASAVersion = "1.0.2" |
| double RDKit::Descriptors::CalcAMW | ( | const ROMol & | mol, | |
| bool | onlyHeavy = false | |||
| ) |
Calculates a molecule's molecular weight
| mol | the molecule of interest | |
| onlyHeavy | (optional) if this is true (the default is false), only heavy atoms will be included in the MW calculation |
| void RDKit::Descriptors::CalcCrippenDescriptors | ( | const ROMol & | mol, | |
| double & | logp, | |||
| double & | mr, | |||
| bool | includeHs = true, |
|||
| bool | force = false | |||
| ) |
generate Wildman-Crippen LogP and MR estimates for a molecule
Uses an atom-based scheme based on the values in the paper: S. A. Wildman and G. M. Crippen JCICS 39 868-873 (1999)
| mol | the molecule of interest | |
| logp | used to return the logp estimate | |
| mr | used to return the MR estimate | |
| includeHs | (optional) if this is true (the default), a copy of mol is made and Hs are added to it. If false, Hs that are not explicitly present in the graph will not be included. | |
| force | forces the value to be recalculated instead of pulled from the cache |
| double RDKit::Descriptors::calcLabuteASA | ( | const ROMol & | mol, | |
| bool | includeHs = true, |
|||
| bool | force = false | |||
| ) |
calculates Labute's Approximate Surface Area (ASA from MOE)
Definition from P. Labute's article in the Journal of the Chemical Computing Group and J. Mol. Graph. Mod. _18_ 464-477 (2000)
| mol | the molecule of interest | |
| includeHs | (optional) if this is true (the default), the contribution of H atoms to the ASA will be included. | |
| force | (optional) calculate the value even if it's cached. |
| void RDKit::Descriptors::getCrippenAtomContribs | ( | const ROMol & | mol, | |
| std::vector< double > & | logpContribs, | |||
| std::vector< double > & | mrContribs, | |||
| bool | force = false | |||
| ) |
generate atomic contributions to the Wildman-Crippen LogP and MR estimates for a molecule
Uses an atom-based scheme based on the values in the paper: S. A. Wildman and G. M. Crippen JCICS 39 868-873 (1999)
| mol | the molecule of interest | |
| logpContribs | used to return the logp contributions, must be equal in length to the number of atoms | |
| mrContribs | used to return the MR contributions, must be equal in length to the number of atoms | |
| force | forces the value to be recalculated instead of pulled from the cache |
| double RDKit::Descriptors::getLabuteAtomContribs | ( | const ROMol & | mol, | |
| std::vector< double > & | Vi, | |||
| double & | hContrib, | |||
| bool | includeHs = true, |
|||
| bool | force = false | |||
| ) |
calculates atomic contributions to Labute's Approximate Surface Area
Definition from P. Labute's article in the Journal of the Chemical Computing Group and J. Mol. Graph. Mod. _18_ 464-477 (2000)
| mol | the molecule of interest | |
| Vi | used to return the explict atom contribs | |
| hContrib | used to return the H contributions (if calculated) | |
| includeHs | (optional) if this is true (the default), the contribution of H atoms to the ASA will be included. | |
| force | (optional) calculate the values even if they are cached. |
| const std::string RDKit::Descriptors::crippenVersion = "1.1.0" |
| const std::string RDKit::Descriptors::labuteASAVersion = "1.0.2" |
1.5.3