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)
75 #define COMBRULE_CHK_TOL 1e-6
76 #define COMBRULE_SIGMA(sig1, sig2) (((sig1) + (sig2))/2)
77 #define COMBRULE_EPS(eps1, eps2) (sqrt((eps1) * (eps2)))
80 * \brief Convert string to integer type.
81 * \param[in] s String to convert from.
82 * \param[in] f Basefield format flag that takes any of the following I/O
83 * manipulators: dec, hex, oct.
84 * \param[out] t Destination variable to convert to.
87 static bool from_string(T
& t
, const string
& s
, ios_base
& (*f
)(ios_base
&))
90 return !(iss
>> f
>> t
).fail();
94 * \brief Split string around a given delimiter.
95 * \param[in] s String to split.
96 * \param[in] delim Delimiter character.
97 * \returns Vector of strings found in \p s.
99 static vector
<string
> split(const string
&s
, char delim
)
101 vector
<string
> elems
;
104 while (getline(ss
, item
, delim
))
106 if (item
.length() != 0)
107 elems
.push_back(item
);
113 * \brief Split a string of the form "option=value" into "option" and "value" strings.
114 * This string corresponds to one option and the associated value from the option list
115 * in the mdrun -device argument.
117 * \param[in] s A string containing an "option=value" pair that needs to be split up.
118 * \param[out] opt The name of the option.
119 * \param[out] val Value of the option.
121 static void splitOptionValue(const string
&s
, string
&opt
, string
&val
)
123 size_t eqPos
= s
.find('=');
124 if (eqPos
!= string::npos
)
126 opt
= s
.substr(0, eqPos
);
127 if (eqPos
!= s
.length()) val
= s
.substr(eqPos
+1);
132 * \brief Compare two strings ignoring case.
133 * This function is in fact a wrapper around the gromacs function gmx_strncasecmp().
134 * \param[in] s1 String.
135 * \param[in] s2 String.
136 * \returns Similarly to the C function strncasecmp(), the return value is an
137 integer less than, equal to, or greater than 0 if \p s1 less than,
138 identical to, or greater than \p s2.
140 static bool isStringEqNCase(const string s1
, const string s2
)
142 return (gmx_strncasecmp(s1
.c_str(), s2
.c_str(), max(s1
.length(), s2
.length())) == 0);
146 * \brief Convert string to upper case.
148 * \param[in] s String to convert to uppercase.
149 * \returns The given string converted to uppercase.
151 static string
toUpper(const string
&s
)
154 std::transform(stmp
.begin(), stmp
.end(), stmp
.begin(), static_cast < int(*)(int) > (toupper
));
159 \name Sizes of constant device option arrays GmxOpenMMPlatformOptions#platforms,
160 GmxOpenMMPlatformOptions#memtests, GmxOpenMMPlatformOptions#deviceid,
161 GmxOpenMMPlatformOptions#force_dev. */
163 #define SIZEOF_PLATFORMS 1 // 2
164 #define SIZEOF_MEMTESTS 3
165 #define SIZEOF_DEVICEIDS 1
166 #define SIZEOF_FORCE_DEV 2
168 #define SIZEOF_CHECK_COMBRULE 2
171 /*! Possible platform options in the mdrun -device option. */
172 static const char *devOptStrings
[] = { "platform", "deviceid", "memtest", "force-device", "check-combrule" };
174 /*! Enumerated platform options in the mdrun -device option. */
184 * \brief Class to extract and manage the platform options in the mdrun -device option.
187 class GmxOpenMMPlatformOptions
190 GmxOpenMMPlatformOptions(const char *opt
);
191 ~GmxOpenMMPlatformOptions() { options
.clear(); }
192 string
getOptionValue(const string
&opt
);
193 void remOption(const string
&opt
);
196 void setOption(const string
&opt
, const string
&val
);
198 map
<string
, string
> options
; /*!< Data structure to store the option (name, value) pairs. */
200 static const char * const platforms
[SIZEOF_PLATFORMS
]; /*!< Available OpenMM platforms; size #SIZEOF_PLATFORMS */
201 static const char * const memtests
[SIZEOF_MEMTESTS
]; /*!< Available types of memory tests, also valid
202 any positive integer >=15; size #SIZEOF_MEMTESTS */
203 static const char * const deviceid
[SIZEOF_DEVICEIDS
]; /*!< Possible values for deviceid option;
204 also valid any positive integer; size #SIZEOF_DEVICEIDS */
205 static const char * const force_dev
[SIZEOF_FORCE_DEV
]; /*!< Possible values for for force-device option;
206 size #SIZEOF_FORCE_DEV */
207 static const char * const check_combrule
[SIZEOF_CHECK_COMBRULE
]; /* XXX temporary debug feature to
208 turn off combination rule check */
211 const char * const GmxOpenMMPlatformOptions::platforms
[SIZEOF_PLATFORMS
]
213 //= { "Reference", "CUDA" /*,"OpenCL"*/ };
214 const char * const GmxOpenMMPlatformOptions::memtests
[SIZEOF_MEMTESTS
]
215 = { "15", "full", "off" };
216 const char * const GmxOpenMMPlatformOptions::deviceid
[SIZEOF_DEVICEIDS
]
218 const char * const GmxOpenMMPlatformOptions::force_dev
[SIZEOF_FORCE_DEV
]
220 const char * const GmxOpenMMPlatformOptions::check_combrule
[SIZEOF_CHECK_COMBRULE
]
225 * Takes the option list, parses it, checks the options and their values for validity.
226 * When certain options are not provided by the user, as default value the first item
227 * of the respective constant array is taken (GmxOpenMMPlatformOptions#platforms,
228 * GmxOpenMMPlatformOptions#memtests, GmxOpenMMPlatformOptions#deviceid,
229 * GmxOpenMMPlatformOptions#force_dev).
230 * \param[in] optionString Option list part of the mdrun -device parameter.
232 GmxOpenMMPlatformOptions::GmxOpenMMPlatformOptions(const char *optionString
)
234 // set default values
235 setOption("platform", platforms
[0]);
236 setOption("memtest", memtests
[0]);
237 setOption("deviceid", deviceid
[0]);
238 setOption("force-device", force_dev
[0]);
239 setOption("check-combrule", check_combrule
[0]);
241 string
opt(optionString
);
243 // remove all whitespaces
244 opt
.erase(remove_if(opt
.begin(), opt
.end(), ::isspace
), opt
.end());
245 // tokenize around ","-s
246 vector
<string
> tokens
= split(opt
, ',');
248 for (vector
<string
>::iterator it
= tokens
.begin(); it
!= tokens
.end(); ++it
)
250 string opt
= "", val
= "";
251 splitOptionValue(*it
, opt
, val
);
253 if (isStringEqNCase(opt
, "platform"))
255 /* no check, this will fail if platform does not exist when we try to set it */
260 if (isStringEqNCase(opt
, "memtest"))
262 /* the value has to be an integer >15(s) or "full" OR "off" */
263 if (!isStringEqNCase(val
, "full") && !isStringEqNCase(val
, "off"))
266 if (!from_string
<int>(secs
, val
, std::dec
))
268 gmx_fatal(FARGS
, "Invalid value for option memtest option: \"%s\"!", val
.c_str());
272 gmx_fatal(FARGS
, "Incorrect value for memtest option (%d). "
273 "Memtest needs to run for at least 15s!", secs
);
280 if (isStringEqNCase(opt
, "deviceid"))
283 if (!from_string
<int>(id
, val
, std::dec
) )
285 gmx_fatal(FARGS
, "Invalid device id: \"%s\"!", val
.c_str());
291 if (isStringEqNCase(opt
, "force-device"))
294 if (!isStringEqNCase(val
, "yes") && !isStringEqNCase(val
, "no"))
296 gmx_fatal(FARGS
, "Invalid OpenMM force option: \"%s\"!", val
.c_str());
302 if (isStringEqNCase(opt
, "check-combrule"))
305 if (!isStringEqNCase(val
, "yes") && !isStringEqNCase(val
, "no"))
307 gmx_fatal(FARGS
, "Invalid OpenMM force option: \"%s\"!", val
.c_str());
314 // if we got till here something went wrong
315 gmx_fatal(FARGS
, "Invalid OpenMM platform option: \"%s\"!", (*it
).c_str());
321 * \brief Getter function.
322 * \param[in] opt Name of the option.
323 * \returns Returns the value associated to an option.
325 string
GmxOpenMMPlatformOptions::getOptionValue(const string
&opt
)
327 map
<string
, string
> :: const_iterator it
= options
.find(toUpper(opt
));
328 if (it
!= options
.end())
339 * \brief Setter function - private, only used from contructor.
340 * \param[in] opt Name of the option.
341 * \param[in] val Value for the option.
343 void GmxOpenMMPlatformOptions::setOption(const string
&opt
, const string
&val
)
345 options
[toUpper(opt
)] = val
;
349 * \brief Removes an option with its value from the map structure. If the option
350 * does not exist, returns without any action.
351 * \param[in] opt Name of the option.
353 void GmxOpenMMPlatformOptions::remOption(const string
&opt
)
355 options
.erase(toUpper(opt
));
359 * \brief Print option-value pairs to a file (debugging function).
361 void GmxOpenMMPlatformOptions::print()
363 cout
<< ">> Platform options: " << endl
364 << ">> platform = " << getOptionValue("platform") << endl
365 << ">> deviceID = " << getOptionValue("deviceid") << endl
366 << ">> memtest = " << getOptionValue("memtest") << endl
367 << ">> force-device = " << getOptionValue("force-device") << endl
;
371 * \brief Container for OpenMM related data structures that represent the bridge
372 * between the Gromacs data-structures and the OpenMM library and is but it's
373 * only passed through the API functions as void to disable direct access.
378 System
* system
; /*! The system to simulate. */
379 Context
* context
; /*! The OpenMM context in which the simulation is carried out. */
380 Integrator
* integrator
; /*! The integrator used in the simulation. */
381 bool removeCM
; /*! If \true remove venter of motion, false otherwise. */
382 GmxOpenMMPlatformOptions
*platformOpt
; /*! Platform options. */
386 * \brief Runs memtest on the GPU that has alreaby been initialized by OpenMM.
387 * \param[in] fplog Pointer to gromacs log file.
388 * \param[in] devId Device id of the GPU to run the test on.
389 Note: as OpenMM previously creates the context,for now this is always -1.
390 * \param[in] pre_post Contains either "Pre" or "Post" just to be able to differentiate in
391 * stdout messages/log between memtest carried out before and after simulation.
392 * \param[in] opt Pointer to platform options object.
394 static void runMemtest(FILE* fplog
, int devId
, const char* pre_post
, GmxOpenMMPlatformOptions
*opt
)
396 char strout_buf
[STRLEN
];
399 string s
= opt
->getOptionValue("memtest");
400 const char *test_type
= s
.c_str();
402 if (!gmx_strcasecmp(test_type
, "off"))
408 if (!gmx_strcasecmp(test_type
, "full"))
414 from_string
<int>(which_test
, test_type
, std::dec
);
420 gmx_fatal(FARGS
, "Amount of seconds for memetest is negative (%d). ", which_test
);
425 case 0: /* no memtest */
426 sprintf(strout_buf
, "%s-simulation GPU memtest skipped. Note, that faulty memory can cause "
427 "incorrect results!", pre_post
);
428 fprintf(fplog
, "%s\n", strout_buf
);
429 gmx_warning(strout_buf
);
432 case 1: /* quick memtest */
433 fprintf(fplog
, "%s-simulation %s GPU memtest in progress...\n", pre_post
, test_type
);
434 fprintf(stdout
, "\n%s-simulation %s GPU memtest in progress...", pre_post
, test_type
);
437 res
= do_quick_memtest(devId
);
440 case 2: /* full memtest */
441 fprintf(fplog
, "%s-simulation %s memtest in progress...\n", pre_post
, test_type
);
442 fprintf(stdout
, "\n%s-simulation %s memtest in progress...", pre_post
, test_type
);
445 res
= do_full_memtest(devId
);
448 default: /* timed memtest */
449 fprintf(fplog
, "%s-simulation ~%ds memtest in progress...\n", pre_post
, which_test
);
450 fprintf(stdout
, "\n%s-simulation ~%ds memtest in progress...", pre_post
, which_test
);
453 res
= do_timed_memtest(devId
, which_test
);
460 gmx_fatal(FARGS
, MEM_ERR_MSG(pre_post
));
464 fprintf(fplog
, "Memory test completed without errors.\n");
466 fprintf(stdout
, "done, no errors detected\n");
473 * \brief Convert Lennard-Jones parameters c12 and c6 to sigma and epsilon.
478 * \param[out] epsilon
480 static void convert_c_12_6(double c12
, double c6
, double *sigma
, double *epsilon
)
482 if (c12
== 0 && c6
== 0)
487 else if (c12
> 0 && c6
> 0)
489 *epsilon
= (c6
*c6
)/(4.0*c12
);
490 *sigma
= pow(c12
/c6
, 1.0/6.0);
494 gmx_fatal(FARGS
,"OpenMM only supports c6 > 0 and c12 > 0 or c6 = c12 = 0.");
499 * \brief Does gromacs option checking.
501 * Checks the gromacs mdp options for features unsupported in OpenMM, case in which
502 * interrupts the execution. It also warns the user about pecularities of OpenMM
504 * \param[in] fplog Gromacs log file pointer.
505 * \param[in] ir Gromacs input parameters, see ::t_inputrec
506 * \param[in] top Gromacs node local topology, \see gmx_localtop_t
507 * \param[in] state Gromacs state structure \see ::t_state
508 * \param[in] mdatoms Gromacs atom parameters, \see ::t_mdatoms
509 * \param[in] fr \see ::t_forcerec
510 * \param[in] state Gromacs systems state, \see ::t_state
512 static void checkGmxOptions(FILE* fplog
, GmxOpenMMPlatformOptions
*opt
,
513 t_inputrec
*ir
, gmx_localtop_t
*top
,
514 t_forcerec
*fr
, t_state
*state
)
516 char warn_buf
[STRLEN
];
518 double c6
, c12
, sigma_ij
, sigma_ji
, sigma_ii
, sigma_jj
, sigma_comb
,
519 eps_ij
, eps_ji
, eps_ii
, eps_jj
, eps_comb
;
521 /* Abort if unsupported critical options are present */
526 gmx_warning( "OpenMM does not support leap-frog, will use velocity-verlet integrator.");
529 if ( (ir
->eI
!= eiMD
) &&
531 (ir
->eI
!= eiVVAK
) &&
536 gmx_fatal(FARGS
, "OpenMM supports only the following integrators: md/md-vv/md-vv-avek, sd/sd1, and bd.");
540 if ( (ir
->coulombtype
!= eelPME
) &&
541 (ir
->coulombtype
!= eelRF
) &&
542 (ir
->coulombtype
!= eelEWALD
) &&
544 ( !(ir
->coulombtype
== eelCUT
&& ir
->rcoulomb
== 0 && ir
->rvdw
== 0)) )
546 gmx_fatal(FARGS
,"OpenMM supports only the following methods for electrostatics: "
547 "NoCutoff (i.e. rcoulomb = rvdw = 0 ),Reaction-Field, Ewald or PME.");
550 if (ir
->etc
!= etcNO
&&
555 gmx_warning("OpenMM supports only Andersen thermostat with the md/md-vv/md-vv-avek integrators.");
558 if (ir
->opts
.ngtc
> 1)
559 gmx_fatal(FARGS
,"OpenMM does not support multiple temperature coupling groups.");
561 if (ir
->epc
!= etcNO
)
562 gmx_fatal(FARGS
,"OpenMM does not support pressure coupling.");
564 if (ir
->opts
.annealing
[0])
565 gmx_fatal(FARGS
,"OpenMM does not support simulated annealing.");
567 if (top
->idef
.il
[F_CONSTR
].nr
> 0 && ir
->eConstrAlg
!= econtSHAKE
)
568 gmx_warning("OpenMM provides contraints as a combination "
569 "of SHAKE, SETTLE and CCMA. Accuracy is based on the SHAKE tolerance set "
570 "by the \"shake_tol\" option.");
573 gmx_fatal(FARGS
,"OpenMM does not support walls.");
575 if (ir
->ePull
!= epullNO
)
576 gmx_fatal(FARGS
,"OpenMM does not support pulling.");
578 /* TODO: DISABLED as it does not work with implicit solvent simulation */
580 /* check for restraints */
581 for (i
= 0; i
< F_EPOT
; i
++)
590 top
->idef
.il
[i
].nr
> 0)
592 /* TODO fix the message */
593 gmx_fatal(FARGS
, "OpenMM does not support (some) of the provided restaint(s).");
598 if (ir
->efep
!= efepNO
)
599 gmx_fatal(FARGS
,"OpenMM does not support free energy calculations.");
601 if (ir
->opts
.ngacc
> 1)
602 gmx_fatal(FARGS
,"OpenMM does not support non-equilibrium MD (accelerated groups).");
604 if (IR_ELEC_FIELD(*ir
))
605 gmx_fatal(FARGS
,"OpenMM does not support electric fields.");
608 gmx_fatal(FARGS
,"OpenMM does not support QMMM calculations.");
610 if (ir
->rcoulomb
!= ir
->rvdw
)
611 gmx_fatal(FARGS
,"OpenMM uses a single cutoff for both Coulomb "
612 "and VdW interactions. Please set rcoulomb equal to rvdw.");
614 if (EEL_FULL(ir
->coulombtype
))
616 if (ir
->ewald_geometry
== eewg3DC
)
617 gmx_fatal(FARGS
,"OpenMM supports only Ewald 3D geometry.");
618 if (ir
->epsilon_surface
!= 0)
619 gmx_fatal(FARGS
,"OpenMM does not support dipole correction in Ewald summation.");
622 if (TRICLINIC(state
->box
))
624 gmx_fatal(FARGS
,"OpenMM does not support triclinic unit cells.");
627 /* XXX this is just debugging code to disable the combination rule check */
628 if ( isStringEqNCase(opt
->getOptionValue("check-combrule"), "yes") )
630 /* As OpenMM by default uses hardcoded combination rules
631 sigma_ij = (sigma_i + sigma_j)/2, eps_ij = sqrt(eps_i * eps_j)
632 we need to check whether the force field params obey this
633 and if not, we can't use this force field so we exit
634 grace-fatal-fully. */
635 real
*nbfp
= fr
->nbfp
;
639 fprintf(debug
, ">> Atom parameters: <<\n%10s%5s %5s %5s %5s COMB\n",
640 "", "i-j", "j-i", "i-i", "j-j");
642 /* loop over all i-j atom pairs and verify if
643 sigma_ij = sigma_ji = sigma_comb and eps_ij = eps_ji = eps_comb */
644 for (i
= 0; i
< natoms
; i
++)
647 c12
= C12(nbfp
, natoms
, i
, i
);
648 c6
= C6(nbfp
, natoms
, i
, i
);
649 convert_c_12_6(c12
, c6
, &sigma_ii
, &eps_ii
);
651 for (j
= 0; j
< i
; j
++)
654 c12
= C12(nbfp
, natoms
, i
, j
);
655 c6
= C6(nbfp
, natoms
, i
, j
);
656 convert_c_12_6(c12
, c6
, &sigma_ij
, &eps_ij
);
658 c12
= C12(nbfp
, natoms
, j
, i
);
659 c6
= C6(nbfp
, natoms
, j
, i
);
660 convert_c_12_6(c12
, c6
, &sigma_ji
, &eps_ji
);
662 c12
= C12(nbfp
, natoms
, j
, j
);
663 c6
= C6(nbfp
, natoms
, j
, j
);
664 convert_c_12_6(c12
, c6
, &sigma_jj
, &eps_jj
);
665 /* OpenMM hardcoded combination rules */
666 sigma_comb
= COMBRULE_SIGMA(sigma_ii
, sigma_jj
);
667 eps_comb
= COMBRULE_EPS(eps_ii
, eps_jj
);
671 fprintf(debug
, "i=%-3d j=%-3d", i
, j
);
672 fprintf(debug
, "%-11s", "sigma");
673 fprintf(debug
, "%5.3f %5.3f %5.3f %5.3f %5.3f\n",
674 sigma_ij
, sigma_ji
, sigma_ii
, sigma_jj
, sigma_comb
);
675 fprintf(debug
, "%11s%-11s", "", "epsilon");
676 fprintf(debug
, "%5.3f %5.3f %5.3f %5.3f %5.3f\n",
677 eps_ij
, eps_ji
, eps_ii
, eps_jj
, eps_comb
);
680 /* check the values against the rule used by omm */
681 if((fabs(eps_ij
) > COMBRULE_CHK_TOL
&&
682 fabs(eps_ji
) > COMBRULE_CHK_TOL
) &&
683 (fabs(sigma_comb
- sigma_ij
) > COMBRULE_CHK_TOL
||
684 fabs(sigma_comb
- sigma_ji
) > COMBRULE_CHK_TOL
||
685 fabs(eps_comb
- eps_ij
) > COMBRULE_CHK_TOL
||
686 fabs(eps_comb
- eps_ji
) > COMBRULE_CHK_TOL
))
689 "The combination rules of the used force-field do not "
690 "match the one supported by OpenMM: "
691 "sigma_ij = (sigma_i + sigma_j)/2, eps_ij = sqrt(eps_i * eps_j). "
692 "Switch to a force-field that uses these rules in order to "
693 "simulate this system using OpenMM.\n");
697 if (debug
) { fprintf(debug
, ">><<\n\n"); }
699 /* if we got here, log that everything is fine */
702 fprintf(debug
, ">> The combination rule of the used force matches the one used by OpenMM.\n");
704 fprintf(fplog
, "The combination rule of the force used field matches the one used by OpenMM.\n");
706 } /* if (are we checking the combination rules) ... */
711 * \brief Initialize OpenMM, run sanity/consistency checks, and return a pointer to
714 * Various gromacs data structures are passed that contain the parameters, state and
715 * other porperties of the system to simulate. These serve as input for initializing
716 * OpenMM. Besides, a set of misc action are taken:
717 * - OpenMM plugins are loaded;
718 * - platform options in \p platformOptStr are parsed and checked;
719 * - Gromacs parameters are checked for OpenMM support and consistency;
720 * - after the OpenMM is initialized memtest executed in the same GPU context.
722 * \param[in] fplog Gromacs log file handler.
723 * \param[in] platformOptStr Platform option string.
724 * \param[in] ir The Gromacs input parameters, see ::t_inputrec
725 * \param[in] top_global Gromacs system toppology, \see ::gmx_mtop_t
726 * \param[in] top Gromacs node local topology, \see gmx_localtop_t
727 * \param[in] mdatoms Gromacs atom parameters, \see ::t_mdatoms
728 * \param[in] fr \see ::t_forcerec
729 * \param[in] state Gromacs systems state, \see ::t_state
730 * \returns Pointer to a
733 void* openmm_init(FILE *fplog
, const char *platformOptStr
,
735 gmx_mtop_t
*top_global
, gmx_localtop_t
*top
,
736 t_mdatoms
*mdatoms
, t_forcerec
*fr
, t_state
*state
)
739 char warn_buf
[STRLEN
];
740 static bool hasLoadedPlugins
= false;
741 string usedPluginDir
;
746 if (!hasLoadedPlugins
)
748 vector
<string
> loadedPlugins
;
749 /* Look for OpenMM plugins at various locations (listed in order of priority):
750 - on the path in OPENMM_PLUGIN_DIR environment variable if this is specified
751 - on the path in the OPENMM_PLUGIN_DIR macro that is set by the build script
752 - at the default location assumed by OpenMM
755 char *pluginDir
= getenv("OPENMM_PLUGIN_DIR");
757 /* no env var or empty */
758 if (pluginDir
!= NULL
&& *pluginDir
!= '\0')
760 loadedPlugins
= Platform::loadPluginsFromDirectory(pluginDir
);
761 if (loadedPlugins
.size() > 0)
763 hasLoadedPlugins
= true;
764 usedPluginDir
= pluginDir
;
768 gmx_fatal(FARGS
, "The directory provided in the OPENMM_PLUGIN_DIR environment variable "
769 "(%s) does not contain valid OpenMM plugins. Check your OpenMM installation!",
774 /* macro set at build time */
775 #ifdef OPENMM_PLUGIN_DIR
776 if (!hasLoadedPlugins
)
778 loadedPlugins
= Platform::loadPluginsFromDirectory(OPENMM_PLUGIN_DIR
);
779 if (loadedPlugins
.size() > 0)
781 hasLoadedPlugins
= true;
782 usedPluginDir
= OPENMM_PLUGIN_DIR
;
786 /* default loocation */
787 if (!hasLoadedPlugins
)
789 loadedPlugins
= Platform::loadPluginsFromDirectory(Platform::getDefaultPluginsDirectory());
790 if (loadedPlugins
.size() > 0)
792 hasLoadedPlugins
= true;
793 usedPluginDir
= Platform::getDefaultPluginsDirectory();
797 /* if there are still no plugins loaded there won't be any */
798 if (!hasLoadedPlugins
)
800 gmx_fatal(FARGS
, "No OpenMM plugins were found! You can provide the"
801 " plugin directory in the OPENMM_PLUGIN_DIR environment variable.", pluginDir
);
804 fprintf(fplog
, "\nPlugins loaded from directory %s:\t", usedPluginDir
.c_str());
805 for (int i
= 0; i
< (int)loadedPlugins
.size(); i
++)
807 fprintf(fplog
, "%s, ", loadedPlugins
[i
].c_str());
809 fprintf(fplog
, "\n");
812 /* parse option string */
813 GmxOpenMMPlatformOptions
*opt
= new GmxOpenMMPlatformOptions(platformOptStr
);
814 devId
= atoi(opt
->getOptionValue("deviceid").c_str());
821 /* check wheter Gromacs options compatibility with OpenMM */
822 checkGmxOptions(fplog
, opt
, ir
, top
, fr
, state
);
824 // Create the system.
825 const t_idef
& idef
= top
->idef
;
826 const int numAtoms
= top_global
->natoms
;
827 const int numConstraints
= idef
.il
[F_CONSTR
].nr
/3;
828 const int numSettle
= idef
.il
[F_SETTLE
].nr
/2;
829 const int numBonds
= idef
.il
[F_BONDS
].nr
/3;
830 const int numUB
= idef
.il
[F_UREY_BRADLEY
].nr
/4;
831 const int numAngles
= idef
.il
[F_ANGLES
].nr
/4;
832 const int numPeriodic
= idef
.il
[F_PDIHS
].nr
/5;
833 const int numRB
= idef
.il
[F_RBDIHS
].nr
/5;
834 const int num14
= idef
.il
[F_LJ14
].nr
/3;
835 System
* sys
= new System();
837 sys
->addForce(new CMMotionRemover(ir
->nstcomm
));
839 // Set bonded force field terms.
840 const int* bondAtoms
= (int*) idef
.il
[F_BONDS
].iatoms
;
841 HarmonicBondForce
* bondForce
= new HarmonicBondForce();
842 sys
->addForce(bondForce
);
844 for (int i
= 0; i
< numBonds
; ++i
)
846 int type
= bondAtoms
[offset
++];
847 int atom1
= bondAtoms
[offset
++];
848 int atom2
= bondAtoms
[offset
++];
849 bondForce
->addBond(atom1
, atom2
,
850 idef
.iparams
[type
].harmonic
.rA
, idef
.iparams
[type
].harmonic
.krA
);
852 // Urey-Bradley includes both the angle and bond potential for 1-3 interactions
853 const int* ubAtoms
= (int*) idef
.il
[F_UREY_BRADLEY
].iatoms
;
854 HarmonicBondForce
* ubBondForce
= new HarmonicBondForce();
855 HarmonicAngleForce
* ubAngleForce
= new HarmonicAngleForce();
856 sys
->addForce(ubBondForce
);
857 sys
->addForce(ubAngleForce
);
859 for (int i
= 0; i
< numUB
; ++i
)
861 int type
= ubAtoms
[offset
++];
862 int atom1
= ubAtoms
[offset
++];
863 int atom2
= ubAtoms
[offset
++];
864 int atom3
= ubAtoms
[offset
++];
865 ubBondForce
->addBond(atom1
, atom3
,
866 idef
.iparams
[type
].u_b
.r13
, idef
.iparams
[type
].u_b
.kUB
);
867 ubAngleForce
->addAngle(atom1
, atom2
, atom3
,
868 idef
.iparams
[type
].u_b
.theta
*M_PI
/180.0, idef
.iparams
[type
].u_b
.ktheta
);
870 const int* angleAtoms
= (int*) idef
.il
[F_ANGLES
].iatoms
;
871 HarmonicAngleForce
* angleForce
= new HarmonicAngleForce();
872 sys
->addForce(angleForce
);
874 for (int i
= 0; i
< numAngles
; ++i
)
876 int type
= angleAtoms
[offset
++];
877 int atom1
= angleAtoms
[offset
++];
878 int atom2
= angleAtoms
[offset
++];
879 int atom3
= angleAtoms
[offset
++];
880 angleForce
->addAngle(atom1
, atom2
, atom3
,
881 idef
.iparams
[type
].harmonic
.rA
*M_PI
/180.0, idef
.iparams
[type
].harmonic
.krA
);
883 const int* periodicAtoms
= (int*) idef
.il
[F_PDIHS
].iatoms
;
884 PeriodicTorsionForce
* periodicForce
= new PeriodicTorsionForce();
885 sys
->addForce(periodicForce
);
887 for (int i
= 0; i
< numPeriodic
; ++i
)
889 int type
= periodicAtoms
[offset
++];
890 int atom1
= periodicAtoms
[offset
++];
891 int atom2
= periodicAtoms
[offset
++];
892 int atom3
= periodicAtoms
[offset
++];
893 int atom4
= periodicAtoms
[offset
++];
894 periodicForce
->addTorsion(atom1
, atom2
, atom3
, atom4
,
895 idef
.iparams
[type
].pdihs
.mult
,
896 idef
.iparams
[type
].pdihs
.phiA
*M_PI
/180.0,
897 idef
.iparams
[type
].pdihs
.cpA
);
899 const int* rbAtoms
= (int*) idef
.il
[F_RBDIHS
].iatoms
;
900 RBTorsionForce
* rbForce
= new RBTorsionForce();
901 sys
->addForce(rbForce
);
903 for (int i
= 0; i
< numRB
; ++i
)
905 int type
= rbAtoms
[offset
++];
906 int atom1
= rbAtoms
[offset
++];
907 int atom2
= rbAtoms
[offset
++];
908 int atom3
= rbAtoms
[offset
++];
909 int atom4
= rbAtoms
[offset
++];
910 rbForce
->addTorsion(atom1
, atom2
, atom3
, atom4
,
911 idef
.iparams
[type
].rbdihs
.rbcA
[0], idef
.iparams
[type
].rbdihs
.rbcA
[1],
912 idef
.iparams
[type
].rbdihs
.rbcA
[2], idef
.iparams
[type
].rbdihs
.rbcA
[3],
913 idef
.iparams
[type
].rbdihs
.rbcA
[4], idef
.iparams
[type
].rbdihs
.rbcA
[5]);
916 // Set nonbonded parameters and masses.
917 int ntypes
= fr
->ntype
;
918 int* types
= mdatoms
->typeA
;
919 real
* nbfp
= fr
->nbfp
;
920 real
* charges
= mdatoms
->chargeA
;
921 real
* masses
= mdatoms
->massT
;
922 NonbondedForce
* nonbondedForce
= new NonbondedForce();
923 sys
->addForce(nonbondedForce
);
928 if (ir
->rcoulomb
== 0)
930 nonbondedForce
->setNonbondedMethod(NonbondedForce::NoCutoff
);
934 nonbondedForce
->setNonbondedMethod(NonbondedForce::CutoffNonPeriodic
);
938 switch (ir
->coulombtype
)
941 nonbondedForce
->setNonbondedMethod(NonbondedForce::CutoffPeriodic
);
945 nonbondedForce
->setNonbondedMethod(NonbondedForce::Ewald
);
949 nonbondedForce
->setNonbondedMethod(NonbondedForce::PME
);
953 gmx_fatal(FARGS
,"Internal error: you should not see this message, it that the"
954 "electrosatics option check failed. Please report this error!");
956 sys
->setPeriodicBoxVectors(Vec3(state
->box
[0][0], 0, 0),
957 Vec3(0, state
->box
[1][1], 0), Vec3(0, 0, state
->box
[2][2]));
958 nonbondedForce
->setCutoffDistance(ir
->rcoulomb
);
962 gmx_fatal(FARGS
,"OpenMM supports only full periodic boundary conditions "
963 "(pbc = xyz), or none (pbc = no).");
967 /* Fix for PME and Ewald error tolerance
969 * OpenMM uses approximate formulas to calculate the Ewald parameter:
970 * alpha = (1.0/cutoff)*sqrt(-log(2.0*tolerlance));
971 * and the grid spacing for PME:
972 * gridX = ceil(alpha*box[0][0]/pow(0.5*tol, 0.2));
973 * gridY = ceil(alpha*box[1][1]/pow(0.5*tol, 0.2));
974 * gridZ = ceil(alpha*box[2][2]/pow(0.5*tol, 0.2));
976 * It overestimates the precision and setting it to
977 * (500 x ewald_rtol) seems to give a reasonable match to the GROMACS settings
979 * If the default ewald_rtol=1e-5 is used we silently adjust the value,
980 * otherwise a warning is issued about the action taken.
982 double corr_ewald_rtol
= 500.0 * ir
->ewald_rtol
;
983 if ((ir
->ePBC
== epbcXYZ
) &&
984 (ir
->coulombtype
== eelEWALD
|| ir
->coulombtype
== eelPME
))
988 fprintf(debug
, ">> ewald_rtol = %e (corrected = %e) \n",
989 ir
->ewald_rtol
, corr_ewald_rtol
);
992 if (fabs(ir
->ewald_rtol
- 1e-5) > 1e-10)
994 gmx_warning("OpenMM uses the ewald_rtol parameter with approximate formulas "
995 "to calculate the alpha and grid spacing parameters of the Ewald "
996 "and PME methods. This tolerance need to be corrected in order to get "
997 "settings close to the ones used in GROMACS. Although the internal correction "
998 "should work for any reasonable value of ewald_rtol, using values other than "
999 "the default 1e-5 might cause incorrect behavior.");
1001 if (corr_ewald_rtol
> 1)
1003 gmx_fatal(FARGS
, "The ewald_rtol accuracy term is >1 after the "
1004 "adjustment for OpenMM (%e)", corr_ewald_rtol
);
1007 nonbondedForce
->setEwaldErrorTolerance(corr_ewald_rtol
);
1010 for (int i
= 0; i
< numAtoms
; ++i
)
1012 double c12
= nbfp
[types
[i
]*2*ntypes
+types
[i
]*2+1];
1013 double c6
= nbfp
[types
[i
]*2*ntypes
+types
[i
]*2];
1014 double sigma
, epsilon
;
1015 convert_c_12_6(c12
, c6
, &sigma
, &epsilon
);
1016 nonbondedForce
->addParticle(charges
[i
], sigma
, epsilon
);
1017 sys
->addParticle(masses
[i
]);
1020 // Build a table of all exclusions.
1021 vector
<set
<int> > exclusions(numAtoms
);
1022 for (int i
= 0; i
< numAtoms
; i
++)
1024 int start
= top
->excls
.index
[i
];
1025 int end
= top
->excls
.index
[i
+1];
1026 for (int j
= start
; j
< end
; j
++)
1027 exclusions
[i
].insert(top
->excls
.a
[j
]);
1030 // Record the 1-4 interactions, and remove them from the list of exclusions.
1031 const int* nb14Atoms
= (int*) idef
.il
[F_LJ14
].iatoms
;
1033 for (int i
= 0; i
< num14
; ++i
)
1035 int type
= nb14Atoms
[offset
++];
1036 int atom1
= nb14Atoms
[offset
++];
1037 int atom2
= nb14Atoms
[offset
++];
1038 double sigma
, epsilon
;
1039 convert_c_12_6(idef
.iparams
[type
].lj14
.c12A
,
1040 idef
.iparams
[type
].lj14
.c6A
,
1042 nonbondedForce
->addException(atom1
, atom2
,
1043 fr
->fudgeQQ
*charges
[atom1
]*charges
[atom2
], sigma
, epsilon
);
1044 exclusions
[atom1
].erase(atom2
);
1045 exclusions
[atom2
].erase(atom1
);
1048 // Record exclusions.
1049 for (int i
= 0; i
< numAtoms
; i
++)
1051 for (set
<int>::const_iterator iter
= exclusions
[i
].begin(); iter
!= exclusions
[i
].end(); ++iter
)
1055 nonbondedForce
->addException(i
, *iter
, 0.0, 1.0, 0.0);
1060 // Add GBSA if needed.
1061 if (ir
->implicit_solvent
== eisGBSA
)
1063 t_atoms atoms
= gmx_mtop_global_atoms(top_global
);
1064 GBSAOBCForce
* gbsa
= new GBSAOBCForce();
1066 sys
->addForce(gbsa
);
1067 gbsa
->setSoluteDielectric(ir
->epsilon_r
);
1068 gbsa
->setSolventDielectric(ir
->gb_epsilon_solvent
);
1069 gbsa
->setCutoffDistance(nonbondedForce
->getCutoffDistance());
1070 if (nonbondedForce
->getNonbondedMethod() == NonbondedForce::NoCutoff
)
1071 gbsa
->setNonbondedMethod(GBSAOBCForce::NoCutoff
);
1072 else if (nonbondedForce
->getNonbondedMethod() == NonbondedForce::CutoffNonPeriodic
)
1073 gbsa
->setNonbondedMethod(GBSAOBCForce::CutoffNonPeriodic
);
1074 else if (nonbondedForce
->getNonbondedMethod() == NonbondedForce::CutoffPeriodic
)
1075 gbsa
->setNonbondedMethod(GBSAOBCForce::CutoffPeriodic
);
1077 gmx_fatal(FARGS
,"OpenMM supports only Reaction-Field electrostatics with OBC/GBSA.");
1079 for (int i
= 0; i
< numAtoms
; ++i
)
1081 gbsa
->addParticle(charges
[i
],
1082 top_global
->atomtypes
.gb_radius
[atoms
.atom
[i
].type
],
1083 top_global
->atomtypes
.S_hct
[atoms
.atom
[i
].type
]);
1085 free_t_atoms(&atoms
, FALSE
);
1089 const int* constraintAtoms
= (int*) idef
.il
[F_CONSTR
].iatoms
;
1091 for (int i
= 0; i
< numConstraints
; ++i
)
1093 int type
= constraintAtoms
[offset
++];
1094 int atom1
= constraintAtoms
[offset
++];
1095 int atom2
= constraintAtoms
[offset
++];
1096 sys
->addConstraint(atom1
, atom2
, idef
.iparams
[type
].constr
.dA
);
1098 const int* settleAtoms
= (int*) idef
.il
[F_SETTLE
].iatoms
;
1100 for (int i
= 0; i
< numSettle
; ++i
)
1102 int type
= settleAtoms
[offset
++];
1103 int oxygen
= settleAtoms
[offset
++];
1104 sys
->addConstraint(oxygen
, oxygen
+1, idef
.iparams
[type
].settle
.doh
);
1105 sys
->addConstraint(oxygen
, oxygen
+2, idef
.iparams
[type
].settle
.doh
);
1106 sys
->addConstraint(oxygen
+1, oxygen
+2, idef
.iparams
[type
].settle
.dhh
);
1109 // Create an integrator for simulating the system.
1110 double friction
= (ir
->opts
.tau_t
[0] == 0.0 ? 0.0 : 1.0/ir
->opts
.tau_t
[0]);
1114 integ
= new BrownianIntegrator(ir
->opts
.ref_t
[0], friction
, ir
->delta_t
);
1115 static_cast<BrownianIntegrator
*>(integ
)->setRandomNumberSeed(ir
->ld_seed
);
1117 else if (EI_SD(ir
->eI
))
1119 integ
= new LangevinIntegrator(ir
->opts
.ref_t
[0], friction
, ir
->delta_t
);
1120 static_cast<LangevinIntegrator
*>(integ
)->setRandomNumberSeed(ir
->ld_seed
);
1124 integ
= new VerletIntegrator(ir
->delta_t
);
1125 if ( ir
->etc
!= etcNO
)
1127 real collisionFreq
= ir
->opts
.tau_t
[0] / 1000; /* tau_t (ps) / 1000 = collisionFreq (fs^-1) */
1128 AndersenThermostat
* thermostat
= new AndersenThermostat(ir
->opts
.ref_t
[0], friction
);
1129 sys
->addForce(thermostat
);
1132 integ
->setConstraintTolerance(ir
->shake_tol
);
1134 // Create a context and initialize it.
1135 Context
* context
= NULL
;
1138 OpenMM could automatically select the "best" GPU, however we're not't
1139 going to let it do that for now, as the current algorithm is very rudimentary
1140 and we anyway support only CUDA.
1141 if (platformOptStr == NULL || platformOptStr == "")
1143 context = new Context(*sys, *integ);
1148 /* which platform should we use */
1149 for (int i
= 0; i
< (int)Platform::getNumPlatforms() && context
== NULL
; i
++)
1151 if (isStringEqNCase(opt
->getOptionValue("platform"), Platform::getPlatform(i
).getName()))
1153 Platform
& platform
= Platform::getPlatform(i
);
1154 // set standard properties
1155 platform
.setPropertyDefaultValue("CudaDevice", opt
->getOptionValue("deviceid"));
1156 // TODO add extra properties
1157 context
= new Context(*sys
, *integ
, platform
);
1160 if (context
== NULL
)
1162 gmx_fatal(FARGS
, "The requested platform \"%s\" could not be found.",
1163 opt
->getOptionValue("platform").c_str());
1167 Platform
& platform
= context
->getPlatform();
1168 fprintf(fplog
, "Gromacs will use the OpenMM platform: %s\n", platform
.getName().c_str());
1170 const vector
<string
>& properties
= platform
.getPropertyNames();
1173 for (int i
= 0; i
< (int)properties
.size(); i
++)
1175 fprintf(debug
, ">> %s: %s\n", properties
[i
].c_str(),
1176 platform
.getPropertyValue(*context
, properties
[i
]).c_str());
1181 if (isStringEqNCase(opt
->getOptionValue("platform"), "CUDA"))
1183 /* For now this is just to double-check if OpenMM selected the GPU we wanted,
1184 but when we'll let OpenMM select the GPU automatically, it will query the devideId.
1187 if (!from_string
<int>(tmp
, platform
.getPropertyValue(*context
, "CudaDevice"), std::dec
))
1189 gmx_fatal(FARGS
, "Internal error: couldn't determine the device selected by OpenMM");
1192 gmx_fatal(FARGS
, "Internal error: OpenMM is using device #%d"
1193 "while initialized for device #%d", tmp
, devId
);
1197 /* check GPU compatibility */
1198 char gpuname
[STRLEN
];
1199 devId
= atoi(opt
->getOptionValue("deviceid").c_str());
1200 if (!is_supported_cuda_gpu(-1, gpuname
))
1202 if (!gmx_strcasecmp(opt
->getOptionValue("force-device").c_str(), "yes"))
1204 sprintf(warn_buf
, "Non-supported GPU selected (#%d, %s), forced continuing."
1205 "Note, that the simulation can be slow or it migth even crash.",
1207 fprintf(fplog
, "%s\n", warn_buf
);
1208 gmx_warning(warn_buf
);
1212 gmx_fatal(FARGS
, "The selected GPU (#%d, %s) is not supported by Gromacs! "
1213 "Most probably you have a low-end GPU which would not perform well, "
1214 "or new hardware that has not been tested with the current release. "
1215 "If you still want to try using the device, use the force-device=yes option.",
1221 fprintf(fplog
, "Gromacs will run on the GPU #%d (%s).\n", devId
, gpuname
);
1226 if (isStringEqNCase(opt
->getOptionValue("platform"), "CUDA"))
1228 /* pre-simulation memtest */
1229 runMemtest(fplog
, -1, "Pre", opt
);
1232 vector
<Vec3
> pos(numAtoms
);
1233 vector
<Vec3
> vel(numAtoms
);
1234 for (int i
= 0; i
< numAtoms
; ++i
)
1236 pos
[i
] = Vec3(state
->x
[i
][0], state
->x
[i
][1], state
->x
[i
][2]);
1237 vel
[i
] = Vec3(state
->v
[i
][0], state
->v
[i
][1], state
->v
[i
][2]);
1239 context
->setPositions(pos
);
1240 context
->setVelocities(vel
);
1242 // Return a structure containing the system, integrator, and context.
1243 OpenMMData
* data
= new OpenMMData();
1245 data
->integrator
= integ
;
1246 data
->context
= context
;
1247 data
->removeCM
= (ir
->nstcomm
> 0);
1248 data
->platformOpt
= opt
;
1251 catch (std::exception
& e
)
1253 gmx_fatal(FARGS
, "OpenMM exception caught while initializating: %s", e
.what());
1255 return NULL
; /* just to avoid warnings */
1259 * \brief Integrate one step.
1261 * \param[in] data OpenMMData object created by openmm_init().
1263 void openmm_take_one_step(void* data
)
1265 // static int step = 0; printf("----> taking step #%d\n", step++);
1268 static_cast<OpenMMData
*>(data
)->integrator
->step(1);
1270 catch (std::exception
& e
)
1272 gmx_fatal(FARGS
, "OpenMM exception caught while taking a step: %s", e
.what());
1277 * \brief Integrate n steps.
1279 * \param[in] data OpenMMData object created by openmm_init().
1281 void openmm_take_steps(void* data
, int nstep
)
1285 static_cast<OpenMMData
*>(data
)->integrator
->step(nstep
);
1287 catch (std::exception
& e
)
1289 gmx_fatal(FARGS
, "OpenMM exception caught while taking a step: %s", e
.what());
1294 * \brief Clean up the data structures cretead for OpenMM.
1296 * \param[in] log Log file pointer.
1297 * \param[in] data OpenMMData object created by openmm_init().
1299 void openmm_cleanup(FILE* fplog
, void* data
)
1301 OpenMMData
* d
= static_cast<OpenMMData
*>(data
);
1303 if (isStringEqNCase(d
->platformOpt
->getOptionValue("platform"), "CUDA"))
1305 /* post-simulation memtest */
1306 runMemtest(fplog
, -1, "Post", d
->platformOpt
);
1309 delete d
->integrator
;
1311 delete d
->platformOpt
;
1316 * \brief Copy the current state information from OpenMM into the Gromacs data structures.
1318 * This function results in the requested proprties to be copied from the
1319 * GPU to host. As this represents a bottleneck, the frequency of pulling data
1320 * should be minimized.
1322 * \param[in] data OpenMMData object created by openmm_init().
1323 * \param[out] time Simulation time for which the state was created.
1324 * \param[out] state State of the system: coordinates and velocities.
1325 * \param[out] f Forces.
1326 * \param[out] enerd Energies.
1327 * \param[in] includePos True if coordinates are requested.
1328 * \param[in] includeVel True if velocities are requested.
1329 * \param[in] includeForce True if forces are requested.
1330 * \param[in] includeEnergy True if energies are requested.
1332 void openmm_copy_state(void *data
,
1333 t_state
*state
, double *time
,
1334 rvec f
[], gmx_enerdata_t
*enerd
,
1335 bool includePos
, bool includeVel
, bool includeForce
, bool includeEnergy
)
1339 types
+= State::Positions
;
1341 types
+= State::Velocities
;
1343 types
+= State::Forces
;
1345 types
+= State::Energy
;
1350 State currentState
= static_cast<OpenMMData
*>(data
)->context
->getState(types
);
1351 int numAtoms
= static_cast<OpenMMData
*>(data
)->system
->getNumParticles();
1354 for (int i
= 0; i
< numAtoms
; i
++)
1356 Vec3 x
= currentState
.getPositions()[i
];
1357 state
->x
[i
][0] = x
[0];
1358 state
->x
[i
][1] = x
[1];
1359 state
->x
[i
][2] = x
[2];
1364 for (int i
= 0; i
< numAtoms
; i
++)
1366 Vec3 v
= currentState
.getVelocities()[i
];
1367 state
->v
[i
][0] = v
[0];
1368 state
->v
[i
][1] = v
[1];
1369 state
->v
[i
][2] = v
[2];
1374 for (int i
= 0; i
< numAtoms
; i
++)
1376 Vec3 force
= currentState
.getForces()[i
];
1384 int numConstraints
= static_cast<OpenMMData
*>(data
)->system
->getNumConstraints();
1385 int dof
= 3*numAtoms
-numConstraints
;
1386 if (static_cast<OpenMMData
*>(data
)->removeCM
)
1388 enerd
->term
[F_EPOT
] = currentState
.getPotentialEnergy();
1389 enerd
->term
[F_EKIN
] = currentState
.getKineticEnergy();
1390 enerd
->term
[F_ETOT
] = enerd
->term
[F_EPOT
] + enerd
->term
[F_EKIN
];
1391 enerd
->term
[F_TEMP
] = 2.0*enerd
->term
[F_EKIN
]/dof
/BOLTZ
;
1393 *time
= currentState
.getTime();
1395 catch (std::exception
& e
)
1397 gmx_fatal(FARGS
, "OpenMM exception caught while retrieving state information: %s", e
.what());