File: rdForceFields.cpp

package info (click to toggle)
rdkit 201809.1%2Bdfsg-6
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 123,688 kB
  • sloc: cpp: 230,509; python: 70,501; java: 6,329; ansic: 5,427; sql: 1,899; yacc: 1,739; lex: 1,243; makefile: 445; xml: 229; fortran: 183; sh: 123; cs: 93
file content (436 lines) | stat: -rw-r--r-- 17,619 bytes parent folder | download
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
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
// $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 <RDBoost/python.h>
#include <GraphMol/GraphMol.h>
#include <RDBoost/Wrap.h>

#include <ForceField/ForceField.h>
#include <ForceField/Wrap/PyForceField.h>
#include <ForceField/UFF/Params.h>
#include <GraphMol/ForceFieldHelpers/UFF/AtomTyper.h>
#include <GraphMol/ForceFieldHelpers/UFF/Builder.h>
#include <GraphMol/ForceFieldHelpers/UFF/UFF.h>
#include <GraphMol/ForceFieldHelpers/MMFF/AtomTyper.h>
#include <GraphMol/ForceFieldHelpers/MMFF/Builder.h>
#include <GraphMol/ForceFieldHelpers/MMFF/MMFF.h>

namespace python = boost::python;

namespace RDKit {
int UFFHelper(ROMol &mol, int maxIters, double vdwThresh, int confId,
              bool ignoreInterfragInteractions) {
  NOGIL gil;
  return UFF::UFFOptimizeMolecule(mol, maxIters, vdwThresh, confId,
                                  ignoreInterfragInteractions)
      .first;
}
python::object UFFConfsHelper(ROMol &mol, int numThreads, int maxIters,
                              double vdwThresh, int confId,
                              bool ignoreInterfragInteractions) {
  RDUNUSED_PARAM(confId);  // XXX FIX ME?
  std::vector<std::pair<int, double>> res;
  {
    NOGIL gil;
    UFF::UFFOptimizeMoleculeConfs(mol, res, numThreads, maxIters, vdwThresh,
                                  ignoreInterfragInteractions);
  }
  python::list pyres;
  for (auto &itm : res) {
    pyres.append(python::make_tuple(itm.first, itm.second));
  }
  return pyres;
}

python::object MMFFConfsHelper(ROMol &mol, int numThreads, int maxIters,
                               std::string mmffVariant, double nonBondedThresh,
                               int confId, bool ignoreInterfragInteractions) {
  RDUNUSED_PARAM(confId);  // Fix me?
  std::vector<std::pair<int, double>> res;
  {
    NOGIL gil;
    MMFF::MMFFOptimizeMoleculeConfs(mol, res, numThreads, maxIters, mmffVariant,
                                    nonBondedThresh,
                                    ignoreInterfragInteractions);
  }
  python::list pyres;
  for (auto &itm : res) {
    pyres.append(python::make_tuple(itm.first, itm.second));
  }
  return pyres;
}

ForceFields::PyForceField *UFFGetMoleculeForceField(
    ROMol &mol, double vdwThresh = 10.0, int confId = -1,
    bool ignoreInterfragInteractions = true) {
  ForceFields::ForceField *ff = UFF::constructForceField(
      mol, vdwThresh, confId, ignoreInterfragInteractions);
  auto *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()) {
    NOGIL gil;
    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) {
  auto *mmffMolProperties =
      new MMFF::MMFFMolProperties(mol, mmffVariant, mmffVerbosity);
  ForceFields::PyMMFFMolProperties *pyMP = nullptr;

  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 = nullptr;
  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(const ROMol &mol) {
  ROMol molCopy(mol);
  MMFF::MMFFMolProperties mmffMolProperties(molCopy);

  return mmffMolProperties.isValid();
}
};

namespace ForceFields {
PyObject *getUFFBondStretchParams(const RDKit::ROMol &mol,
                                  const unsigned int idx1,
                                  const unsigned int idx2) {
  PyObject *res = nullptr;
  ForceFields::UFF::UFFBond uffBondStretchParams;
  if (RDKit::UFF::getUFFBondStretchParams(mol, idx1, idx2,
                                          uffBondStretchParams)) {
    res = PyTuple_New(2);
    PyTuple_SetItem(res, 0, PyFloat_FromDouble(uffBondStretchParams.kb));
    PyTuple_SetItem(res, 1, PyFloat_FromDouble(uffBondStretchParams.r0));
  }
  return res;
};

PyObject *getUFFAngleBendParams(const RDKit::ROMol &mol,
                                const unsigned int idx1,
                                const unsigned int idx2,
                                const unsigned int idx3) {
  PyObject *res = nullptr;
  ForceFields::UFF::UFFAngle uffAngleBendParams;
  if (RDKit::UFF::getUFFAngleBendParams(mol, idx1, idx2, idx3,
                                        uffAngleBendParams)) {
    res = PyTuple_New(2);
    PyTuple_SetItem(res, 0, PyFloat_FromDouble(uffAngleBendParams.ka));
    PyTuple_SetItem(res, 1, PyFloat_FromDouble(uffAngleBendParams.theta0));
  }
  return res;
};

PyObject *getUFFTorsionParams(const RDKit::ROMol &mol, const unsigned int idx1,
                              const unsigned int idx2, const unsigned int idx3,
                              const unsigned int idx4) {
  PyObject *res = nullptr;
  ForceFields::UFF::UFFTor uffTorsionParams;
  if (RDKit::UFF::getUFFTorsionParams(mol, idx1, idx2, idx3, idx4,
                                      uffTorsionParams)) {
    res = PyFloat_FromDouble(uffTorsionParams.V);
  }
  return res;
};

PyObject *getUFFInversionParams(const RDKit::ROMol &mol,
                                const unsigned int idx1,
                                const unsigned int idx2,
                                const unsigned int idx3,
                                const unsigned int idx4) {
  PyObject *res = nullptr;
  ForceFields::UFF::UFFInv uffInversionParams;
  if (RDKit::UFF::getUFFInversionParams(mol, idx1, idx2, idx3, idx4,
                                        uffInversionParams)) {
    res = PyFloat_FromDouble(uffInversionParams.K);
  }
  return res;
};

PyObject *getUFFVdWParams(const RDKit::ROMol &mol, const unsigned int idx1,
                          const unsigned int idx2) {
  PyObject *res = nullptr;
  ForceFields::UFF::UFFVdW uffVdWParams;
  if (RDKit::UFF::getUFFVdWParams(mol, idx1, idx2, uffVdWParams)) {
    res = PyTuple_New(2);
    PyTuple_SetItem(res, 0, PyFloat_FromDouble(uffVdWParams.x_ij));
    PyTuple_SetItem(res, 1, PyFloat_FromDouble(uffVdWParams.D_ij));
  }
  return res;
};
}

BOOST_PYTHON_MODULE(rdForceFieldHelpers) {
  python::scope().attr("__doc__") =
      "Module containing functions to handle force fields";

  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::UFFHelper,
              (python::arg("self"), python::arg("maxIters") = 200,
               python::arg("vdwThresh") = 10.0, python::arg("confId") = -1,
               python::arg("ignoreInterfragInteractions") = true),
              docString.c_str());

