1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255
|
// $Id$
//
// Copyright (C) 2004-2008 Greg Landrum and Rational Discovery LLC
//
// @@ All Rights Reserved @@
// This file is part of the RDKit.
// The contents are covered by the terms of the BSD license
// which is included in the file license.txt, found at the root
// of the RDKit source tree.
//
#include <boost/python.hpp>
#include <GraphMol/GraphMol.h>
#include <RDBoost/Wrap.h>
#include <ForceField/ForceField.h>
#include <ForceField/Wrap/PyForceField.h>
#include <GraphMol/ForceFieldHelpers/UFF/AtomTyper.h>
#include <GraphMol/ForceFieldHelpers/UFF/Builder.h>
#include <GraphMol/ForceFieldHelpers/MMFF/AtomTyper.h>
#include <GraphMol/ForceFieldHelpers/MMFF/Builder.h>
namespace python = boost::python;
namespace RDKit {
int UFFOptimizeMolecule(ROMol &mol, int maxIters=200,
double vdwThresh=10.0, int confId=-1,
bool ignoreInterfragInteractions=true ){
ForceFields::ForceField *ff=UFF::constructForceField(mol,vdwThresh, confId,
ignoreInterfragInteractions);
ff->initialize();
int res=ff->minimize(maxIters);
delete ff;
return res;
}
ForceFields::PyForceField *UFFGetMoleculeForceField(ROMol &mol,
double vdwThresh=10.0,
int confId=-1,
bool ignoreInterfragInteractions=true ){
ForceFields::ForceField *ff=UFF::constructForceField(mol,vdwThresh, confId,
ignoreInterfragInteractions);
ForceFields::PyForceField *res=new ForceFields::PyForceField(ff);
res->initialize();
return res;
}
bool UFFHasAllMoleculeParams(const ROMol &mol){
UFF::AtomicParamVect types;
bool foundAll;
boost::tie(types,foundAll)=UFF::getAtomTypes(mol);
return foundAll;
}
int MMFFOptimizeMolecule(ROMol &mol, std::string mmffVariant = "MMFF94",
int maxIters = 200, double nonBondedThresh = 100.0, int confId = -1,
bool ignoreInterfragInteractions = true)
{
int res = -1;
MMFF::MMFFMolProperties mmffMolProperties(mol, mmffVariant);
if (mmffMolProperties.isValid()) {
ForceFields::ForceField *ff = MMFF::constructForceField(mol,
&mmffMolProperties, nonBondedThresh, confId, ignoreInterfragInteractions);
ff->initialize();
res = ff->minimize(maxIters);
delete ff;
}
return res;
}
unsigned int SanitizeMMFFMol(ROMol &mol)
{
return MMFF::sanitizeMMFFMol((RWMol &)mol);
};
ForceFields::PyMMFFMolProperties *GetMMFFMolProperties
(ROMol &mol, std::string mmffVariant = "MMFF94",
unsigned int mmffVerbosity = MMFF::MMFF_VERBOSITY_NONE)
{
MMFF::MMFFMolProperties *mmffMolProperties =
new MMFF::MMFFMolProperties(mol, mmffVariant, mmffVerbosity);
ForceFields::PyMMFFMolProperties *pyMP = NULL;
if (mmffMolProperties && mmffMolProperties->isValid()) {
pyMP = new ForceFields::PyMMFFMolProperties(mmffMolProperties);
}
return pyMP;
}
ForceFields::PyForceField *MMFFGetMoleculeForceField(ROMol &mol,
ForceFields::PyMMFFMolProperties *pyMMFFMolProperties,
double nonBondedThresh = 100.0, int confId = -1,
bool ignoreInterfragInteractions = true)
{
ForceFields::PyForceField *pyFF = NULL;
boost::python::list res;
if (pyMMFFMolProperties) {
MMFF::MMFFMolProperties *mmffMolProperties =
&(*(pyMMFFMolProperties->mmffMolProperties));
ForceFields::ForceField *ff = MMFF::constructForceField(mol,
mmffMolProperties, nonBondedThresh,
confId, ignoreInterfragInteractions);
pyFF = new ForceFields::PyForceField(ff);
if (pyFF) {
pyFF->initialize();
}
}
return pyFF;
}
bool MMFFHasAllMoleculeParams(ROMol &mol)
{
MMFF::MMFFMolProperties mmffMolProperties(mol);
return mmffMolProperties.isValid();
}
}
BOOST_PYTHON_MODULE(rdForceFieldHelpers) {
python::scope().attr("__doc__") =
"Module containing functions to handle UFF force field"
;
std::string docString = "uses UFF to optimize a molecule's structure\n\n\
\n\
ARGUMENTS:\n\n\
- mol : the molecule of interest\n\
- maxIters : the maximum number of iterations (defaults to 200)\n\
- vdwThresh : used to exclude long-range van der Waals interactions\n\
(defaults to 10.0)\n\
- confId : indicates which conformer to optimize\n\
- ignoreInterfragInteractions : if true, nonbonded terms between\n\
fragments will not be added to the forcefield.\n\
\n\
RETURNS: 0 if the optimization converged, 1 if more iterations are required.\n\
\n";
python::def("UFFOptimizeMolecule", RDKit::UFFOptimizeMolecule,
(python::arg("self"),python::arg("maxIters")=200,
python::arg("vdwThresh")=10.0,python::arg("confId")=-1,
python::arg("ignoreInterfragInteractions")=true),
docString.c_str());
docString = "returns a UFF force field for a molecule\n\n\
\n\
ARGUMENTS:\n\n\
- mol : the molecule of interest\n\
- vdwThresh : used to exclude long-range van der Waals interactions\n\
(defaults to 10.0)\n\
- confId : indicates which conformer to optimize\n\
- ignoreInterfragInteractions : if true, nonbonded terms between\n\
fragments will not be added to the forcefield.\n\
\n";
python::def("UFFGetMoleculeForceField", RDKit::UFFGetMoleculeForceField,
(python::arg("mol"),python::arg("vdwThresh")=10.0,python::arg("confId")=-1,
python::arg("ignoreInterfragInteractions")=true),
python::return_value_policy<python::manage_new_object>(),
docString.c_str());
docString = "checks if UFF parameters are available for all of a molecule's atoms\n\n\
\n\
ARGUMENTS:\n\n\
- mol : the molecule of interest.\n\
\n";
python::def("UFFHasAllMoleculeParams", RDKit::UFFHasAllMoleculeParams,
(python::arg("mol")),
docString.c_str());
python::scope().attr("__doc__") =
"Module containing functions to handle MMFF force field"
;
docString = "uses MMFF to optimize a molecule's structure\n\n\
\n\
ARGUMENTS:\n\n\
- mol : the molecule of interest\n\
- mmffVariant : \"MMFF94\" or \"MMFF94s\"\n\
- maxIters : the maximum number of iterations (defaults to 200)\n\
- nonBondedThresh : used to exclude long-range non-bonded\n\
interactions (defaults to 100.0)\n\
- confId : indicates which conformer to optimize\n\
- ignoreInterfragInteractions : if true, nonbonded terms between\n\
fragments will not be added to the forcefield\n\
\n\
RETURNS: 0 if the optimization converged, -1 if the forcefield could\n\
not be set up, 1 if more iterations are required.\n\
\n";
python::def("MMFFOptimizeMolecule", RDKit::MMFFOptimizeMolecule,
(python::arg("mol"), python::arg("mmffVariant") = "MMFF94",
python::arg("maxIters") = 200,
python::arg("nonBondedThresh") = 100.0,
python::arg("confId") = -1,
python::arg("ignoreInterfragInteractions") = true),
docString.c_str());
docString = "sanitizes a molecule according to MMFF requirements.\n\n\
- mol : the molecule of interest.\n\
\n";
python::def("MMFFSanitizeMolecule", RDKit::SanitizeMMFFMol,
(python::arg("mol")),
docString.c_str());
docString = "returns a PyMMFFMolProperties object for a\n\
molecule, which is required by MMFFGetMoleculeForceField()\n\
and can be used to get/set MMFF properties\n\n\
\n\
ARGUMENTS:\n\n\
- mol : the molecule of interest\n\
- mmffVariant : \"MMFF94\" or \"MMFF94s\"\n\
(defaults to \"MMFF94\")\n\
- mmffVerbosity : 0: none; 1: low; 2: high (defaults to 0).\n\
\n";
python::def("MMFFGetMoleculeProperties", RDKit::GetMMFFMolProperties,
(python::arg("mol"), python::arg("mmffVariant") = "MMFF94",
python::arg("mmffVerbosity") = 0),
python::return_value_policy<python::manage_new_object>(),
docString.c_str());
docString = "returns a MMFF force field for a molecule\n\n\
\n\
ARGUMENTS:\n\n\
- mol : the molecule of interest\n\
- pyMMFFMolProperties : PyMMFFMolProperties object as returned\n\
by MMFFGetMoleculeProperties()\n\
- nonBondedThresh : used to exclude long-range non-bonded\n\
interactions (defaults to 100.0)\n\
- confId : indicates which conformer to optimize\n\
- ignoreInterfragInteractions : if true, nonbonded terms between\n\
fragments will not be added to the forcefield\n\
\n";
python::def("MMFFGetMoleculeForceField", RDKit::MMFFGetMoleculeForceField,
(python::arg("mol"),
python::arg("pyMMFFMolProperties"),
python::arg("nonBondedThresh") = 100.0,
python::arg("confId") = -1,
python::arg("ignoreInterfragInteractions") = true),
python::return_value_policy<python::manage_new_object>(),
docString.c_str());
docString = "checks if MMFF parameters are available for all of a molecule's atoms\n\n\
\n\
ARGUMENTS:\n\n\
- mol : the molecule of interest\n\
\n";
python::def("MMFFHasAllMoleculeParams", RDKit::MMFFHasAllMoleculeParams,
(python::arg("mol")),
docString.c_str());
}
|