Fixed improper error message in openmm_wrapper.
[gromacs.git] / src / kernel / openmm_wrapper.cpp
blob39a6642687b94156a6d8ce33c016c56050560dac
1 #include <cmath>
2 #include <set>
3 #include <iostream>
4 #include <sstream>
5 #include <fstream>
6 #include <map>
7 #include <vector>
8 #include <cctype>
9 #include <algorithm>
11 using namespace std;
13 #include "OpenMM.h"
15 #include "gmx_fatal.h"
16 #include "typedefs.h"
17 #include "mdrun.h"
18 #include "physics.h"
19 #include "string2.h"
20 #include "gmx_gpu_utils.h"
21 #include "mtop_util.h"
22 #include "warninp.h"
24 #include "openmm_wrapper.h"
26 using namespace OpenMM;
28 /*! \cond */
29 #define MEM_ERR_MSG(str) \
30 "The %s-simulation GPU memory test detected errors. As memory errors would cause incorrect " \
31 "simulation results, gromacs has aborted execution.\n Make sure that your GPU's memory is not " \
32 "overclocked and that the device is properly cooled.\n", (str)
33 /*! \endcond */
35 /*!
36 * \brief Convert string to integer type.
37 * \param[in] s String to convert from.
38 * \param[in] ios_base Basefield format flag that takes any of the following I/O
39 * manipulators: dec, hex, oct.
40 * \param[out] Destination variable to convert to.
42 template <class T>
43 static bool from_string(T& t, const string& s, ios_base& (*f)(ios_base&))
45 istringstream iss(s);
46 return !(iss >> f >> t).fail();
49 /*!
50 * \brief Split string around a given delimiter.
51 * \param[in] s String to split.
52 * \param[in] delim Delimiter character that defines the boundaries of substring in \p s.
53 * \returns Vector of strings found in \p s.
55 static vector<string> split(const string &s, char delim)
57 vector<string> elems;
58 stringstream ss(s);
59 string item;
60 while (getline(ss, item, delim))
62 if (item.length() != 0)
63 elems.push_back(item);
65 return elems;
68 /*!
69 * \brief Split a string of the form "option=value" into "option" and "value" strings.
70 * This string corresponds to one option and the associated value from the option list
71 * in the mdrun -device argument.
73 * \param[in] s A string containing an "option=value" pair that needs to be split up.
74 * \param[out] opt The name of the option.
75 * \param[out] val Value of the option.
77 static void splitOptionValue(const string &s, string &opt, string &val)
79 size_t eqPos = s.find('=');
80 if (eqPos != string::npos)
82 opt = s.substr(0, eqPos);
83 if (eqPos != s.length()) val = s.substr(eqPos+1);
87 /*!
88 * \brief Compare two strings ignoring case.
89 * This function is in fact a wrapper around the gromacs function gmx_strncasecmp().
90 * \param[in] s1 String.
91 * \param[in] s2 String.
92 * \returns Similarly to the C function strncasecmp(), the return value is an
93 integer less than, equal to, or greater than 0 if \p s1 less than,
94 identical to, or greater than \p s2.
96 static bool isStringEqNCase(const string s1, const string s2)
98 return (gmx_strncasecmp(s1.c_str(), s2.c_str(), max(s1.length(), s2.length())) == 0);
102 * \brief Convert string to upper case.
104 * \param[in] s String to convert to uppercase.
105 * \returns The given string converted to uppercase.
107 static string toUpper(const string &s)
109 string stmp(s);
110 std::transform(stmp.begin(), stmp.end(), stmp.begin(), static_cast < int(*)(int) > (toupper));
111 return stmp;
114 /*!
115 \name Sizes of constant device option arrays GmxOpenMMPlatformOptions#platforms,
116 GmxOpenMMPlatformOptions#memtests, GmxOpenMMPlatformOptions#deviceid,
117 GmxOpenMMPlatformOptions#force_dev. */
118 /* {@ */
119 #define SIZEOF_PLATFORMS 1
120 #define SIZEOF_MEMTESTS 3
121 #define SIZEOF_DEVICEIDS 1
122 #define SIZEOF_FORCE_DEV 2
123 /* @} */
125 /*! Possible platform options in the mdrun -device option. */
126 static const char *devOptStrings[] = { "platform", "deviceid", "memtest", "force-device" };
128 /*! Enumerated platform options in the mdrun -device option. */
129 enum devOpt
131 PLATFORM = 0,
132 DEVICEID = 1,
133 MEMTEST = 2,
134 FORCE_DEVICE = 3
138 * \brief Class to extract and manage the platform options in the mdrun -device option.
141 class GmxOpenMMPlatformOptions
143 public:
144 GmxOpenMMPlatformOptions(const char *opt);
145 ~GmxOpenMMPlatformOptions() { options.clear(); }
146 string getOptionValue(const string &opt);
147 void remOption(const string &opt);
148 private:
149 void setOption(const string &opt, const string &val);
151 map<string, string> options; /*!< Data structure to store the option (name, value) pairs. */
153 static const char * const platforms[SIZEOF_PLATFORMS]; /*!< Available OpenMM platforms; size #SIZEOF_PLATFORMS */
154 static const char * const memtests[SIZEOF_MEMTESTS]; /*!< Available types of memory tests, also valid
155 any positive integer >=15; size #SIZEOF_MEMTESTS */
156 static const char * const deviceid[SIZEOF_DEVICEIDS]; /*!< Possible values for deviceid option;
157 also valid any positive integer; size #SIZEOF_DEVICEIDS */
158 static const char * const force_dev[SIZEOF_FORCE_DEV]; /*!< Possible values for for force-device option;
159 size #SIZEOF_FORCE_DEV */
162 const char * const GmxOpenMMPlatformOptions::platforms[SIZEOF_PLATFORMS] = { "Cuda" /*,"OpenCL"*/ };
163 const char * const GmxOpenMMPlatformOptions::memtests[SIZEOF_MEMTESTS] = { "15", "full", "off" };
164 const char * const GmxOpenMMPlatformOptions::deviceid[SIZEOF_DEVICEIDS] = { "0" };
165 const char * const GmxOpenMMPlatformOptions::force_dev[SIZEOF_FORCE_DEV] = { "no", "yes" };
168 * \brief Contructor.
169 * Takes the option list, parses it, checks the options and their values for validity.
170 * When certain options are not provided by the user, as default value the first item
171 * of the respective constant array is taken (GmxOpenMMPlatformOptions#platforms,
172 * GmxOpenMMPlatformOptions#memtests, GmxOpenMMPlatformOptions#deviceid,
173 * GmxOpenMMPlatformOptions#force_dev).
174 * \param[in] optionString Option string part of the mdrun -deviceoption parameter.
176 GmxOpenMMPlatformOptions::GmxOpenMMPlatformOptions(const char *optionString)
178 // set default values
179 setOption("platform", platforms[0]);
180 setOption("memtest", memtests[0]);
181 setOption("deviceid", deviceid[0]);
182 setOption("force-device", force_dev[0]);
184 string opt(optionString);
186 // remove all whitespaces
187 opt.erase(remove_if(opt.begin(), opt.end(), ::isspace), opt.end());
188 // tokenize around ","-s
189 vector<string> tokens = split(opt, ',');
191 for (vector<string>::iterator it = tokens.begin(); it != tokens.end(); ++it)
193 string opt = "", val = "";
194 splitOptionValue(*it, opt, val);
196 if (isStringEqNCase(opt, "platform"))
198 setOption(opt, val);
199 continue;
202 if (isStringEqNCase(opt, "memtest"))
204 if (!isStringEqNCase(val, "full") && !isStringEqNCase(val, "off")) /* the value has to be an integer >15 */
206 int secs;
207 if (!from_string<int>(secs, val, std::dec))
209 gmx_fatal(FARGS, "Invalid value for option memtestoption: \"%s\"!\n", val.c_str());
211 if (secs < 15)
213 gmx_fatal(FARGS, "Incorrect value for memtest option (%d). "
214 "Memtest needs to run for at least 15s!\n", secs);
217 setOption(opt, val);
218 continue;
221 if (isStringEqNCase(opt, "deviceid"))
223 int id;
224 if (!from_string<int>(id, val, std::dec) )
226 gmx_fatal(FARGS, "Invalid device id: \"%s\"!\n", val.c_str());
228 setOption(opt, val);
229 continue;
232 if (isStringEqNCase(opt, "force-device"))
234 if (!isStringEqNCase(val, "yes") && !isStringEqNCase(val, "no"))
236 gmx_fatal(FARGS, "Invalid OpenMM force option: \"%s\"!\n", val.c_str());
238 setOption(opt, val);
239 continue;
242 // if we got till here something went wrong
243 gmx_fatal(FARGS, "Invalid OpenMM platform option: \"%s\"!\n", (*it).c_str());
249 * \brief Returns the value of an option.
250 * \param[in] opt Name of the option.
252 string GmxOpenMMPlatformOptions::getOptionValue(const string &opt)
254 if (options.find(toUpper(opt)) != options.end())
256 return options[toUpper(opt)];
258 else
260 return NULL;
265 * \brief Setter function - private, only used from contructor.
266 * \param[in] opt Name of the option.
267 * \param[in] val Value for the option.
269 void GmxOpenMMPlatformOptions::setOption(const string &opt, const string &val)
271 options[toUpper(opt)] = val;
275 * \brief Removes an option with its value from the map structure. If the option
276 * does not exist, returns without any action.
277 * \param[in] opt Name of the option.
279 void GmxOpenMMPlatformOptions::remOption(const string &opt)
281 options.erase(toUpper(opt));
286 * \brief Container for OpenMM related data structures that represent the bridge
287 * between the Gromacs data-structures and the OpenMM library and is but it's
288 * only passed through the API functions as void to disable direct access.
290 class OpenMMData
292 public:
293 System* system; /*! The system to simulate. */
294 Context* context; /*! The OpenMM context in which the simulation is carried out. */
295 Integrator* integrator; /*! The integrator used in the simulation. */
296 bool removeCM; /*! If \true remove venter of motion, false otherwise. */
297 GmxOpenMMPlatformOptions *platformOpt; /*! Platform options. */
301 * \brief Runs memtest on the GPU that has alreaby been initialized by OpenMM.
302 * \param[in] fplog Pointer to gromacs log file.
303 * \param[in] devId Device id of the GPU to run the test on. TODO: this can be removed!
304 * \param[in] pre_post Contains either "Pre" or "Post" just to be able to differentiate in
305 * stdout messages/log between memtest carried out before and after simulation.
306 * \param[in] opt Pointer to platform options object.
308 void runMemtest(FILE* fplog, int devId, const char* pre_post, GmxOpenMMPlatformOptions *opt)
311 char warn_buf[STRLEN];
313 int which_test;
314 int res;
315 const char * test_type = opt->getOptionValue("memtest").c_str();
316 if (!gmx_strcasecmp(test_type, "off"))
318 which_test = 0;
320 else
322 if (!gmx_strcasecmp(test_type, "full"))
324 which_test = 2;
326 else
328 from_string<int>(which_test, test_type, std::dec);
332 switch (which_test)
334 case 0: /* no memtest */
335 sprintf(warn_buf, "%s-simulation GPU memtest skipped. Note, that faulty memory can cause "
336 "incorrect results!\n", pre_post);
337 fprintf(fplog, "%s", warn_buf);
338 gmx_warning(warn_buf);
339 break; /* case 0 */
341 case 1: /* quick memtest */
342 fprintf(fplog, "%s-simulation %s GPU memtest in progress...\n", pre_post, test_type);
343 fprintf(stdout, "\n%s-simulation %s GPU memtest in progress...", pre_post, test_type);
344 fflush(fplog);
345 fflush(stdout);
346 if (do_quick_memtest(-1) != 0)
348 gmx_fatal(FARGS,MEM_ERR_MSG(pre_post));
350 else
352 fprintf(fplog, "Memory test completed without errors.\n");
353 fprintf(stdout, "done, no errors detected\n");
355 break; /* case 1 */
357 case 2: /* full memtest */
358 fprintf(fplog, "%s-simulation %s memtest in progress...\n", pre_post, test_type);
359 fprintf(stdout, "\n%s-simulation %s memtest in progress...", pre_post, test_type);
360 fflush(fplog);
361 fflush(stdout);
362 if (do_full_memtest(-1) != 0)
364 gmx_fatal(FARGS, MEM_ERR_MSG(pre_post) );
367 else
369 fprintf(fplog, "Memory test completed without errors.\n");
370 fprintf(stdout, "done, no errors detected\n");
372 break; /* case 2 */
374 default: /* timed memtest */
375 fprintf(fplog, "%s-simulation memtest for ~%ds in progress...\n", pre_post, which_test);
376 fprintf(stdout, "\n%s-simulation memtest for ~%ds in progress...", pre_post, which_test);
377 fflush(fplog);
378 fflush(stdout);
379 if (do_timed_memtest(-1, which_test) != 0)
381 gmx_fatal(FARGS, MEM_ERR_MSG(pre_post));
384 else
386 fprintf(fplog, "Memory test completed without errors.\n");
387 fprintf(stdout, "done, no errors detected.\n");
390 fflush(fplog);
391 fflush(stdout);
395 * \brief Does gromacs option checking.
397 * Checks the gromacs mdp options for features unsupported in OpenMM, case in which
398 * interrupts the execution. It also warns the user about pecularities of OpenMM
399 * implementations.
400 * \param[in] ir Gromacs structure for input options, \see ::t_inputrec
401 * \param[in] top Gromacs topology, \see ::gmx_localtop_t
403 void checkGmxOptions(t_inputrec *ir, gmx_localtop_t *top)
406 char warn_buf[STRLEN];
408 // Abort if unsupported critical options are present
410 /* Integrator */
411 if (ir->eI == eiMD)
412 gmx_warning( "OpenMM does not support leap-frog, will use velocity-verlet integrator.\n");
414 if ( (ir->eI != eiMD) &&
415 (ir->eI != eiVV) &&
416 (ir->eI != eiVVAK) &&
417 (ir->eI != eiSD1) &&
418 (ir->eI != eiSD2) &&
419 (ir->eI != eiBD) )
421 gmx_fatal(FARGS, "OpenMM supports only the following integrators: md/md-vv/md-vv-avek, sd/sd1, and bd.\n");
424 /* Electroctstics */
425 if ( (ir->coulombtype != eelPME) &&
426 (ir->coulombtype != eelRF) &&
427 (ir->coulombtype != eelEWALD) &&
428 // no-cutoff
429 ( !(ir->coulombtype == eelCUT && ir->rcoulomb == 0 && ir->rvdw == 0)) )
431 gmx_fatal(FARGS,"OpenMM supports only the following methods for electrostatics: "
432 "NoCutoff (i.e. rcoulomb = rvdw = 0 ),Reaction-Field, Ewald or PME.\n");
435 if ( (ir->etc != etcNO) &&
436 (ir->eI != eiSD1) &&
437 (ir->eI != eiSD2) &&
438 (ir->eI != eiBD) )
440 gmx_warning("OpenMM supports only Andersen thermostat with the md/md-vv/md-vv-avek integrators.\n");
443 if (ir->opts.ngtc > 1)
444 gmx_fatal(FARGS,"OpenMM does not support multiple temperature coupling groups.\n");
446 if (ir->epc != etcNO)
447 gmx_fatal(FARGS,"OpenMM does not support pressure coupling.\n");
449 if (ir->opts.annealing[0])
450 gmx_fatal(FARGS,"OpenMM does not support simulated annealing.\n");
452 if (ir->eConstrAlg != econtSHAKE)
453 gmx_warning("Constraints in OpenMM are done by a combination "
454 "of SHAKE, SETTLE and CCMA. Accuracy is based on the SHAKE tolerance set "
455 "by the \"shake_tol\" option.\n");
457 if (ir->nwall != 0)
458 gmx_fatal(FARGS,"OpenMM does not support walls.\n");
460 if (ir->ePull != epullNO)
461 gmx_fatal(FARGS,"OpenMM does not support pulling.\n");
463 if (top->idef.il[F_DISRES].nr > 0)
464 gmx_fatal(FARGS,"OpenMM does not support distant restraints.\n");
466 if (top->idef.il[F_ORIRES].nr > 0)
467 gmx_fatal(FARGS,"OpenMM does not support orientation restraints.\n");
469 if (top->idef.il[F_ANGRES].nr > 0)
470 gmx_fatal(FARGS,"OpenMM does not support angle restraints.\n");
472 if (top->idef.il[F_DIHRES].nr > 0)
473 gmx_fatal(FARGS,"OpenMM does not support dihedral restraints.\n");
475 if (ir->efep != efepNO)
476 gmx_fatal(FARGS,"OpenMM does not support free energy calculations.\n");
478 if (ir->opts.ngacc > 1)
479 gmx_fatal(FARGS,"OpenMM does not support non-equilibrium MD (accelerated groups).\n");
481 if (IR_ELEC_FIELD(*ir))
482 gmx_fatal(FARGS,"OpenMM does not support electric fields.\n");
484 if (ir->bQMMM)
485 gmx_fatal(FARGS,"OpenMM does not support QMMM calculations.\n");
487 if (ir->rcoulomb != ir->rvdw)
488 gmx_fatal(FARGS,"OpenMM uses a single cutoff for both Coulomb "
489 "and VdW interactions. Please set rcoulomb equal to rvdw.\n");
494 * \brief Initialize OpenMM, run sanity/consistency checks, and return a pointer to
495 * the OpenMMData.
497 * Various gromacs data structures are passed that contain the parameters, state and
498 * other porperties of the system to simulate. These serve as input for initializing
499 * OpenMM. Besides, a set of misc action are taken:
500 * - OpenMM plugins are loaded;
501 * - platform options in \p platformOptStr are parsed and checked;
502 * - Gromacs parameters are checked for OpenMM support and consistency;
503 * - after the OpenMM is initialized memtest executed in the same GPU context.
505 * \param[in] fplog Gromacs log file handler.
506 * \param[in] platformOptStr Platform option string.
507 * \param[in] cr TODO remove!
508 * \param[in] ir The Gromacs input parameters.
509 * \param[in] top_global TODO
510 * \param[in] top TODO
511 * \param[in] mdatoms TODO
512 * \param[in] fr TODO
513 * \param[in] state TODO
516 void* openmm_init(FILE *fplog, const char *platformOptStr,
517 t_commrec *cr,t_inputrec *ir,
518 gmx_mtop_t *top_global, gmx_localtop_t *top,
519 t_mdatoms *mdatoms, t_forcerec *fr,t_state *state)
522 char warn_buf[STRLEN];
523 static bool hasLoadedPlugins = false;
524 string usedPluginDir;
525 int devId;
529 if (!hasLoadedPlugins)
531 vector<string> loadedPlugins;
532 /* Look for OpenMM plugins at various locations (listed in order of priority):
533 - on the path in OPENMM_PLUGIN_DIR environment variable if this is specified
534 - on the path in the OPENMM_PLUGIN_DIR macro that is set by the build script
535 - at the default location assumed by OpenMM
537 /* env var */
538 char *pluginDir = getenv("OPENMM_PLUGIN_DIR");
539 trim(pluginDir);
540 /* no env var or empty */
541 if (pluginDir != NULL && *pluginDir != '\0')
543 loadedPlugins = Platform::loadPluginsFromDirectory(pluginDir);
544 if (loadedPlugins.size() > 0)
546 hasLoadedPlugins = true;
547 usedPluginDir = pluginDir;
549 else
551 gmx_fatal(FARGS, "The directory provided in the OPENMM_PLUGIN_DIR environment variable "
552 "(%s) does not contain valid OpenMM plugins. Check your OpenMM installation!",
553 pluginDir);
557 /* macro set at build time */
558 #ifdef OPENMM_PLUGIN_DIR
559 if (!hasLoadedPlugins)
561 loadedPlugins = Platform::loadPluginsFromDirectory(OPENMM_PLUGIN_DIR);
562 if (loadedPlugins.size() > 0)
564 hasLoadedPlugins = true;
565 usedPluginDir = OPENMM_PLUGIN_DIR;
568 #endif
569 /* default loocation */
570 if (!hasLoadedPlugins)
572 loadedPlugins = Platform::loadPluginsFromDirectory(Platform::getDefaultPluginsDirectory());
573 if (loadedPlugins.size() > 0)
575 hasLoadedPlugins = true;
576 usedPluginDir = Platform::getDefaultPluginsDirectory();
580 /* if there are still no plugins loaded there won't be any */
581 if (!hasLoadedPlugins)
583 gmx_fatal(FARGS, "No OpenMM plugins were found! You can provide the"
584 " plugin directory in the OPENMM_PLUGIN_DIR environment variable.", pluginDir);
587 fprintf(fplog, "\nPlugins loaded from directory %s:\t", usedPluginDir.c_str());
588 for (int i = 0; i < loadedPlugins.size(); i++)
590 fprintf(fplog, "%s, ", loadedPlugins[i].c_str());
592 fprintf(fplog, "\n");
595 /* parse option string */
596 GmxOpenMMPlatformOptions *opt = new GmxOpenMMPlatformOptions(platformOptStr);
598 /* check wheter Gromacs options compatibility with OpenMM */
599 checkGmxOptions(ir, top);
601 // Create the system.
602 const t_idef& idef = top->idef;
603 const int numAtoms = top_global->natoms;
604 const int numConstraints = idef.il[F_CONSTR].nr/3;
605 const int numSettle = idef.il[F_SETTLE].nr/2;
606 const int numBonds = idef.il[F_BONDS].nr/3;
607 const int numAngles = idef.il[F_ANGLES].nr/4;
608 const int numPeriodic = idef.il[F_PDIHS].nr/5;
609 const int numRB = idef.il[F_RBDIHS].nr/5;
610 const int num14 = idef.il[F_LJ14].nr/3;
611 System* sys = new System();
612 if (ir->nstcomm > 0)
613 sys->addForce(new CMMotionRemover(ir->nstcomm));
615 // Set bonded force field terms.
616 const int* bondAtoms = (int*) idef.il[F_BONDS].iatoms;
617 HarmonicBondForce* bondForce = new HarmonicBondForce();
618 sys->addForce(bondForce);
619 int offset = 0;
620 for (int i = 0; i < numBonds; ++i)
622 int type = bondAtoms[offset++];
623 int atom1 = bondAtoms[offset++];
624 int atom2 = bondAtoms[offset++];
625 bondForce->addBond(atom1, atom2,
626 idef.iparams[type].harmonic.rA, idef.iparams[type].harmonic.krA);
628 const int* angleAtoms = (int*) idef.il[F_ANGLES].iatoms;
629 HarmonicAngleForce* angleForce = new HarmonicAngleForce();
630 sys->addForce(angleForce);
631 offset = 0;
632 for (int i = 0; i < numAngles; ++i)
634 int type = angleAtoms[offset++];
635 int atom1 = angleAtoms[offset++];
636 int atom2 = angleAtoms[offset++];
637 int atom3 = angleAtoms[offset++];
638 angleForce->addAngle(atom1, atom2, atom3,
639 idef.iparams[type].harmonic.rA*M_PI/180.0, idef.iparams[type].harmonic.krA);
641 const int* periodicAtoms = (int*) idef.il[F_PDIHS].iatoms;
642 PeriodicTorsionForce* periodicForce = new PeriodicTorsionForce();
643 sys->addForce(periodicForce);
644 offset = 0;
645 for (int i = 0; i < numPeriodic; ++i)
647 int type = periodicAtoms[offset++];
648 int atom1 = periodicAtoms[offset++];
649 int atom2 = periodicAtoms[offset++];
650 int atom3 = periodicAtoms[offset++];
651 int atom4 = periodicAtoms[offset++];
652 periodicForce->addTorsion(atom1, atom2, atom3, atom4,
653 idef.iparams[type].pdihs.mult,
654 idef.iparams[type].pdihs.phiA*M_PI/180.0,
655 idef.iparams[type].pdihs.cpA);
657 const int* rbAtoms = (int*) idef.il[F_RBDIHS].iatoms;
658 RBTorsionForce* rbForce = new RBTorsionForce();
659 sys->addForce(rbForce);
660 offset = 0;
661 for (int i = 0; i < numRB; ++i)
663 int type = rbAtoms[offset++];
664 int atom1 = rbAtoms[offset++];
665 int atom2 = rbAtoms[offset++];
666 int atom3 = rbAtoms[offset++];
667 int atom4 = rbAtoms[offset++];
668 rbForce->addTorsion(atom1, atom2, atom3, atom4,
669 idef.iparams[type].rbdihs.rbcA[0], idef.iparams[type].rbdihs.rbcA[1],
670 idef.iparams[type].rbdihs.rbcA[2], idef.iparams[type].rbdihs.rbcA[3],
671 idef.iparams[type].rbdihs.rbcA[4], idef.iparams[type].rbdihs.rbcA[5]);
674 // Set nonbonded parameters and masses.
675 int ntypes = fr->ntype;
676 int* types = mdatoms->typeA;
677 real* nbfp = fr->nbfp;
678 real* charges = mdatoms->chargeA;
679 real* masses = mdatoms->massT;
680 NonbondedForce* nonbondedForce = new NonbondedForce();
681 sys->addForce(nonbondedForce);
683 if (ir->rcoulomb == 0)
685 nonbondedForce->setNonbondedMethod(NonbondedForce::NoCutoff);
687 else
689 switch (ir->coulombtype)
691 case eelRF: // TODO what is the correct condition?
692 if (ir->ePBC == epbcXYZ)
694 nonbondedForce->setNonbondedMethod(NonbondedForce::CutoffPeriodic);
695 sys->setPeriodicBoxVectors(Vec3(state->box[0][0], 0, 0),
696 Vec3(0, state->box[1][1], 0), Vec3(0, 0, state->box[2][2]));
698 else if (ir->ePBC == epbcNONE)
699 nonbondedForce->setNonbondedMethod(NonbondedForce::CutoffNonPeriodic);
700 else
701 gmx_fatal(FARGS,"OpenMM supports only full periodic boundary conditions "
702 "(pbc = xyz), or none (pbc = no).\n");
703 nonbondedForce->setCutoffDistance(ir->rcoulomb);
704 break;
706 case eelEWALD:
707 if (ir->ewald_geometry == eewg3DC)
708 gmx_fatal(FARGS,"OpenMM supports only Ewald 3D geometry.\n");
709 if (ir->epsilon_surface != 0)
710 gmx_fatal(FARGS,"OpenMM does not support dipole correction in Ewald summation.\n");
711 nonbondedForce->setNonbondedMethod(NonbondedForce::Ewald);
712 nonbondedForce->setCutoffDistance(ir->rcoulomb);
713 sys->setPeriodicBoxVectors(Vec3(state->box[0][0], 0, 0),
714 Vec3(0, state->box[1][1], 0), Vec3(0, 0, state->box[2][2]));
715 break;
717 case eelPME:
718 nonbondedForce->setNonbondedMethod(NonbondedForce::PME);
719 nonbondedForce->setCutoffDistance(ir->rcoulomb);
720 sys->setPeriodicBoxVectors(Vec3(state->box[0][0], 0, 0),
721 Vec3(0, state->box[1][1], 0), Vec3(0, 0, state->box[2][2]));
722 break;
724 default:
725 gmx_fatal(FARGS,"Internal error: you should not see this message, it that the"
726 "electrosatics option check failed. Please report this error!");
730 for (int i = 0; i < numAtoms; ++i)
732 real c6 = nbfp[types[i]*2*ntypes+types[i]*2];
733 real c12 = nbfp[types[i]*2*ntypes+types[i]*2+1];
734 if (c12 <= 0)
735 nonbondedForce->addParticle(charges[i], 1.0, 0.0);
736 else
737 nonbondedForce->addParticle(charges[i], pow(c12/c6, (real) (1.0/6.0)), c6*c6/(4.0*c12));
738 sys->addParticle(masses[i]);
741 // Build a table of all exclusions.
742 vector<set<int> > exclusions(numAtoms);
743 for (int i = 0; i < numAtoms; i++)
745 int start = top->excls.index[i];
746 int end = top->excls.index[i+1];
747 for (int j = start; j < end; j++)
748 exclusions[i].insert(top->excls.a[j]);
751 // Record the 1-4 interactions, and remove them from the list of exclusions.
752 const int* nb14Atoms = (int*) idef.il[F_LJ14].iatoms;
753 offset = 0;
754 for (int i = 0; i < num14; ++i)
756 int type = nb14Atoms[offset++];
757 int atom1 = nb14Atoms[offset++];
758 int atom2 = nb14Atoms[offset++];
759 real c6 = idef.iparams[type].lj14.c6A;
760 real c12 = idef.iparams[type].lj14.c12A;
761 real sigma, epsilon;
762 if (c12 <= 0)
764 epsilon = (real) 0.0;
765 sigma = (real) 1.0;
767 else
769 epsilon = (real) ((c6*c6)/(4.0*c12));
770 sigma = (real) pow(c12/c6, (real) (1.0/6.0));
772 nonbondedForce->addException(atom1, atom2,
773 fr->fudgeQQ*charges[atom1]*charges[atom2], sigma, epsilon);
774 exclusions[atom1].erase(atom2);
775 exclusions[atom2].erase(atom1);
778 // Record exclusions.
779 for (int i = 0; i < numAtoms; i++)
781 for (set<int>::const_iterator iter = exclusions[i].begin(); iter != exclusions[i].end(); ++iter)
783 if (i < *iter)
784 nonbondedForce->addException(i, *iter, 0.0, 1.0, 0.0);
788 // Add GBSA if needed.
789 t_atoms atoms;
790 atoms = gmx_mtop_global_atoms(top_global);
792 if (ir->implicit_solvent == eisGBSA)
794 GBSAOBCForce* gbsa = new GBSAOBCForce();
795 sys->addForce(gbsa);
796 gbsa->setSoluteDielectric(ir->epsilon_r);
797 gbsa->setSolventDielectric(ir->gb_epsilon_solvent);
798 gbsa->setCutoffDistance(nonbondedForce->getCutoffDistance());
799 if (nonbondedForce->getNonbondedMethod() == NonbondedForce::NoCutoff)
800 gbsa->setNonbondedMethod(GBSAOBCForce::NoCutoff);
801 else if (nonbondedForce->getNonbondedMethod() == NonbondedForce::CutoffNonPeriodic)
802 gbsa->setNonbondedMethod(GBSAOBCForce::CutoffNonPeriodic);
803 else if (nonbondedForce->getNonbondedMethod() == NonbondedForce::CutoffPeriodic)
804 gbsa->setNonbondedMethod(GBSAOBCForce::CutoffPeriodic);
805 else
806 gmx_fatal(FARGS,"OpenMM supports only Reaction-Field electrostatics with OBC/GBSA.\n");
808 for (int i = 0; i < numAtoms; ++i)
809 gbsa->addParticle(charges[i],
810 top_global->atomtypes.gb_radius[atoms.atom[i].type],
811 top_global->atomtypes.S_hct[atoms.atom[i].type]);
814 // Set constraints.
815 const int* constraintAtoms = (int*) idef.il[F_CONSTR].iatoms;
816 offset = 0;
817 for (int i = 0; i < numConstraints; ++i)
819 int type = constraintAtoms[offset++];
820 int atom1 = constraintAtoms[offset++];
821 int atom2 = constraintAtoms[offset++];
822 sys->addConstraint(atom1, atom2, idef.iparams[type].constr.dA);
824 const int* settleAtoms = (int*) idef.il[F_SETTLE].iatoms;
825 offset = 0;
826 for (int i = 0; i < numSettle; ++i)
828 int type = settleAtoms[offset++];
829 int oxygen = settleAtoms[offset++];
830 sys->addConstraint(oxygen, oxygen+1, idef.iparams[type].settle.doh);
831 sys->addConstraint(oxygen, oxygen+2, idef.iparams[type].settle.doh);
832 sys->addConstraint(oxygen+1, oxygen+2, idef.iparams[type].settle.dhh);
835 // Create an integrator for simulating the system.
836 real friction = (ir->opts.tau_t[0] == 0.0 ? 0.0 : 1.0/ir->opts.tau_t[0]);
837 Integrator* integ;
838 if (ir->eI == eiMD || ir->eI == eiVV || ir->eI == eiVVAK)
840 integ = new VerletIntegrator(ir->delta_t);
841 if ( ir->etc != etcNO)
843 real collisionFreq = ir->opts.tau_t[0] / 1000; /* tau_t (ps) / 1000 = collisionFreq (fs^-1) */
844 AndersenThermostat* thermostat = new AndersenThermostat(ir->opts.ref_t[0], friction);
845 sys->addForce(thermostat);
848 else if (ir->eI == eiBD)
850 integ = new BrownianIntegrator(ir->opts.ref_t[0], friction, ir->delta_t);
851 static_cast<BrownianIntegrator*>(integ)->setRandomNumberSeed(ir->ld_seed);
853 else if (EI_SD(ir->eI))
855 integ = new LangevinIntegrator(ir->opts.ref_t[0], friction, ir->delta_t);
856 static_cast<LangevinIntegrator*>(integ)->setRandomNumberSeed(ir->ld_seed);
858 integ->setConstraintTolerance(ir->shake_tol);
860 // Create a context and initialize it.
861 Context* context = NULL;
864 OpenMM could automatically select the "best" GPU, however we're not't
865 going to let it do that for now, as the current algorithm is very rudimentary
866 and we anyway support only CUDA.
867 if (platformOptStr == NULL || platformOptStr == "")
869 context = new Context(*sys, *integ);
871 else
874 // Find which platform is it.
875 for (int i = 0; i < Platform::getNumPlatforms() && context == NULL; i++)
877 if (isStringEqNCase(opt->getOptionValue("platform"), Platform::getPlatform(i).getName()))
879 Platform& platform = Platform::getPlatform(i);
880 // set standard properties
881 platform.setPropertyDefaultValue("CudaDevice", opt->getOptionValue("deviceid"));
882 // TODO add extra properties
883 context = new Context(*sys, *integ, platform);
886 if (context == NULL)
888 gmx_fatal(FARGS, "The requested platform \"%s\" could not be found.\n",
889 opt->getOptionValue("platform").c_str());
893 Platform& platform = context->getPlatform();
894 fprintf(fplog, "Gromacs will use the OpenMM platform: %s\n", platform.getName().c_str());
896 const vector<string>& properties = platform.getPropertyNames();
897 if (debug)
899 for (int i = 0; i < properties.size(); i++)
901 printf(">> %s: %s\n", properties[i].c_str(),
902 platform.getPropertyValue(*context, properties[i]).c_str());
903 fprintf(fplog, ">> %s: %s\n", properties[i].c_str(),
904 platform.getPropertyValue(*context, properties[i]).c_str());
908 /* For now this is just to double-check if OpenMM selected the GPU we wanted,
909 but when we'll let OpenMM select the GPU automatically it will query the devideId.
911 int tmp;
912 if (!from_string<int>(tmp, platform.getPropertyValue(*context, "CudaDevice"), std::dec))
914 gmx_fatal(FARGS, "Internal error: couldn't determine the device selected by OpenMM");
915 if (tmp != devId)
917 gmx_fatal(FARGS, "Internal error: OpenMM is using device #%d while initialized for device #%d",
918 tmp, devId);
922 /* check GPU compatibility */
923 char gpuname[STRLEN];
924 devId = atoi(opt->getOptionValue("deviceid").c_str());
925 if (!is_supported_cuda_gpu(-1, gpuname))
927 if (!gmx_strcasecmp(opt->getOptionValue("force-device").c_str(), "yes"))
929 sprintf(warn_buf, "Non-supported GPU selected (#%d, %s), forced continuing.\n"
930 "Note, that the simulation can be slow or it migth even crash.",
931 devId, gpuname);
932 fprintf(fplog, "%s", warn_buf);
933 gmx_warning(warn_buf);
935 else
937 gmx_fatal(FARGS, "The selected GPU (#%d, %s) is not supported by Gromacs! "
938 "Most probably you have a low-end GPU which would not perform well, "
939 "or new hardware that has not been tested yet with Gromacs-OpenMM. "
940 "If you still want to try using the device, use the force-device=yes option.",
941 devId, gpuname);
944 else
946 fprintf(fplog, "Gromacs will run on the GPU #%d (%s).\n", devId, gpuname);
949 /* do the pre-simulation memtest */
950 runMemtest(fplog, -1, "Pre", opt);
952 vector<Vec3> pos(numAtoms);
953 vector<Vec3> vel(numAtoms);
954 for (int i = 0; i < numAtoms; ++i)
956 pos[i] = Vec3(state->x[i][0], state->x[i][1], state->x[i][2]);
957 vel[i] = Vec3(state->v[i][0], state->v[i][1], state->v[i][2]);
959 context->setPositions(pos);
960 context->setVelocities(vel);
962 // Return a structure containing the system, integrator, and context.
963 OpenMMData* data = new OpenMMData();
964 data->system = sys;
965 data->integrator = integ;
966 data->context = context;
967 data->removeCM = (ir->nstcomm > 0);
968 data->platformOpt = opt;
969 return data;
972 catch (std::exception& e)
974 gmx_fatal(FARGS, "OpenMM exception caught while initializating: %s\n", e.what());
979 * \brief Integrate one step.
981 * \param[in] data OpenMMData object created by openmm_init().
983 void openmm_take_one_step(void* data)
985 // static int step = 0; printf("----> taking step #%d\n", step++);
988 static_cast<OpenMMData*>(data)->integrator->step(1);
990 catch (std::exception& e)
992 gmx_fatal(FARGS, "OpenMM exception caught while taking a step: %s\n", e.what());
997 * \brief Integrate n steps.
999 * \param[in] data OpenMMData object created by openmm_init().
1001 void openmm_take_steps(void* data, int nstep)
1005 static_cast<OpenMMData*>(data)->integrator->step(nstep);
1007 catch (std::exception& e)
1009 gmx_fatal(FARGS, "OpenMM exception caught while taking a step: %s\n", e.what());
1014 * \brief Clean up the data structures cretead for OpenMM.
1016 * \param[in] log Log file pointer.
1017 * \param[in] data OpenMMData object created by openmm_init().
1019 void openmm_cleanup(FILE* fplog, void* data)
1021 OpenMMData* d = static_cast<OpenMMData*>(data);
1022 runMemtest(fplog, -1, "Post", d->platformOpt);
1023 delete d->system;
1024 delete d->integrator;
1025 delete d->context;
1026 delete d->platformOpt;
1027 delete d;
1031 * \brief Copy the current state information from OpenMM into the Gromacs data structures.
1033 * This function results in the requested proprties to be copied from the
1034 * GPU to host. As this represents a bottleneck, the frequency of pulling data
1035 * should be minimized.
1037 * \param[in] data OpenMMData object created by openmm_init().
1038 * \param[out] time
1039 * \param[out] state State of the system: coordinates and velocities.
1040 * \param[out] f Forces.
1041 * \param[out] enerd Energies.
1042 * \param[in] includePos True if coordinates are requested.
1043 * \param[in] includeVel True if velocities are requested.
1044 * \param[in] includeForce True if forces are requested.
1045 * \param[in] includeEnergy True if energies are requested.
1047 void openmm_copy_state(void *data,
1048 t_state *state, double *time,
1049 rvec f[], gmx_enerdata_t *enerd,
1050 bool includePos, bool includeVel, bool includeForce, bool includeEnergy)
1052 int types = 0;
1053 if (includePos)
1054 types += State::Positions;
1055 if (includeVel)
1056 types += State::Velocities;
1057 if (includeForce)
1058 types += State::Forces;
1059 if (includeEnergy)
1060 types += State::Energy;
1061 if (types == 0)
1062 return;
1065 State currentState = static_cast<OpenMMData*>(data)->context->getState(types);
1066 int numAtoms = static_cast<OpenMMData*>(data)->system->getNumParticles();
1067 if (includePos)
1069 for (int i = 0; i < numAtoms; i++)
1071 Vec3 x = currentState.getPositions()[i];
1072 state->x[i][0] = x[0];
1073 state->x[i][1] = x[1];
1074 state->x[i][2] = x[2];
1077 if (includeVel)
1079 for (int i = 0; i < numAtoms; i++)
1081 Vec3 v = currentState.getVelocities()[i];
1082 state->v[i][0] = v[0];
1083 state->v[i][1] = v[1];
1084 state->v[i][2] = v[2];
1087 if (includeForce)
1089 for (int i = 0; i < numAtoms; i++)
1091 Vec3 force = currentState.getForces()[i];
1092 f[i][0] = force[0];
1093 f[i][1] = force[1];
1094 f[i][2] = force[2];
1097 if (includeEnergy)
1099 int numConstraints = static_cast<OpenMMData*>(data)->system->getNumConstraints();
1100 int dof = 3*numAtoms-numConstraints;
1101 if (static_cast<OpenMMData*>(data)->removeCM)
1102 dof -= 3;
1103 enerd->term[F_EPOT] = currentState.getPotentialEnergy();
1104 enerd->term[F_EKIN] = currentState.getKineticEnergy();
1105 enerd->term[F_ETOT] = enerd->term[F_EPOT] + enerd->term[F_EKIN];
1106 enerd->term[F_TEMP] = 2.0*enerd->term[F_EKIN]/dof/BOLTZ;
1108 *time = currentState.getTime();
1110 catch (std::exception& e)
1112 gmx_fatal(FARGS, "OpenMM exception caught while retrieving state information: %s\n", e.what());