1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2010, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
37 * Note, that parts of this source code originate from the Simtk release
38 * of OpenMM accelerated Gromacs, for more details see:
39 * https://simtk.org/project/xml/downloads.xml?group_id=161#package_id600
56 #include "gmx_fatal.h"
61 #include "gmx_gpu_utils.h"
62 #include "mtop_util.h"
64 #include "openmm_wrapper.h"
66 using namespace OpenMM
;
69 #define MEM_ERR_MSG(str) \
70 "The %s-simulation GPU memory test detected errors. As memory errors would cause incorrect " \
71 "simulation results, gromacs has aborted execution.\n Make sure that your GPU's memory is not " \
72 "overclocked and that the device is properly cooled.\n", (str)
76 * \brief Convert string to integer type.
77 * \param[in] s String to convert from.
78 * \param[in] f Basefield format flag that takes any of the following I/O
79 * manipulators: dec, hex, oct.
80 * \param[out] t Destination variable to convert to.
83 static bool from_string(T
& t
, const string
& s
, ios_base
& (*f
)(ios_base
&))
86 return !(iss
>> f
>> t
).fail();
90 * \brief Split string around a given delimiter.
91 * \param[in] s String to split.
92 * \param[in] delim Delimiter character.
93 * \returns Vector of strings found in \p s.
95 static vector
<string
> split(const string
&s
, char delim
)
100 while (getline(ss
, item
, delim
))
102 if (item
.length() != 0)
103 elems
.push_back(item
);
109 * \brief Split a string of the form "option=value" into "option" and "value" strings.
110 * This string corresponds to one option and the associated value from the option list
111 * in the mdrun -device argument.
113 * \param[in] s A string containing an "option=value" pair that needs to be split up.
114 * \param[out] opt The name of the option.
115 * \param[out] val Value of the option.
117 static void splitOptionValue(const string
&s
, string
&opt
, string
&val
)
119 size_t eqPos
= s
.find('=');
120 if (eqPos
!= string::npos
)
122 opt
= s
.substr(0, eqPos
);
123 if (eqPos
!= s
.length()) val
= s
.substr(eqPos
+1);
128 * \brief Compare two strings ignoring case.
129 * This function is in fact a wrapper around the gromacs function gmx_strncasecmp().
130 * \param[in] s1 String.
131 * \param[in] s2 String.
132 * \returns Similarly to the C function strncasecmp(), the return value is an
133 integer less than, equal to, or greater than 0 if \p s1 less than,
134 identical to, or greater than \p s2.
136 static bool isStringEqNCase(const string s1
, const string s2
)
138 return (gmx_strncasecmp(s1
.c_str(), s2
.c_str(), max(s1
.length(), s2
.length())) == 0);
142 * \brief Convert string to upper case.
144 * \param[in] s String to convert to uppercase.
145 * \returns The given string converted to uppercase.
147 static string
toUpper(const string
&s
)
150 std::transform(stmp
.begin(), stmp
.end(), stmp
.begin(), static_cast < int(*)(int) > (toupper
));
155 \name Sizes of constant device option arrays GmxOpenMMPlatformOptions#platforms,
156 GmxOpenMMPlatformOptions#memtests, GmxOpenMMPlatformOptions#deviceid,
157 GmxOpenMMPlatformOptions#force_dev. */
159 #define SIZEOF_PLATFORMS 1 // 2
160 #define SIZEOF_MEMTESTS 3
161 #define SIZEOF_DEVICEIDS 1
162 #define SIZEOF_FORCE_DEV 2
165 /*! Possible platform options in the mdrun -device option. */
166 static const char *devOptStrings
[] = { "platform", "deviceid", "memtest", "force-device" };
168 /*! Enumerated platform options in the mdrun -device option. */
178 * \brief Class to extract and manage the platform options in the mdrun -device option.
181 class GmxOpenMMPlatformOptions
184 GmxOpenMMPlatformOptions(const char *opt
);
185 ~GmxOpenMMPlatformOptions() { options
.clear(); }
186 string
getOptionValue(const string
&opt
);
187 void remOption(const string
&opt
);
190 void setOption(const string
&opt
, const string
&val
);
192 map
<string
, string
> options
; /*!< Data structure to store the option (name, value) pairs. */
194 static const char * const platforms
[SIZEOF_PLATFORMS
]; /*!< Available OpenMM platforms; size #SIZEOF_PLATFORMS */
195 static const char * const memtests
[SIZEOF_MEMTESTS
]; /*!< Available types of memory tests, also valid
196 any positive integer >=15; size #SIZEOF_MEMTESTS */
197 static const char * const deviceid
[SIZEOF_DEVICEIDS
]; /*!< Possible values for deviceid option;
198 also valid any positive integer; size #SIZEOF_DEVICEIDS */
199 static const char * const force_dev
[SIZEOF_FORCE_DEV
]; /*!< Possible values for for force-device option;
200 size #SIZEOF_FORCE_DEV */
203 const char * const GmxOpenMMPlatformOptions::platforms
[SIZEOF_PLATFORMS
]
205 //= { "Reference", "CUDA" /*,"OpenCL"*/ };
206 const char * const GmxOpenMMPlatformOptions::memtests
[SIZEOF_MEMTESTS
]
207 = { "15", "full", "off" };
208 const char * const GmxOpenMMPlatformOptions::deviceid
[SIZEOF_DEVICEIDS
]
210 const char * const GmxOpenMMPlatformOptions::force_dev
[SIZEOF_FORCE_DEV
]
215 * Takes the option list, parses it, checks the options and their values for validity.
216 * When certain options are not provided by the user, as default value the first item
217 * of the respective constant array is taken (GmxOpenMMPlatformOptions#platforms,
218 * GmxOpenMMPlatformOptions#memtests, GmxOpenMMPlatformOptions#deviceid,
219 * GmxOpenMMPlatformOptions#force_dev).
220 * \param[in] optionString Option list part of the mdrun -device parameter.
222 GmxOpenMMPlatformOptions::GmxOpenMMPlatformOptions(const char *optionString
)
224 // set default values
225 setOption("platform", platforms
[0]);
226 setOption("memtest", memtests
[0]);
227 setOption("deviceid", deviceid
[0]);
228 setOption("force-device", force_dev
[0]);
230 string
opt(optionString
);
232 // remove all whitespaces
233 opt
.erase(remove_if(opt
.begin(), opt
.end(), ::isspace
), opt
.end());
234 // tokenize around ","-s
235 vector
<string
> tokens
= split(opt
, ',');
237 for (vector
<string
>::iterator it
= tokens
.begin(); it
!= tokens
.end(); ++it
)
239 string opt
= "", val
= "";
240 splitOptionValue(*it
, opt
, val
);
242 if (isStringEqNCase(opt
, "platform"))
244 /* no check, this will fail if platform does not exist when we try to set it */
249 if (isStringEqNCase(opt
, "memtest"))
251 /* the value has to be an integer >15(s) or "full" OR "off" */
252 if (!isStringEqNCase(val
, "full") && !isStringEqNCase(val
, "off"))
255 if (!from_string
<int>(secs
, val
, std::dec
))
257 gmx_fatal(FARGS
, "Invalid value for option memtest option: \"%s\"!", val
.c_str());
261 gmx_fatal(FARGS
, "Incorrect value for memtest option (%d). "
262 "Memtest needs to run for at least 15s!", secs
);
269 if (isStringEqNCase(opt
, "deviceid"))
272 if (!from_string
<int>(id
, val
, std::dec
) )
274 gmx_fatal(FARGS
, "Invalid device id: \"%s\"!", val
.c_str());
280 if (isStringEqNCase(opt
, "force-device"))
283 if (!isStringEqNCase(val
, "yes") && !isStringEqNCase(val
, "no"))
285 gmx_fatal(FARGS
, "Invalid OpenMM force option: \"%s\"!", val
.c_str());
291 // if we got till here something went wrong
292 gmx_fatal(FARGS
, "Invalid OpenMM platform option: \"%s\"!", (*it
).c_str());
298 * \brief Getter function.
299 * \param[in] opt Name of the option.
300 * \returns Returns the value associated to an option.
302 string
GmxOpenMMPlatformOptions::getOptionValue(const string
&opt
)
304 map
<string
, string
> :: const_iterator it
= options
.find(toUpper(opt
));
305 if (it
!= options
.end())
307 // cout << "@@@>> " << it->first << " : " << it->second << endl;
317 * \brief Setter function - private, only used from contructor.
318 * \param[in] opt Name of the option.
319 * \param[in] val Value for the option.
321 void GmxOpenMMPlatformOptions::setOption(const string
&opt
, const string
&val
)
323 options
[toUpper(opt
)] = val
;
327 * \brief Removes an option with its value from the map structure. If the option
328 * does not exist, returns without any action.
329 * \param[in] opt Name of the option.
331 void GmxOpenMMPlatformOptions::remOption(const string
&opt
)
333 options
.erase(toUpper(opt
));
337 * \brief Print option-value pairs to a file (debugging function).
339 void GmxOpenMMPlatformOptions::print()
341 cout
<< ">> Platform options: " << endl
342 << ">> platform = " << getOptionValue("platform") << endl
343 << ">> deviceID = " << getOptionValue("deviceid") << endl
344 << ">> memtest = " << getOptionValue("memtest") << endl
345 << ">> force-device = " << getOptionValue("force-device") << endl
;
349 * \brief Container for OpenMM related data structures that represent the bridge
350 * between the Gromacs data-structures and the OpenMM library and is but it's
351 * only passed through the API functions as void to disable direct access.
356 System
* system
; /*! The system to simulate. */
357 Context
* context
; /*! The OpenMM context in which the simulation is carried out. */
358 Integrator
* integrator
; /*! The integrator used in the simulation. */
359 bool removeCM
; /*! If \true remove venter of motion, false otherwise. */
360 GmxOpenMMPlatformOptions
*platformOpt
; /*! Platform options. */
364 * \brief Runs memtest on the GPU that has alreaby been initialized by OpenMM.
365 * \param[in] fplog Pointer to gromacs log file.
366 * \param[in] devId Device id of the GPU to run the test on.
367 Note: as OpenMM previously creates the context,for now this is always -1.
368 * \param[in] pre_post Contains either "Pre" or "Post" just to be able to differentiate in
369 * stdout messages/log between memtest carried out before and after simulation.
370 * \param[in] opt Pointer to platform options object.
372 static void runMemtest(FILE* fplog
, int devId
, const char* pre_post
, GmxOpenMMPlatformOptions
*opt
)
374 char strout_buf
[STRLEN
];
377 string s
= opt
->getOptionValue("memtest");
378 const char *test_type
= s
.c_str();
380 if (!gmx_strcasecmp(test_type
, "off"))
386 if (!gmx_strcasecmp(test_type
, "full"))
392 from_string
<int>(which_test
, test_type
, std::dec
);
398 gmx_fatal(FARGS
, "Amount of seconds for memetest is negative (%d). ", which_test
);
403 case 0: /* no memtest */
404 sprintf(strout_buf
, "%s-simulation GPU memtest skipped. Note, that faulty memory can cause "
405 "incorrect results!", pre_post
);
406 fprintf(fplog
, "%s\n", strout_buf
);
407 gmx_warning(strout_buf
);
410 case 1: /* quick memtest */
411 fprintf(fplog
, "%s-simulation %s GPU memtest in progress...\n", pre_post
, test_type
);
412 fprintf(stdout
, "\n%s-simulation %s GPU memtest in progress...", pre_post
, test_type
);
415 res
= do_quick_memtest(devId
);
418 case 2: /* full memtest */
419 fprintf(fplog
, "%s-simulation %s memtest in progress...\n", pre_post
, test_type
);
420 fprintf(stdout
, "\n%s-simulation %s memtest in progress...", pre_post
, test_type
);
423 res
= do_full_memtest(devId
);
426 default: /* timed memtest */
427 fprintf(fplog
, "%s-simulation ~%ds memtest in progress...\n", pre_post
, which_test
);
428 fprintf(stdout
, "\n%s-simulation ~%ds memtest in progress...", pre_post
, which_test
);
431 res
= do_timed_memtest(devId
, which_test
);
438 gmx_fatal(FARGS
, MEM_ERR_MSG(pre_post
));
442 fprintf(fplog
, "Memory test completed without errors.\n");
444 fprintf(stdout
, "done, no errors detected\n");
451 * \brief Convert Lennard-Jones parameters c12 and c6 to sigma and epsilon.
456 * \param[out] epsilon
458 static void convert_c_12_6(double c12
, double c6
, double *sigma
, double *epsilon
)
460 if (c12
== 0 && c6
== 0)
465 else if (c12
> 0 && c6
> 0)
467 *epsilon
= (c6
*c6
)/(4.0*c12
);
468 *sigma
= pow(c12
/c6
, 1.0/6.0);
472 gmx_fatal(FARGS
,"OpenMM only supports c6 > 0 and c12 > 0 or c6 = c12 = 0.");
477 * \brief Does gromacs option checking.
479 * Checks the gromacs mdp options for features unsupported in OpenMM, case in which
480 * interrupts the execution. It also warns the user about pecularities of OpenMM
482 * \param[in] fplog Gromacs log file pointer.
483 * \param[in] ir Gromacs input parameters, see ::t_inputrec
484 * \param[in] top Gromacs node local topology, \see gmx_localtop_t
485 * \param[in] state Gromacs state structure \see ::t_state
486 * \param[in] mdatoms Gromacs atom parameters, \see ::t_mdatoms
487 * \param[in] fr \see ::t_forcerec
488 * \param[in] state Gromacs systems state, \see ::t_state
490 static void checkGmxOptions(FILE* fplog
,
491 t_inputrec
*ir
, gmx_localtop_t
*top
,
492 t_forcerec
*fr
, t_state
*state
)
494 char warn_buf
[STRLEN
];
496 double c6
, c12
, sigma_ij
, sigma_ji
, sigma_ii
, sigma_jj
, sigma_comb
,
497 eps_ij
, eps_ji
, eps_ii
, eps_jj
, eps_comb
;
499 /* Abort if unsupported critical options are present */
504 gmx_warning( "OpenMM does not support leap-frog, will use velocity-verlet integrator.");
507 if ( (ir
->eI
!= eiMD
) &&
509 (ir
->eI
!= eiVVAK
) &&
514 gmx_fatal(FARGS
, "OpenMM supports only the following integrators: md/md-vv/md-vv-avek, sd/sd1, and bd.");
518 if ( (ir
->coulombtype
!= eelPME
) &&
519 (ir
->coulombtype
!= eelRF
) &&
520 (ir
->coulombtype
!= eelEWALD
) &&
522 ( !(ir
->coulombtype
== eelCUT
&& ir
->rcoulomb
== 0 && ir
->rvdw
== 0)) )
524 gmx_fatal(FARGS
,"OpenMM supports only the following methods for electrostatics: "
525 "NoCutoff (i.e. rcoulomb = rvdw = 0 ),Reaction-Field, Ewald or PME.");
528 if (ir
->etc
!= etcNO
&&
533 gmx_warning("OpenMM supports only Andersen thermostat with the md/md-vv/md-vv-avek integrators.");
536 if (ir
->opts
.ngtc
> 1)
537 gmx_fatal(FARGS
,"OpenMM does not support multiple temperature coupling groups.");
539 if (ir
->epc
!= etcNO
)
540 gmx_fatal(FARGS
,"OpenMM does not support pressure coupling.");
542 if (ir
->opts
.annealing
[0])
543 gmx_fatal(FARGS
,"OpenMM does not support simulated annealing.");
545 if (top
->idef
.il
[F_CONSTR
].nr
> 0 && ir
->eConstrAlg
!= econtSHAKE
)
546 gmx_warning("OpenMM provides contraints as a combination "
547 "of SHAKE, SETTLE and CCMA. Accuracy is based on the SHAKE tolerance set "
548 "by the \"shake_tol\" option.");
551 gmx_fatal(FARGS
,"OpenMM does not support walls.");
553 if (ir
->ePull
!= epullNO
)
554 gmx_fatal(FARGS
,"OpenMM does not support pulling.");
556 /* check for restraints */
557 for (i
= 0; i
< F_EPOT
; i
++)
566 top
->idef
.il
[i
].nr
> 0)
568 /* TODO fix the message */
569 gmx_fatal(FARGS
, "OpenMM does not support (some) of the provided restaint(s).");
573 if (ir
->efep
!= efepNO
)
574 gmx_fatal(FARGS
,"OpenMM does not support free energy calculations.");
576 if (ir
->opts
.ngacc
> 1)
577 gmx_fatal(FARGS
,"OpenMM does not support non-equilibrium MD (accelerated groups).");
579 if (IR_ELEC_FIELD(*ir
))
580 gmx_fatal(FARGS
,"OpenMM does not support electric fields.");
583 gmx_fatal(FARGS
,"OpenMM does not support QMMM calculations.");
585 if (ir
->rcoulomb
!= ir
->rvdw
)
586 gmx_fatal(FARGS
,"OpenMM uses a single cutoff for both Coulomb "
587 "and VdW interactions. Please set rcoulomb equal to rvdw.");
589 if (EEL_FULL(ir
->coulombtype
))
591 if (ir
->ewald_geometry
== eewg3DC
)
592 gmx_fatal(FARGS
,"OpenMM supports only Ewald 3D geometry.");
593 if (ir
->epsilon_surface
!= 0)
594 gmx_fatal(FARGS
,"OpenMM does not support dipole correction in Ewald summation.");
597 if (TRICLINIC(state
->box
))
599 gmx_fatal(FARGS
,"OpenMM does not support triclinic unit cells.");
605 * \brief Initialize OpenMM, run sanity/consistency checks, and return a pointer to
608 * Various gromacs data structures are passed that contain the parameters, state and
609 * other porperties of the system to simulate. These serve as input for initializing
610 * OpenMM. Besides, a set of misc action are taken:
611 * - OpenMM plugins are loaded;
612 * - platform options in \p platformOptStr are parsed and checked;
613 * - Gromacs parameters are checked for OpenMM support and consistency;
614 * - after the OpenMM is initialized memtest executed in the same GPU context.
616 * \param[in] fplog Gromacs log file handler.
617 * \param[in] platformOptStr Platform option string.
618 * \param[in] ir The Gromacs input parameters, see ::t_inputrec
619 * \param[in] top_global Gromacs system toppology, \see ::gmx_mtop_t
620 * \param[in] top Gromacs node local topology, \see gmx_localtop_t
621 * \param[in] mdatoms Gromacs atom parameters, \see ::t_mdatoms
622 * \param[in] fr \see ::t_forcerec
623 * \param[in] state Gromacs systems state, \see ::t_state
624 * \returns Pointer to a
627 void* openmm_init(FILE *fplog
, const char *platformOptStr
,
629 gmx_mtop_t
*top_global
, gmx_localtop_t
*top
,
630 t_mdatoms
*mdatoms
, t_forcerec
*fr
, t_state
*state
)
633 char warn_buf
[STRLEN
];
634 static bool hasLoadedPlugins
= false;
635 string usedPluginDir
;
640 if (!hasLoadedPlugins
)
642 vector
<string
> loadedPlugins
;
643 /* Look for OpenMM plugins at various locations (listed in order of priority):
644 - on the path in OPENMM_PLUGIN_DIR environment variable if this is specified
645 - on the path in the OPENMM_PLUGIN_DIR macro that is set by the build script
646 - at the default location assumed by OpenMM
649 char *pluginDir
= getenv("OPENMM_PLUGIN_DIR");
651 /* no env var or empty */
652 if (pluginDir
!= NULL
&& *pluginDir
!= '\0')
654 loadedPlugins
= Platform::loadPluginsFromDirectory(pluginDir
);
655 if (loadedPlugins
.size() > 0)
657 hasLoadedPlugins
= true;
658 usedPluginDir
= pluginDir
;
662 gmx_fatal(FARGS
, "The directory provided in the OPENMM_PLUGIN_DIR environment variable "
663 "(%s) does not contain valid OpenMM plugins. Check your OpenMM installation!",
668 /* macro set at build time */
669 #ifdef OPENMM_PLUGIN_DIR
670 if (!hasLoadedPlugins
)
672 loadedPlugins
= Platform::loadPluginsFromDirectory(OPENMM_PLUGIN_DIR
);
673 if (loadedPlugins
.size() > 0)
675 hasLoadedPlugins
= true;
676 usedPluginDir
= OPENMM_PLUGIN_DIR
;
680 /* default loocation */
681 if (!hasLoadedPlugins
)
683 loadedPlugins
= Platform::loadPluginsFromDirectory(Platform::getDefaultPluginsDirectory());
684 if (loadedPlugins
.size() > 0)
686 hasLoadedPlugins
= true;
687 usedPluginDir
= Platform::getDefaultPluginsDirectory();
691 /* if there are still no plugins loaded there won't be any */
692 if (!hasLoadedPlugins
)
694 gmx_fatal(FARGS
, "No OpenMM plugins were found! You can provide the"
695 " plugin directory in the OPENMM_PLUGIN_DIR environment variable.", pluginDir
);
698 fprintf(fplog
, "\nPlugins loaded from directory %s:\t", usedPluginDir
.c_str());
699 for (int i
= 0; i
< (int)loadedPlugins
.size(); i
++)
701 fprintf(fplog
, "%s, ", loadedPlugins
[i
].c_str());
703 fprintf(fplog
, "\n");
706 /* parse option string */
707 GmxOpenMMPlatformOptions
*opt
= new GmxOpenMMPlatformOptions(platformOptStr
);
708 devId
= atoi(opt
->getOptionValue("deviceid").c_str());
715 /* check wheter Gromacs options compatibility with OpenMM */
716 checkGmxOptions(fplog
, ir
, top
, fr
, state
);
718 // Create the system.
719 const t_idef
& idef
= top
->idef
;
720 const int numAtoms
= top_global
->natoms
;
721 const int numConstraints
= idef
.il
[F_CONSTR
].nr
/3;
722 const int numSettle
= idef
.il
[F_SETTLE
].nr
/2;
723 const int numBonds
= idef
.il
[F_BONDS
].nr
/3;
724 const int numUB
= idef
.il
[F_UREY_BRADLEY
].nr
/4;
725 const int numAngles
= idef
.il
[F_ANGLES
].nr
/4;
726 const int numPeriodic
= idef
.il
[F_PDIHS
].nr
/5;
727 const int numRB
= idef
.il
[F_RBDIHS
].nr
/5;
728 const int num14
= idef
.il
[F_LJ14
].nr
/3;
729 System
* sys
= new System();
731 sys
->addForce(new CMMotionRemover(ir
->nstcomm
));
733 // Set bonded force field terms.
734 const int* bondAtoms
= (int*) idef
.il
[F_BONDS
].iatoms
;
735 HarmonicBondForce
* bondForce
= new HarmonicBondForce();
736 sys
->addForce(bondForce
);
738 for (int i
= 0; i
< numBonds
; ++i
)
740 int type
= bondAtoms
[offset
++];
741 int atom1
= bondAtoms
[offset
++];
742 int atom2
= bondAtoms
[offset
++];
743 bondForce
->addBond(atom1
, atom2
,
744 idef
.iparams
[type
].harmonic
.rA
, idef
.iparams
[type
].harmonic
.krA
);
746 // Urey-Bradley includes both the angle and bond potential for 1-3 interactions
747 const int* ubAtoms
= (int*) idef
.il
[F_UREY_BRADLEY
].iatoms
;
748 HarmonicBondForce
* ubBondForce
= new HarmonicBondForce();
749 HarmonicAngleForce
* ubAngleForce
= new HarmonicAngleForce();
750 sys
->addForce(ubBondForce
);
751 sys
->addForce(ubAngleForce
);
753 for (int i
= 0; i
< numUB
; ++i
)
755 int type
= ubAtoms
[offset
++];
756 int atom1
= ubAtoms
[offset
++];
757 int atom2
= ubAtoms
[offset
++];
758 int atom3
= ubAtoms
[offset
++];
759 ubBondForce
->addBond(atom1
, atom3
,
760 idef
.iparams
[type
].u_b
.r13
, idef
.iparams
[type
].u_b
.kUB
);
761 ubAngleForce
->addAngle(atom1
, atom2
, atom3
,
762 idef
.iparams
[type
].u_b
.theta
*M_PI
/180.0, idef
.iparams
[type
].u_b
.ktheta
);
764 const int* angleAtoms
= (int*) idef
.il
[F_ANGLES
].iatoms
;
765 HarmonicAngleForce
* angleForce
= new HarmonicAngleForce();
766 sys
->addForce(angleForce
);
768 for (int i
= 0; i
< numAngles
; ++i
)
770 int type
= angleAtoms
[offset
++];
771 int atom1
= angleAtoms
[offset
++];
772 int atom2
= angleAtoms
[offset
++];
773 int atom3
= angleAtoms
[offset
++];
774 angleForce
->addAngle(atom1
, atom2
, atom3
,
775 idef
.iparams
[type
].harmonic
.rA
*M_PI
/180.0, idef
.iparams
[type
].harmonic
.krA
);
777 const int* periodicAtoms
= (int*) idef
.il
[F_PDIHS
].iatoms
;
778 PeriodicTorsionForce
* periodicForce
= new PeriodicTorsionForce();
779 sys
->addForce(periodicForce
);
781 for (int i
= 0; i
< numPeriodic
; ++i
)
783 int type
= periodicAtoms
[offset
++];
784 int atom1
= periodicAtoms
[offset
++];
785 int atom2
= periodicAtoms
[offset
++];
786 int atom3
= periodicAtoms
[offset
++];
787 int atom4
= periodicAtoms
[offset
++];
788 periodicForce
->addTorsion(atom1
, atom2
, atom3
, atom4
,
789 idef
.iparams
[type
].pdihs
.mult
,
790 idef
.iparams
[type
].pdihs
.phiA
*M_PI
/180.0,
791 idef
.iparams
[type
].pdihs
.cpA
);
793 const int* rbAtoms
= (int*) idef
.il
[F_RBDIHS
].iatoms
;
794 RBTorsionForce
* rbForce
= new RBTorsionForce();
795 sys
->addForce(rbForce
);
797 for (int i
= 0; i
< numRB
; ++i
)
799 int type
= rbAtoms
[offset
++];
800 int atom1
= rbAtoms
[offset
++];
801 int atom2
= rbAtoms
[offset
++];
802 int atom3
= rbAtoms
[offset
++];
803 int atom4
= rbAtoms
[offset
++];
804 rbForce
->addTorsion(atom1
, atom2
, atom3
, atom4
,
805 idef
.iparams
[type
].rbdihs
.rbcA
[0], idef
.iparams
[type
].rbdihs
.rbcA
[1],
806 idef
.iparams
[type
].rbdihs
.rbcA
[2], idef
.iparams
[type
].rbdihs
.rbcA
[3],
807 idef
.iparams
[type
].rbdihs
.rbcA
[4], idef
.iparams
[type
].rbdihs
.rbcA
[5]);
810 // Set nonbonded parameters and masses.
811 int ntypes
= fr
->ntype
;
812 int* types
= mdatoms
->typeA
;
813 real
* nbfp
= fr
->nbfp
;
814 real
* charges
= mdatoms
->chargeA
;
815 real
* masses
= mdatoms
->massT
;
816 NonbondedForce
* nonbondedForce
= new NonbondedForce();
817 sys
->addForce(nonbondedForce
);
822 if (ir
->rcoulomb
== 0)
824 nonbondedForce
->setNonbondedMethod(NonbondedForce::NoCutoff
);
828 nonbondedForce
->setNonbondedMethod(NonbondedForce::CutoffNonPeriodic
);
832 switch (ir
->coulombtype
)
835 nonbondedForce
->setNonbondedMethod(NonbondedForce::CutoffPeriodic
);
839 nonbondedForce
->setNonbondedMethod(NonbondedForce::Ewald
);
843 nonbondedForce
->setNonbondedMethod(NonbondedForce::PME
);
847 gmx_fatal(FARGS
,"Internal error: you should not see this message, it that the"
848 "electrosatics option check failed. Please report this error!");
850 sys
->setPeriodicBoxVectors(Vec3(state
->box
[0][0], 0, 0),
851 Vec3(0, state
->box
[1][1], 0), Vec3(0, 0, state
->box
[2][2]));
852 nonbondedForce
->setCutoffDistance(ir
->rcoulomb
);
856 gmx_fatal(FARGS
,"OpenMM supports only full periodic boundary conditions "
857 "(pbc = xyz), or none (pbc = no).");
861 /* Fix for PME and Ewald error tolerance
863 * OpenMM uses approximate formulas to calculate the Ewald parameter:
864 * alpha = (1.0/cutoff)*sqrt(-log(2.0*tolerlance));
865 * and the grid spacing for PME:
866 * gridX = ceil(alpha*box[0][0]/pow(0.5*tol, 0.2));
867 * gridY = ceil(alpha*box[1][1]/pow(0.5*tol, 0.2));
868 * gridZ = ceil(alpha*box[2][2]/pow(0.5*tol, 0.2));
870 * It overestimates the precision and setting it to
871 * (500 x ewald_rtol) seems to give a reasonable match to the GROMACS settings
873 * If the default ewald_rtol=1e-5 is used we silently adjust the value,
874 * otherwise a warning is issued about the action taken.
876 double corr_ewald_rtol
= 500.0 * ir
->ewald_rtol
;
877 if ((ir
->ePBC
== epbcXYZ
) &&
878 (ir
->coulombtype
== eelEWALD
|| ir
->coulombtype
== eelPME
))
882 fprintf(debug
, ">> ewald_rtol = %e (corrected = %e) \n",
883 ir
->ewald_rtol
, corr_ewald_rtol
);
886 if (fabs(ir
->ewald_rtol
- 1e-5) > 1e-10)
888 gmx_warning("OpenMM uses the ewald_rtol parameter with approximate formulas "
889 "to calculate the alpha and grid spacing parameters of the Ewald "
890 "and PME methods. This tolerance need to be corrected in order to get "
891 "settings close to the ones used in GROMACS. Although the internal correction "
892 "should work for any reasonable value of ewald_rtol, using values other than "
893 "the default 1e-5 might cause incorrect behavior.");
895 if (corr_ewald_rtol
> 1)
897 gmx_fatal(FARGS
, "The ewald_rtol accuracy term is >1 after the "
898 "adjustment for OpenMM (%e)", corr_ewald_rtol
);
901 nonbondedForce
->setEwaldErrorTolerance(corr_ewald_rtol
);
904 for (int i
= 0; i
< numAtoms
; ++i
)
906 double c12
= nbfp
[types
[i
]*2*ntypes
+types
[i
]*2+1];
907 double c6
= nbfp
[types
[i
]*2*ntypes
+types
[i
]*2];
908 double sigma
, epsilon
;
909 convert_c_12_6(c12
, c6
, &sigma
, &epsilon
);
910 nonbondedForce
->addParticle(charges
[i
], sigma
, epsilon
);
911 sys
->addParticle(masses
[i
]);
914 // Build a table of all exclusions.
915 vector
<set
<int> > exclusions(numAtoms
);
916 for (int i
= 0; i
< numAtoms
; i
++)
918 int start
= top
->excls
.index
[i
];
919 int end
= top
->excls
.index
[i
+1];
920 for (int j
= start
; j
< end
; j
++)
921 exclusions
[i
].insert(top
->excls
.a
[j
]);
924 // Record the 1-4 interactions, and remove them from the list of exclusions.
925 const int* nb14Atoms
= (int*) idef
.il
[F_LJ14
].iatoms
;
927 for (int i
= 0; i
< num14
; ++i
)
929 int type
= nb14Atoms
[offset
++];
930 int atom1
= nb14Atoms
[offset
++];
931 int atom2
= nb14Atoms
[offset
++];
932 double sigma
, epsilon
;
933 convert_c_12_6(idef
.iparams
[type
].lj14
.c12A
,
934 idef
.iparams
[type
].lj14
.c6A
,
936 nonbondedForce
->addException(atom1
, atom2
,
937 fr
->fudgeQQ
*charges
[atom1
]*charges
[atom2
], sigma
, epsilon
);
938 exclusions
[atom1
].erase(atom2
);
939 exclusions
[atom2
].erase(atom1
);
942 // Record exclusions.
943 for (int i
= 0; i
< numAtoms
; i
++)
945 for (set
<int>::const_iterator iter
= exclusions
[i
].begin(); iter
!= exclusions
[i
].end(); ++iter
)
949 nonbondedForce
->addException(i
, *iter
, 0.0, 1.0, 0.0);
954 // Add GBSA if needed.
955 if (ir
->implicit_solvent
== eisGBSA
)
957 t_atoms atoms
= gmx_mtop_global_atoms(top_global
);
958 GBSAOBCForce
* gbsa
= new GBSAOBCForce();
961 gbsa
->setSoluteDielectric(ir
->epsilon_r
);
962 gbsa
->setSolventDielectric(ir
->gb_epsilon_solvent
);
963 gbsa
->setCutoffDistance(nonbondedForce
->getCutoffDistance());
964 if (nonbondedForce
->getNonbondedMethod() == NonbondedForce::NoCutoff
)
965 gbsa
->setNonbondedMethod(GBSAOBCForce::NoCutoff
);
966 else if (nonbondedForce
->getNonbondedMethod() == NonbondedForce::CutoffNonPeriodic
)
967 gbsa
->setNonbondedMethod(GBSAOBCForce::CutoffNonPeriodic
);
968 else if (nonbondedForce
->getNonbondedMethod() == NonbondedForce::CutoffPeriodic
)
969 gbsa
->setNonbondedMethod(GBSAOBCForce::CutoffPeriodic
);
971 gmx_fatal(FARGS
,"OpenMM supports only Reaction-Field electrostatics with OBC/GBSA.");
973 for (int i
= 0; i
< numAtoms
; ++i
)
975 gbsa
->addParticle(charges
[i
],
976 top_global
->atomtypes
.gb_radius
[atoms
.atom
[i
].type
],
977 top_global
->atomtypes
.S_hct
[atoms
.atom
[i
].type
]);
979 free_t_atoms(&atoms
, FALSE
);
983 const int* constraintAtoms
= (int*) idef
.il
[F_CONSTR
].iatoms
;
985 for (int i
= 0; i
< numConstraints
; ++i
)
987 int type
= constraintAtoms
[offset
++];
988 int atom1
= constraintAtoms
[offset
++];
989 int atom2
= constraintAtoms
[offset
++];
990 sys
->addConstraint(atom1
, atom2
, idef
.iparams
[type
].constr
.dA
);
992 const int* settleAtoms
= (int*) idef
.il
[F_SETTLE
].iatoms
;
994 for (int i
= 0; i
< numSettle
; ++i
)
996 int type
= settleAtoms
[offset
++];
997 int oxygen
= settleAtoms
[offset
++];
998 sys
->addConstraint(oxygen
, oxygen
+1, idef
.iparams
[type
].settle
.doh
);
999 sys
->addConstraint(oxygen
, oxygen
+2, idef
.iparams
[type
].settle
.doh
);
1000 sys
->addConstraint(oxygen
+1, oxygen
+2, idef
.iparams
[type
].settle
.dhh
);
1003 // Create an integrator for simulating the system.
1004 double friction
= (ir
->opts
.tau_t
[0] == 0.0 ? 0.0 : 1.0/ir
->opts
.tau_t
[0]);
1008 integ
= new BrownianIntegrator(ir
->opts
.ref_t
[0], friction
, ir
->delta_t
);
1009 static_cast<BrownianIntegrator
*>(integ
)->setRandomNumberSeed(ir
->ld_seed
);
1011 else if (EI_SD(ir
->eI
))
1013 integ
= new LangevinIntegrator(ir
->opts
.ref_t
[0], friction
, ir
->delta_t
);
1014 static_cast<LangevinIntegrator
*>(integ
)->setRandomNumberSeed(ir
->ld_seed
);
1018 integ
= new VerletIntegrator(ir
->delta_t
);
1019 if ( ir
->etc
!= etcNO
)
1021 real collisionFreq
= ir
->opts
.tau_t
[0] / 1000; /* tau_t (ps) / 1000 = collisionFreq (fs^-1) */
1022 AndersenThermostat
* thermostat
= new AndersenThermostat(ir
->opts
.ref_t
[0], friction
);
1023 sys
->addForce(thermostat
);
1026 integ
->setConstraintTolerance(ir
->shake_tol
);
1028 // Create a context and initialize it.
1029 Context
* context
= NULL
;
1032 OpenMM could automatically select the "best" GPU, however we're not't
1033 going to let it do that for now, as the current algorithm is very rudimentary
1034 and we anyway support only CUDA.
1035 if (platformOptStr == NULL || platformOptStr == "")
1037 context = new Context(*sys, *integ);
1042 /* which platform should we use */
1043 for (int i
= 0; i
< (int)Platform::getNumPlatforms() && context
== NULL
; i
++)
1045 if (isStringEqNCase(opt
->getOptionValue("platform"), Platform::getPlatform(i
).getName()))
1047 Platform
& platform
= Platform::getPlatform(i
);
1048 // set standard properties
1049 platform
.setPropertyDefaultValue("CudaDevice", opt
->getOptionValue("deviceid"));
1050 // TODO add extra properties
1051 context
= new Context(*sys
, *integ
, platform
);
1054 if (context
== NULL
)
1056 gmx_fatal(FARGS
, "The requested platform \"%s\" could not be found.",
1057 opt
->getOptionValue("platform").c_str());
1061 Platform
& platform
= context
->getPlatform();
1062 fprintf(fplog
, "Gromacs will use the OpenMM platform: %s\n", platform
.getName().c_str());
1064 const vector
<string
>& properties
= platform
.getPropertyNames();
1067 for (int i
= 0; i
< (int)properties
.size(); i
++)
1069 fprintf(debug
, ">> %s: %s\n", properties
[i
].c_str(),
1070 platform
.getPropertyValue(*context
, properties
[i
]).c_str());
1075 if (isStringEqNCase(opt
->getOptionValue("platform"), "CUDA"))
1077 /* For now this is just to double-check if OpenMM selected the GPU we wanted,
1078 but when we'll let OpenMM select the GPU automatically, it will query the devideId.
1081 if (!from_string
<int>(tmp
, platform
.getPropertyValue(*context
, "CudaDevice"), std::dec
))
1083 gmx_fatal(FARGS
, "Internal error: couldn't determine the device selected by OpenMM");
1086 gmx_fatal(FARGS
, "Internal error: OpenMM is using device #%d"
1087 "while initialized for device #%d", tmp
, devId
);
1091 /* check GPU compatibility */
1092 char gpuname
[STRLEN
];
1093 devId
= atoi(opt
->getOptionValue("deviceid").c_str());
1094 if (!is_supported_cuda_gpu(-1, gpuname
))
1096 if (!gmx_strcasecmp(opt
->getOptionValue("force-device").c_str(), "yes"))
1098 sprintf(warn_buf
, "Non-supported GPU selected (#%d, %s), forced continuing."
1099 "Note, that the simulation can be slow or it migth even crash.",
1101 fprintf(fplog
, "%s\n", warn_buf
);
1102 gmx_warning(warn_buf
);
1106 gmx_fatal(FARGS
, "The selected GPU (#%d, %s) is not supported by Gromacs! "
1107 "Most probably you have a low-end GPU which would not perform well, "
1108 "or new hardware that has not been tested with the current release. "
1109 "If you still want to try using the device, use the force-device=yes option.",
1115 fprintf(fplog
, "Gromacs will run on the GPU #%d (%s).\n", devId
, gpuname
);
1120 if (isStringEqNCase(opt
->getOptionValue("platform"), "CUDA"))
1122 /* pre-simulation memtest */
1123 runMemtest(fplog
, -1, "Pre", opt
);
1126 vector
<Vec3
> pos(numAtoms
);
1127 vector
<Vec3
> vel(numAtoms
);
1128 for (int i
= 0; i
< numAtoms
; ++i
)
1130 pos
[i
] = Vec3(state
->x
[i
][0], state
->x
[i
][1], state
->x
[i
][2]);
1131 vel
[i
] = Vec3(state
->v
[i
][0], state
->v
[i
][1], state
->v
[i
][2]);
1133 context
->setPositions(pos
);
1134 context
->setVelocities(vel
);
1136 // Return a structure containing the system, integrator, and context.
1137 OpenMMData
* data
= new OpenMMData();
1139 data
->integrator
= integ
;
1140 data
->context
= context
;
1141 data
->removeCM
= (ir
->nstcomm
> 0);
1142 data
->platformOpt
= opt
;
1145 catch (std::exception
& e
)
1147 gmx_fatal(FARGS
, "OpenMM exception caught while initializating: %s", e
.what());
1149 return NULL
; /* just to avoid warnings */
1153 * \brief Integrate one step.
1155 * \param[in] data OpenMMData object created by openmm_init().
1157 void openmm_take_one_step(void* data
)
1159 // static int step = 0; printf("----> taking step #%d\n", step++);
1162 static_cast<OpenMMData
*>(data
)->integrator
->step(1);
1164 catch (std::exception
& e
)
1166 gmx_fatal(FARGS
, "OpenMM exception caught while taking a step: %s", e
.what());
1171 * \brief Integrate n steps.
1173 * \param[in] data OpenMMData object created by openmm_init().
1175 void openmm_take_steps(void* data
, int nstep
)
1179 static_cast<OpenMMData
*>(data
)->integrator
->step(nstep
);
1181 catch (std::exception
& e
)
1183 gmx_fatal(FARGS
, "OpenMM exception caught while taking a step: %s", e
.what());
1188 * \brief Clean up the data structures cretead for OpenMM.
1190 * \param[in] log Log file pointer.
1191 * \param[in] data OpenMMData object created by openmm_init().
1193 void openmm_cleanup(FILE* fplog
, void* data
)
1195 OpenMMData
* d
= static_cast<OpenMMData
*>(data
);
1197 if (isStringEqNCase(d
->platformOpt
->getOptionValue("platform"), "CUDA"))
1199 /* post-simulation memtest */
1200 runMemtest(fplog
, -1, "Post", d
->platformOpt
);
1203 delete d
->integrator
;
1205 delete d
->platformOpt
;
1210 * \brief Copy the current state information from OpenMM into the Gromacs data structures.
1212 * This function results in the requested proprties to be copied from the
1213 * GPU to host. As this represents a bottleneck, the frequency of pulling data
1214 * should be minimized.
1216 * \param[in] data OpenMMData object created by openmm_init().
1217 * \param[out] time Simulation time for which the state was created.
1218 * \param[out] state State of the system: coordinates and velocities.
1219 * \param[out] f Forces.
1220 * \param[out] enerd Energies.
1221 * \param[in] includePos True if coordinates are requested.
1222 * \param[in] includeVel True if velocities are requested.
1223 * \param[in] includeForce True if forces are requested.
1224 * \param[in] includeEnergy True if energies are requested.
1226 void openmm_copy_state(void *data
,
1227 t_state
*state
, double *time
,
1228 rvec f
[], gmx_enerdata_t
*enerd
,
1229 bool includePos
, bool includeVel
, bool includeForce
, bool includeEnergy
)
1233 types
+= State::Positions
;
1235 types
+= State::Velocities
;
1237 types
+= State::Forces
;
1239 types
+= State::Energy
;
1244 State currentState
= static_cast<OpenMMData
*>(data
)->context
->getState(types
);
1245 int numAtoms
= static_cast<OpenMMData
*>(data
)->system
->getNumParticles();
1248 for (int i
= 0; i
< numAtoms
; i
++)
1250 Vec3 x
= currentState
.getPositions()[i
];
1251 state
->x
[i
][0] = x
[0];
1252 state
->x
[i
][1] = x
[1];
1253 state
->x
[i
][2] = x
[2];
1258 for (int i
= 0; i
< numAtoms
; i
++)
1260 Vec3 v
= currentState
.getVelocities()[i
];
1261 state
->v
[i
][0] = v
[0];
1262 state
->v
[i
][1] = v
[1];
1263 state
->v
[i
][2] = v
[2];
1268 for (int i
= 0; i
< numAtoms
; i
++)
1270 Vec3 force
= currentState
.getForces()[i
];
1278 int numConstraints
= static_cast<OpenMMData
*>(data
)->system
->getNumConstraints();
1279 int dof
= 3*numAtoms
-numConstraints
;
1280 if (static_cast<OpenMMData
*>(data
)->removeCM
)
1282 enerd
->term
[F_EPOT
] = currentState
.getPotentialEnergy();
1283 enerd
->term
[F_EKIN
] = currentState
.getKineticEnergy();
1284 enerd
->term
[F_ETOT
] = enerd
->term
[F_EPOT
] + enerd
->term
[F_EKIN
];
1285 enerd
->term
[F_TEMP
] = 2.0*enerd
->term
[F_EKIN
]/dof
/BOLTZ
;
1287 *time
= currentState
.getTime();
1289 catch (std::exception
& e
)
1291 gmx_fatal(FARGS
, "OpenMM exception caught while retrieving state information: %s", e
.what());