  docString =
      "uses UFF to optimize all of a molecule's conformations\n\n\
 \n\
 ARGUMENTS:\n\n\
    - mol : the molecule of interest\n\
    - numThreads : the number of threads to use, only has an effect if the RDKit\n\
                   was built with thread support (defaults to 1)\n\
                   If set to zero, the max supported by the system will be used.\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: a list of (not_converged, energy) 2-tuples. \n\
     If not_converged is 0 the optimization converged for that conformer.\n\
\n";
  python::def("UFFOptimizeMoleculeConfs", RDKit::UFFConfsHelper,
              (python::arg("self"), python::arg("numThreads") = 1,
               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());

  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());

  docString =
      "uses MMFF to optimize all of a molecule's conformations\n\n\
 \n\
 ARGUMENTS:\n\n\
    - mol : the molecule of interest\n\
    - numThreads : the number of threads to use, only has an effect if the RDKit\n\
                   was built with thread support (defaults to 1)\n\
                   If set to zero, the max supported by the system will be used.\n\
    - maxIters : the maximum number of iterations (defaults to 200)\n\
    - mmffVariant : \"MMFF94\" or \"MMFF94s\"\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: a list of (not_converged, energy) 2-tuples. \n\
    If not_converged is 0 the optimization converged for that conformer.\n\
\n";
  python::def(
      "MMFFOptimizeMoleculeConfs", RDKit::MMFFConfsHelper,
      (python::arg("self"), python::arg("numThreads") = 1,
       python::arg("maxIters") = 200, python::arg("mmffVariant") = "MMFF94",
       python::arg("nonBondedThresh") = 10.0, python::arg("confId") = -1,
       python::arg("ignoreInterfragInteractions") = true),
      docString.c_str());

  python::def(
      "GetUFFBondStretchParams", ForceFields::getUFFBondStretchParams,
      (python::arg("mol"), python::arg("idx1"), python::arg("idx2")),
      "Retrieves UFF bond stretch parameters for atoms with indexes idx1, idx2 "
      "as a (kb, r0) tuple, or None if no parameters could be found");
  python::def("GetUFFAngleBendParams", ForceFields::getUFFAngleBendParams,
              (python::arg("mol"), python::arg("idx1"), python::arg("idx2"),
               python::arg("idx3")),
              "Retrieves UFF angle bend parameters for atoms with indexes "
              "idx1, idx2, idx3 as a (ka, theta0) tuple, or None if no "
              "parameters could be found");
  python::def("GetUFFTorsionParams", ForceFields::getUFFTorsionParams,
              (python::arg("mol"), python::arg("idx1"), python::arg("idx2"),
               python::arg("idx3"), python::arg("idx4")),
              "Retrieves UFF torsion parameters for atoms "
              "with indexes idx1, idx2, idx3, idx4 as a V float value, or None "
              "if no parameters "
              "could be found");
  python::def("GetUFFInversionParams", ForceFields::getUFFInversionParams,
              (python::arg("mol"), python::arg("idx1"), python::arg("idx2"),
               python::arg("idx3"), python::arg("idx4")),
              "Retrieves UFF inversion parameters for atoms "
              "with indexes idx1, idx2, idx3, idx4 as a K float value, or None "
              "if no parameters "
              "could be found");
  python::def(
      "GetUFFVdWParams", ForceFields::getUFFVdWParams,
      (python::arg("mol"), python::arg("idx1"), python::arg("idx2")),
      "Retrieves UFF van der Waals parameters for atoms with indexes idx1, "
      "idx2 "
      "as a (x_ij, D_ij) tuple, or None if no parameters could be found");
}