2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
39 * Implements Interactive Molecular Dynamics
41 * Re-implementation of basic IMD functions to work with VMD,
42 * see imdsocket.h for references to the IMD API.
44 * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
57 #include "gromacs/commandline/filenm.h"
58 #include "gromacs/domdec/domdec_struct.h"
59 #include "gromacs/domdec/ga2la.h"
60 #include "gromacs/fileio/confio.h"
61 #include "gromacs/fileio/gmxfio.h"
62 #include "gromacs/fileio/xvgr.h"
63 #include "gromacs/gmxlib/network.h"
64 #include "gromacs/imd/imdsocket.h"
65 #include "gromacs/math/units.h"
66 #include "gromacs/math/vec.h"
67 #include "gromacs/mdlib/broadcaststructs.h"
68 #include "gromacs/mdlib/groupcoord.h"
69 #include "gromacs/mdlib/sighandler.h"
70 #include "gromacs/mdlib/stat.h"
71 #include "gromacs/mdrunutility/multisim.h"
72 #include "gromacs/mdtypes/enerdata.h"
73 #include "gromacs/mdtypes/imdmodule.h"
74 #include "gromacs/mdtypes/inputrec.h"
75 #include "gromacs/mdtypes/md_enums.h"
76 #include "gromacs/mdtypes/mdrunoptions.h"
77 #include "gromacs/mdtypes/state.h"
78 #include "gromacs/pbcutil/pbc.h"
79 #include "gromacs/timing/wallcycle.h"
80 #include "gromacs/topology/mtop_util.h"
81 #include "gromacs/topology/topology.h"
82 #include "gromacs/utility/fatalerror.h"
83 #include "gromacs/utility/logger.h"
84 #include "gromacs/utility/smalloc.h"
85 #include "gromacs/utility/stringutil.h"
90 /*! \brief How long shall we wait in seconds until we check for a connection again? */
91 constexpr int c_loopWait = 1;
93 /*! \brief How long shall we check for the IMD_GO? */
94 constexpr int c_connectWait = 1;
96 /*! \brief IMD Header Size. */
97 constexpr int c_headerSize = 8;
99 /*! \brief IMD Protocol Version. */
100 constexpr int c_protocolVersion = 2;
104 * IMD (interactive molecular dynamics) energy record.
106 * As in the original IMD implementation. Energies in kcal/mol.
107 * NOTE: We return the energies in GROMACS / SI units,
108 * so they also show up as SI in VMD.
113 int32_t tstep; /**< time step */
114 float T_abs; /**< absolute temperature */
115 float E_tot; /**< total energy */
116 float E_pot; /**< potential energy */
117 float E_vdw; /**< van der Waals energy */
118 float E_coul; /**< Coulomb interaction energy */
119 float E_bond; /**< bonds energy */
120 float E_angle; /**< angles energy */
121 float E_dihe; /**< dihedrals energy */
122 float E_impr; /**< improper dihedrals energy */
127 * \brief IMD (interactive molecular dynamics) communication structure.
129 * This structure defines the IMD communication message header & protocol version.
133 int32_t type; /**< Type of IMD message, see IMDType_t above */
134 int32_t length; /**< Length */
139 * \brief Implementation type for the IMD session
141 * \todo Make this class (or one extracted from it) model
143 * \todo Use RAII for files and allocations
145 class ImdSession::Impl
149 Impl(const MDLogger &mdlog);
152 /*! \brief Prepare the socket on the MASTER. */
153 void prepareMasterSocket();
154 /*! \brief Disconnect the client. */
155 void disconnectClient();
156 /*! \brief Prints an error message and disconnects the client.
158 * Does not terminate mdrun!
160 void issueFatalError(const char *msg);
161 /*! \brief Check whether we got an incoming connection. */
163 /*! \brief Wrap imd_tryconnect in order to make it blocking.
165 * Used when the simulation should wait for an incoming connection.
168 /*! \brief Make sure that our array holding the forces received via IMD is large enough. */
169 void prepareVmdForces();
170 /*! \brief Reads forces received via IMD. */
171 void readVmdForces();
172 /*! \brief Prepares the MD force arrays. */
173 void prepareMDForces();
174 /*! \brief Copy IMD forces to MD forces.
176 * Do conversion from Cal->Joule and from
177 * Angstrom -> nm and from a pointer array to arrays to 3*N array.
179 void copyToMDForces();
180 /*! \brief Return true if any of the forces or indices changed. */
181 bool bForcesChanged() const;
182 /*! \brief Update the old_f_ind and old_forces arrays to contain the current values. */
183 void keepOldValues();
184 /*! \brief Write the applied pull forces to logfile.
186 * Call on master only!
188 void outputForces(double time);
189 /*! \brief Synchronize the nodes. */
190 void syncNodes(const t_commrec *cr, double t);
191 /*! \brief Reads header from the client and decides what to do. */
193 /*! \brief Open IMD output file and write header information.
195 * Call on master only.
197 void openOutputFile(const char *fn,
199 const gmx_output_env_t *oenv,
200 const ContinuationOptions &continuationOptions);
201 /*! \brief Creates the molecule start-end position array of molecules in the IMD group. */
202 void prepareMoleculesInImdGroup(const gmx_mtop_t *top_global);
203 /*! \brief Removes shifts of molecules diffused outside of the box. */
204 void removeMolecularShifts(const matrix box);
205 /*! \brief Initialize arrays used to assemble the positions from the other nodes. */
206 void prepareForPositionAssembly(const t_commrec *cr, const rvec x[]);
207 /*! \brief Interact with any connected VMD session */
208 bool run(int64_t step,
214 // TODO rename all the data members to have underscore suffixes
216 //! True if tpr and mdrun input combine to permit IMD sessions
217 bool sessionPossible = false;
218 //! Output file for IMD data, mainly forces.
219 FILE *outf = nullptr;
221 //! Number of atoms that can be pulled via IMD.
223 //! Part of the atoms that are local.
225 //! Global indices of the IMD atoms.
227 //! Local indices of the IMD atoms.
228 int *ind_loc = nullptr;
229 //! Allocation size for ind_loc.
231 //! Positions for all IMD atoms assembled on the master node.
233 //! Shifts for all IMD atoms, to make molecule(s) whole.
234 ivec *xa_shifts = nullptr;
235 //! Extra shifts since last DD step.
236 ivec *xa_eshifts = nullptr;
237 //! Old positions for all IMD atoms on master.
238 rvec *xa_old = nullptr;
239 //! Position of each local atom in the collective array.
240 int *xa_ind = nullptr;
242 //! Global IMD frequency, known to all ranks.
244 //! New frequency from IMD client, master only.
246 //! Default IMD frequency when disconnected.
247 int defaultNstImd = -1;
249 //! Port to use for network socket.
251 //! The IMD socket on the master node.
252 IMDSocket *socket = nullptr;
253 //! The IMD socket on the client.
254 IMDSocket *clientsocket = nullptr;
255 //! Length we got with last header.
258 //! Shall we block and wait for connection?
259 bool bWConnect = false;
260 //! Set if MD is terminated.
261 bool bTerminated = false;
262 //! Set if MD can be terminated.
263 bool bTerminatable = false;
264 //! Set if connection is present.
265 bool bConnected = false;
266 //! Set if we received new forces.
267 bool bNewForces = false;
268 //! Set if pulling from VMD is allowed.
269 bool bForceActivated = false;
271 //! Pointer to energies we send back.
272 IMDEnergyBlock *energies = nullptr;
274 //! Number of VMD forces.
275 int32_t vmd_nforces = 0;
276 //! VMD forces indices.
277 int32_t *vmd_f_ind = nullptr;
278 //! The VMD forces flat in memory.
279 float *vmd_forces = nullptr;
280 //! Number of actual MD forces; this gets communicated to the clients.
283 int *f_ind = nullptr;
284 //! The IMD pulling forces.
287 //! Buffer for force sending.
288 char *forcesendbuf = nullptr;
289 //! Buffer for coordinate sending.
290 char *coordsendbuf = nullptr;
291 //! Send buffer for energies.
292 char *energysendbuf = nullptr;
293 //! Buffer to make molecules whole before sending.
294 rvec *sendxbuf = nullptr;
296 //! Molecules block in IMD group.
299 /* The next block is used on the master node only to reduce the output
300 * without sacrificing information. If any of these values changes,
301 * we need to write output */
302 //! Old value for nforces.
304 //! Old values for force indices.
305 int *old_f_ind = nullptr;
306 //! Old values for IMD pulling forces.
307 rvec *old_forces = nullptr;
310 const MDLogger &mdlog;
311 //! Commmunication object
312 const t_commrec *cr = nullptr;
313 //! Wallcycle counting manager.
314 gmx_wallcycle *wcycle = nullptr;
315 //! Energy output handler
316 gmx_enerdata_t *enerd = nullptr;
320 * \brief Implement interactive molecular dynamics.
322 * \todo Some aspects of this module provides forces (when the user
323 * pulls on things in VMD), so in future it should have a class that
324 * models IForceProvider and is contributed to the collection of such
327 class InteractiveMolecularDynamics final : public IMDModule
330 IMdpOptionProvider *mdpOptionProvider() override { return nullptr; }
331 IMDOutputProvider *outputProvider() override { return nullptr; }
332 void initForceProviders(ForceProviders * /* forceProviders */) override {}
335 std::unique_ptr<IMDModule> createInteractiveMolecularDynamicsModule()
337 return std::make_unique<InteractiveMolecularDynamics>();
340 /*! \brief Enum for types of IMD messages.
342 * We use the same records as the NAMD/VMD IMD implementation.
344 typedef enum IMDType_t
346 IMD_DISCONNECT, /**< client disconnect */
347 IMD_ENERGIES, /**< energy data */
348 IMD_FCOORDS, /**< atomic coordinates */
349 IMD_GO, /**< start command for the simulation */
350 IMD_HANDSHAKE, /**< handshake to determine little/big endianness */
351 IMD_KILL, /**< terminates the simulation */
352 IMD_MDCOMM, /**< force data */
353 IMD_PAUSE, /**< pauses the simulation */
354 IMD_TRATE, /**< sets the IMD transmission and processing rate */
355 IMD_IOERROR, /**< I/O error */
356 IMD_NR /**< number of entries */
360 /*! \brief Names of the IMDType for error messages.
362 static const char *eIMDType_names[IMD_NR + 1] = {
378 /*! \brief Fills the header with message and the length argument. */
379 static void fill_header(IMDHeader *header, IMDMessageType type, int32_t length)
381 /* We (ab-)use htonl network function for the correct endianness */
382 header->type = imd_htonl(static_cast<int32_t>(type));
383 header->length = imd_htonl(length);
387 /*! \brief Swaps the endianess of the header. */
388 static void swap_header(IMDHeader *header)
390 /* and vice versa... */
391 header->type = imd_ntohl(header->type);
392 header->length = imd_ntohl(header->length);
396 /*! \brief Reads multiple bytes from socket. */
397 static int32_t imd_read_multiple(IMDSocket *socket, char *datptr, int32_t toread)
399 int32_t leftcount, countread;
403 /* Try to read while we haven't reached toread */
404 while (leftcount != 0)
406 if ((countread = imdsock_read(socket, datptr, leftcount)) < 0)
408 /* interrupted function call, try again... */
413 /* this is an unexpected error, return what we got */
416 return toread - leftcount;
419 /* nothing read, finished */
421 else if (countread == 0)
425 leftcount -= countread;
429 /* return nr of bytes read */
430 return toread - leftcount;
434 /*! \brief Writes multiple bytes to socket in analogy to imd_read_multiple. */
435 static int32_t imd_write_multiple(IMDSocket *socket, const char *datptr, int32_t towrite)
437 int32_t leftcount, countwritten;
441 while (leftcount != 0)
443 if ((countwritten = imdsock_write(socket, datptr, leftcount)) <= 0)
451 return towrite - leftcount;
454 leftcount -= countwritten;
455 datptr += countwritten;
458 return towrite - leftcount;
462 /*! \brief Handshake with IMD client. */
463 static int imd_handshake(IMDSocket *socket)
468 fill_header(&header, IMD_HANDSHAKE, 1);
469 header.length = c_protocolVersion; /* client wants unswapped version */
471 return static_cast<int>(imd_write_multiple(socket, reinterpret_cast<char *>(&header), c_headerSize) != c_headerSize);
475 /*! \brief Send energies using the energy block and the send buffer. */
476 static int imd_send_energies(IMDSocket *socket, const IMDEnergyBlock *energies, char *buffer)
481 recsize = c_headerSize + sizeof(IMDEnergyBlock);
482 fill_header(reinterpret_cast<IMDHeader *>(buffer), IMD_ENERGIES, 1);
483 memcpy(buffer + c_headerSize, energies, sizeof(IMDEnergyBlock));
485 return static_cast<int>(imd_write_multiple(socket, buffer, recsize) != recsize);
489 /*! \brief Receive IMD header from socket, sets the length and returns the IMD message. */
490 static IMDMessageType imd_recv_header(IMDSocket *socket, int32_t *length)
495 if (imd_read_multiple(socket, reinterpret_cast<char *>(&header), c_headerSize) != c_headerSize)
499 swap_header(&header);
500 *length = header.length;
502 return static_cast<IMDMessageType>(header.type);
506 /*! \brief Receive force indices and forces.
508 * The number of forces was previously communicated via the header.
510 static bool imd_recv_mdcomm(IMDSocket *socket, int32_t nforces, int32_t *forcendx, float *forces)
512 /* reading indices */
513 int retsize = sizeof(int32_t) * nforces;
514 int retbytes = imd_read_multiple(socket, reinterpret_cast<char *>(forcendx), retsize);
515 if (retbytes != retsize)
520 /* reading forces as float array */
521 retsize = 3 * sizeof(float) * nforces;
522 retbytes = imd_read_multiple(socket, reinterpret_cast<char *>(forces), retsize);
523 return (retbytes == retsize);
526 /* GROMACS specific functions for the IMD implementation */
527 void write_IMDgroup_to_file(bool bIMD, t_inputrec *ir, const t_state *state,
528 const gmx_mtop_t *sys, int nfile, const t_filenm fnm[])
535 IMDatoms = gmx_mtop_global_atoms(sys);
536 write_sto_conf_indexed(opt2fn("-imd", nfile, fnm), "IMDgroup", &IMDatoms,
537 state->x.rvec_array(), state->v.rvec_array(), ir->ePBC, state->box, ir->imd->nat, ir->imd->ind);
543 ImdSession::dd_make_local_IMD_atoms(const gmx_domdec_t *dd)
545 if (!impl_->sessionPossible)
550 dd_make_local_group_indices(dd->ga2la, impl_->nat, impl_->ind, &impl_->nat_loc,
551 &impl_->ind_loc, &impl_->nalloc_loc, impl_->xa_ind);
555 /*! \brief Send positions from rvec.
557 * We need a separate send buffer and conversion to Angstrom.
559 static int imd_send_rvecs(IMDSocket *socket, int nat, rvec *x, char *buffer)
564 int tuplesize = 3 * sizeof(float);
567 /* Required size for the send buffer */
568 size = c_headerSize + 3 * sizeof(float) * nat;
571 fill_header(reinterpret_cast<IMDHeader *>(buffer), IMD_FCOORDS, static_cast<int32_t>(nat));
572 for (i = 0; i < nat; i++)
574 sendx[0] = static_cast<float>(x[i][0]) * NM2A;
575 sendx[1] = static_cast<float>(x[i][1]) * NM2A;
576 sendx[2] = static_cast<float>(x[i][2]) * NM2A;
577 memcpy(buffer + c_headerSize + i * tuplesize, sendx, tuplesize);
580 return static_cast<int>(imd_write_multiple(socket, buffer, size) != size);
585 ImdSession::Impl::prepareMasterSocket()
587 if (imdsock_winsockinit() == -1)
589 gmx_fatal(FARGS, "%s Failed to initialize winsock.\n", IMDstr);
592 /* The rest is identical, first create and bind a socket and set to listen then. */
593 GMX_LOG(mdlog.warning).appendTextFormatted("%s Setting up incoming socket.", IMDstr);
594 socket = imdsock_create();
597 gmx_fatal(FARGS, "%s Failed to create socket.", IMDstr);
601 int ret = imdsock_bind(socket, port);
604 gmx_fatal(FARGS, "%s binding socket to port %d failed with error %d.\n", IMDstr, port, ret);
607 if (imd_sock_listen(socket))
609 gmx_fatal(FARGS, "%s socket listen failed with error %d.\n", IMDstr, ret);
612 if (imdsock_getport(socket, &port))
614 gmx_fatal(FARGS, "%s Could not determine port number.\n", IMDstr);
617 GMX_LOG(mdlog.warning).appendTextFormatted("%s Listening for IMD connection on port %d.", IMDstr, port);
622 ImdSession::Impl::disconnectClient()
624 /* Write out any buffered pulling data */
627 /* we first try to shut down the clientsocket */
628 imdsock_shutdown(clientsocket);
629 if (!imdsock_destroy(clientsocket))
631 GMX_LOG(mdlog.warning).appendTextFormatted("%s Failed to destroy socket.", IMDstr);
634 /* then we reset the IMD step to its default, and reset the connection boolean */
635 nstimd_new = defaultNstImd;
636 clientsocket = nullptr;
642 ImdSession::Impl::issueFatalError(const char *msg)
644 GMX_LOG(mdlog.warning).appendTextFormatted("%s %s", IMDstr, msg);
646 GMX_LOG(mdlog.warning).appendTextFormatted("%s disconnected.", IMDstr);
651 ImdSession::Impl::tryConnect()
653 if (imdsock_tryread(socket, 0, 0) > 0)
655 /* yes, we got something, accept on clientsocket */
656 clientsocket = imdsock_accept(socket);
659 GMX_LOG(mdlog.warning).appendTextFormatted("%s Accepting the connection on the socket failed.", IMDstr);
663 /* handshake with client */
664 if (imd_handshake(clientsocket))
666 issueFatalError("Connection failed.");
670 GMX_LOG(mdlog.warning).appendTextFormatted("%s Connection established, checking if I got IMD_GO orders.", IMDstr);
672 /* Check if we get the proper "GO" command from client. */
673 if (imdsock_tryread(clientsocket, c_connectWait, 0) != 1 || imd_recv_header(clientsocket, &(length)) != IMD_GO)
675 issueFatalError("No IMD_GO order received. IMD connection failed.");
689 ImdSession::Impl::blockConnect()
691 /* do not wait for connection, when e.g. ctrl+c is pressed and we will terminate anyways. */
692 if (!(static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none))
697 GMX_LOG(mdlog.warning).appendTextFormatted("%s Will wait until I have a connection and IMD_GO orders.", IMDstr);
699 /* while we have no clientsocket... 2nd part: we should still react on ctrl+c */
700 while ((!clientsocket) && (static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none))
703 imd_sleep(c_loopWait);
709 ImdSession::Impl::prepareVmdForces()
711 srenew((vmd_f_ind), vmd_nforces);
712 srenew((vmd_forces), 3*vmd_nforces);
717 ImdSession::Impl::readVmdForces()
719 /* the length of the previously received header tells us the nr of forces we will receive */
720 vmd_nforces = length;
721 /* prepare the arrays */
723 /* Now we read the forces... */
724 if (!(imd_recv_mdcomm(clientsocket, vmd_nforces, vmd_f_ind, vmd_forces)))
726 issueFatalError("Error while reading forces from remote. Disconnecting");
732 ImdSession::Impl::prepareMDForces()
734 srenew((f_ind), nforces);
735 srenew((f ), nforces);
740 ImdSession::Impl::copyToMDForces()
743 real conversion = CAL2JOULE * NM2A;
746 for (i = 0; i < nforces; i++)
748 /* Copy the indices, a copy is important because we may update the incoming forces
749 * whenever we receive new forces while the MD forces are only communicated upon
750 * IMD communication.*/
751 f_ind[i] = vmd_f_ind[i];
753 /* Convert to rvecs and do a proper unit conversion */
754 f[i][0] = vmd_forces[3*i ] * conversion;
755 f[i][1] = vmd_forces[3*i + 1] * conversion;
756 f[i][2] = vmd_forces[3*i + 2] * conversion;
761 /*! \brief Returns true if any component of the two rvecs differs. */
762 static inline bool rvecs_differ(const rvec v1, const rvec v2)
764 for (int i = 0; i < DIM; i++)
776 ImdSession::Impl::bForcesChanged() const
778 /* First, check whether the number of pulled atoms changed */
779 if (nforces != old_nforces)
784 /* Second, check whether any of the involved atoms changed */
785 for (int i = 0; i < nforces; i++)
787 if (f_ind[i] != old_f_ind[i])
793 /* Third, check whether all forces are the same */
794 for (int i = 0; i < nforces; i++)
796 if (rvecs_differ(f[i], old_forces[i]))
802 /* All old and new forces are identical! */
808 ImdSession::Impl::keepOldValues()
810 old_nforces = nforces;
812 for (int i = 0; i < nforces; i++)
814 old_f_ind[i] = f_ind[i];
815 copy_rvec(f[i], old_forces[i]);
821 ImdSession::Impl::outputForces(double time)
823 if (!bForcesChanged())
828 /* Write time and total number of applied IMD forces */
829 fprintf(outf, "%14.6e%6d", time, nforces);
831 /* Write out the global atom indices of the pulled atoms and the forces itself,
832 * write out a force only if it has changed since the last output */
833 for (int i = 0; i < nforces; i++)
835 if (rvecs_differ(f[i], old_forces[i]))
837 fprintf(outf, "%9d", ind[f_ind[i]] + 1);
838 fprintf(outf, "%12.4e%12.4e%12.4e", f[i][0], f[i][1], f[i][2]);
848 ImdSession::Impl::syncNodes(const t_commrec *cr, double t)
850 /* Notify the other nodes whether we are still connected. */
853 block_bc(cr, bConnected);
856 /* ...if not connected, the job is done here. */
862 /* Let the other nodes know whether we got a new IMD synchronization frequency. */
865 block_bc(cr, nstimd_new);
868 /* Now we all set the (new) nstimd communication time step */
871 /* We're done if we don't allow pulling at all */
872 if (!(bForceActivated))
877 /* OK, let's check if we have received forces which we need to communicate
878 * to the other nodes */
884 new_nforces = vmd_nforces;
886 /* make the "new_forces" negative, when we did not receive new ones */
889 new_nforces = vmd_nforces * -1;
893 /* make new_forces known to the clients */
896 block_bc(cr, new_nforces);
899 /* When new_natoms < 0 then we know that these are still the same forces
900 * so we don't communicate them, otherwise... */
906 /* set local VMD and nforces */
907 vmd_nforces = new_nforces;
908 nforces = vmd_nforces;
910 /* now everybody knows the number of forces in f_ind, so we can prepare
911 * the target arrays for indices and forces */
914 /* we first update the MD forces on the master by converting the VMD forces */
918 /* We also write out forces on every update, so that we know which
919 * forces are applied for every step */
926 /* In parallel mode we communicate the to-be-applied forces to the other nodes */
929 nblock_bc(cr, nforces, f_ind);
930 nblock_bc(cr, nforces, f );
933 /* done communicating the forces, reset bNewForces */
939 ImdSession::Impl::readCommand()
941 bool IMDpaused = false;
943 while (clientsocket && (imdsock_tryread(clientsocket, 0, 0) > 0 || IMDpaused))
945 IMDMessageType itype = imd_recv_header(clientsocket, &(length));
946 /* let's see what we got: */
949 /* IMD asks us to terminate the simulation, check if the user allowed this */
953 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Terminating connection and running simulation (if supported by integrator).", IMDstr);
956 gmx_set_stop_condition(gmx_stop_cond_next);
960 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Set -imdterm command line switch to allow mdrun termination from within IMD.", IMDstr);
965 /* the client doen't want to talk to us anymore */
967 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Disconnecting client.", IMDstr);
971 /* we got new forces, read them and set bNewForces flag */
977 /* the client asks us to (un)pause the simulation. So we toggle the IMDpaused state */
981 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Un-pause command received.", IMDstr);
986 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Pause command received.", IMDstr);
992 /* the client sets a new transfer rate, if we get 0, we reset the rate
993 * to the default. VMD filters 0 however */
995 nstimd_new = (length > 0) ? length : defaultNstImd;
996 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Update frequency will be set to %d.", IMDstr, nstimd_new);
999 /* Catch all rule for the remaining IMD types which we don't expect */
1001 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Received unexpected %s.", IMDstr, enum_name(static_cast<int>(itype), IMD_NR, eIMDType_names));
1002 issueFatalError("Terminating connection");
1010 ImdSession::Impl::openOutputFile(const char *fn,
1012 const gmx_output_env_t *oenv,
1013 const ContinuationOptions &continuationOptions)
1015 /* Open log file of applied IMD forces if requested */
1018 fprintf(stdout, "%s For a log of the IMD pull forces explicitly specify '-if' on the command line.\n"
1019 "%s (Not possible with energy minimization.)\n", IMDstr, IMDstr);
1023 /* If we append to an existing file, all the header information is already there */
1024 if (continuationOptions.appendFiles)
1026 outf = gmx_fio_fopen(fn, "a+");
1030 outf = gmx_fio_fopen(fn, "w+");
1031 if (nat == nat_total)
1033 fprintf(outf, "# Note that you can select an IMD index group in the .mdp file if a subset of the atoms suffices.\n");
1036 xvgr_header(outf, "IMD Pull Forces", "Time (ps)", "# of Forces / Atom IDs / Forces (kJ/mol)", exvggtNONE, oenv);
1038 fprintf(outf, "# Can display and manipulate %d (of a total of %d) atoms via IMD.\n", nat, nat_total);
1039 fprintf(outf, "# column 1 : time (ps)\n");
1040 fprintf(outf, "# column 2 : total number of atoms feeling an IMD pulling force at that time\n");
1041 fprintf(outf, "# cols. 3.-6 : global atom number of pulled atom, x-force, y-force, z-force (kJ/mol)\n");
1042 fprintf(outf, "# then follow : atom-ID, f[x], f[y], f[z] for more atoms in case the force on multiple atoms is changed simultaneously.\n");
1043 fprintf(outf, "# Note that the force on any atom is always equal to the last value for that atom-ID found in the data.\n");
1047 /* To reduce the output file size we remember the old values and output only
1048 * when something changed */
1049 snew(old_f_ind, nat); /* One can never pull on more atoms */
1050 snew(old_forces, nat);
1054 ImdSession::Impl::Impl(const MDLogger &mdlog)
1060 ImdSession::Impl::~Impl()
1064 gmx_fio_fclose(outf);
1071 ImdSession::Impl::prepareMoleculesInImdGroup(const gmx_mtop_t *top_global)
1073 /* check whether index is sorted */
1074 for (int i = 0; i < nat-1; i++)
1076 if (ind[i] > ind[i+1])
1078 gmx_fatal(FARGS, "%s IMD index is not sorted. This is currently not supported.\n", IMDstr);
1082 RangePartitioning gmols = gmx_mtop_molecules(*top_global);
1085 snew(lmols.index, gmols.numBlocks() + 1);
1088 for (int i = 0; i < gmols.numBlocks(); i++)
1090 auto mol = gmols.block(i);
1092 for (int ii = 0; ii < nat; ii++)
1094 if (mol.inRange(ind[ii]))
1101 lmols.index[lmols.nr+1] = lmols.index[lmols.nr]+count;
1105 srenew(lmols.index, lmols.nr+1);
1106 lmols.nalloc_index = lmols.nr+1;
1111 /*! \brief Copied and modified from groupcoord.c shift_positions_group(). */
1112 static void shift_positions(
1114 rvec x[], /* The positions [0..nr] */
1115 const ivec is, /* The shift [0..nr] */
1116 int nr) /* The number of positions */
1120 /* Loop over the group's atoms */
1123 for (i = 0; i < nr; i++)
1129 x[i][XX] = x[i][XX]-tx*box[XX][XX]-ty*box[YY][XX]-tz*box[ZZ][XX];
1130 x[i][YY] = x[i][YY]-ty*box[YY][YY]-tz*box[ZZ][YY];
1131 x[i][ZZ] = x[i][ZZ]-tz*box[ZZ][ZZ];
1136 for (i = 0; i < nr; i++)
1142 x[i][XX] = x[i][XX]-tx*box[XX][XX];
1143 x[i][YY] = x[i][YY]-ty*box[YY][YY];
1144 x[i][ZZ] = x[i][ZZ]-tz*box[ZZ][ZZ];
1150 void ImdSession::Impl::removeMolecularShifts(const matrix box)
1152 /* for each molecule also present in IMD group */
1153 for (int i = 0; i < mols.nr; i++)
1155 /* first we determine the minimum and maximum shifts for each molecule */
1157 ivec largest, smallest, shift;
1158 clear_ivec(largest);
1159 clear_ivec(smallest);
1162 copy_ivec(xa_shifts[mols.index[i]], largest);
1163 copy_ivec(xa_shifts[mols.index[i]], smallest);
1165 for (int ii = mols.index[i]+1; ii < mols.index[i+1]; ii++)
1167 if (xa_shifts[ii][XX] > largest[XX])
1169 largest[XX] = xa_shifts[ii][XX];
1171 if (xa_shifts[ii][XX] < smallest[XX])
1173 smallest[XX] = xa_shifts[ii][XX];
1176 if (xa_shifts[ii][YY] > largest[YY])
1178 largest[YY] = xa_shifts[ii][YY];
1180 if (xa_shifts[ii][YY] < smallest[YY])
1182 smallest[YY] = xa_shifts[ii][YY];
1185 if (xa_shifts[ii][ZZ] > largest[ZZ])
1187 largest[ZZ] = xa_shifts[ii][ZZ];
1189 if (xa_shifts[ii][ZZ] < smallest[ZZ])
1191 smallest[ZZ] = xa_shifts[ii][ZZ];
1196 /* check if we what we can subtract/add to the positions
1197 * to put them back into the central box */
1198 if (smallest[XX] > 0)
1200 shift[XX] = smallest[XX];
1202 if (smallest[YY] > 0)
1204 shift[YY] = smallest[YY];
1206 if (smallest[ZZ] > 0)
1208 shift[ZZ] = smallest[ZZ];
1211 if (largest[XX] < 0)
1213 shift[XX] = largest[XX];
1215 if (largest[YY] < 0)
1217 shift[YY] = largest[YY];
1219 if (largest[ZZ] < 0)
1221 shift[ZZ] = largest[ZZ];
1224 /* is there a shift at all? */
1225 if ((shift[XX]) || (shift[YY]) || (shift[ZZ]))
1227 int molsize = mols.index[i+1]-mols.index[i];
1228 /* shift the positions */
1229 shift_positions(box, &(xa[mols.index[i]]), shift, molsize);
1237 ImdSession::Impl::prepareForPositionAssembly(const t_commrec *cr, const rvec x[])
1241 snew(xa_shifts, nat);
1242 snew(xa_eshifts, nat);
1245 /* Save the original (whole) set of positions such that later the
1246 * molecule can always be made whole again */
1249 for (int i = 0; i < nat; i++)
1252 copy_rvec(x[ii], xa_old[i]);
1261 /* xa_ind[i] needs to be set to i for serial runs */
1262 for (int i = 0; i < nat; i++)
1268 /* Communicate initial coordinates xa_old to all processes */
1271 gmx_bcast(nat * sizeof(xa_old[0]), xa_old, cr);
1276 /*! \brief Check for non-working integrator / parallel options. */
1277 static void imd_check_integrator_parallel(const t_inputrec *ir, const t_commrec *cr)
1281 if (((ir->eI) == eiSteep) || ((ir->eI) == eiCG) || ((ir->eI) == eiLBFGS) || ((ir->eI) == eiNM))
1283 gmx_fatal(FARGS, "%s Energy minimization via steep, CG, lbfgs and nm in parallel is currently not supported by IMD.\n", IMDstr);
1288 std::unique_ptr<ImdSession>
1289 makeImdSession(const t_inputrec *ir,
1290 const t_commrec *cr,
1291 gmx_wallcycle *wcycle,
1292 gmx_enerdata_t *enerd,
1293 const gmx_multisim_t *ms,
1294 const gmx_mtop_t *top_global,
1295 const MDLogger &mdlog,
1298 const t_filenm fnm[],
1299 const gmx_output_env_t *oenv,
1300 const MdrunOptions &mdrunOptions)
1302 std::unique_ptr<ImdSession> session(new ImdSession(mdlog));
1303 auto impl = session->impl_.get();
1305 /* We will allow IMD sessions only if supported by the binary and
1306 explicitly enabled in the .tpr file */
1307 if (!GMX_IMD || !ir->bIMD)
1312 // TODO As IMD is intended for interactivity, and the .tpr file
1313 // opted in for that, it is acceptable to write more terminal
1314 // output than in a typical simulation. However, all the GMX_LOG
1315 // statements below should go to both the log file and to the
1316 // terminal. This is probably be implemented by adding a logging
1317 // stream named like ImdInfo, to separate it from warning and to
1318 // send it to both destinations.
1320 if (EI_DYNAMICS(ir->eI))
1322 impl->defaultNstImd = ir->nstcalcenergy;
1324 else if (EI_ENERGY_MINIMIZATION(ir->eI))
1326 impl->defaultNstImd = 1;
1330 GMX_LOG(mdlog.warning).appendTextFormatted("%s Integrator '%s' is not supported for Interactive Molecular Dynamics, running normally instead", IMDstr, ei_names[ir->eI]);
1335 GMX_LOG(mdlog.warning).appendTextFormatted("%s Cannot use IMD for multiple simulations or replica exchange, running normally instead", IMDstr);
1339 const auto &options = mdrunOptions.imdOptions;
1341 bool createSession = false;
1342 /* It seems we have a .tpr file that defines an IMD group and thus allows IMD connections.
1343 * Check whether we can actually provide the IMD functionality for this setting: */
1346 /* Check whether IMD was enabled by one of the command line switches: */
1347 if (options.wait || options.terminatable || options.pull)
1349 GMX_LOG(mdlog.warning).appendTextFormatted("%s Enabled. This simulation will accept incoming IMD connections.", IMDstr);
1350 createSession = true;
1354 GMX_LOG(mdlog.warning).appendTextFormatted("%s None of the -imd switches was used.\n"
1355 "%s This run will not accept incoming IMD connections", IMDstr, IMDstr);
1357 } /* end master only */
1359 /* Let the other nodes know whether we want IMD */
1362 block_bc(cr, createSession);
1365 /*... if not we are done.*/
1372 /* check if we're using a sane integrator / parallel combination */
1373 imd_check_integrator_parallel(ir, cr);
1377 *****************************************************
1378 * From here on we assume that IMD is turned on *
1379 *****************************************************
1382 int nat_total = top_global->natoms;
1384 /* Initialize IMD session. If we read in a pre-IMD .tpr file, ir->imd->nat
1385 * will be zero. For those cases we transfer _all_ atomic positions */
1386 impl->sessionPossible = true;
1387 impl->nat = ir->imd->nat > 0 ? ir->imd->nat : nat_total;
1388 if (options.port >= 1)
1390 impl->port = options.port;
1393 impl->wcycle = wcycle;
1394 impl->enerd = enerd;
1396 /* We might need to open an output file for IMD forces data */
1399 impl->openOutputFile(opt2fn("-if", nfile, fnm), nat_total, oenv, mdrunOptions.continuationOptions);
1402 /* Make sure that we operate with a valid atom index array for the IMD atoms */
1403 if (ir->imd->nat > 0)
1405 /* Point to the user-supplied array of atom numbers */
1406 impl->ind = ir->imd->ind;
1410 /* Make a dummy (ind[i] = i) array of all atoms */
1411 snew(impl->ind, nat_total);
1412 for (int i = 0; i < nat_total; i++)
1418 /* read environment on master and prepare socket for incoming connections */
1421 /* we allocate memory for our IMD energy structure */
1422 int32_t recsize = c_headerSize + sizeof(IMDEnergyBlock);
1423 snew(impl->energysendbuf, recsize);
1425 /* Shall we wait for a connection? */
1428 impl->bWConnect = true;
1429 GMX_LOG(mdlog.warning).appendTextFormatted("%s Pausing simulation while no IMD connection present (-imdwait).", IMDstr);
1432 /* Will the IMD clients be able to terminate the simulation? */
1433 if (options.terminatable)
1435 impl->bTerminatable = true;
1436 GMX_LOG(mdlog.warning).appendTextFormatted("%s Allow termination of the simulation from IMD client (-imdterm).", IMDstr);
1439 /* Is pulling from IMD client allowed? */
1442 impl->bForceActivated = true;
1443 GMX_LOG(mdlog.warning).appendTextFormatted("%s Pulling from IMD remote is enabled (-imdpull).", IMDstr);
1446 /* Initialize send buffers with constant size */
1447 snew(impl->sendxbuf, impl->nat);
1448 snew(impl->energies, 1);
1449 int32_t bufxsize = c_headerSize + 3 * sizeof(float) * impl->nat;
1450 snew(impl->coordsendbuf, bufxsize);
1453 /* do we allow interactive pulling? If so let the other nodes know. */
1456 block_bc(cr, impl->bForceActivated);
1459 /* setup the listening socket on master process */
1462 GMX_LOG(mdlog.warning).appendTextFormatted("%s Setting port for connection requests to %d.", IMDstr, impl->port);
1463 impl->prepareMasterSocket();
1464 /* Wait until we have a connection if specified before */
1465 if (impl->bWConnect)
1467 impl->blockConnect();
1471 GMX_LOG(mdlog.warning).appendTextFormatted("%s -imdwait not set, starting simulation.", IMDstr);
1474 /* Let the other nodes know whether we are connected */
1475 impl->syncNodes(cr, 0);
1477 /* Initialize arrays used to assemble the positions from the other nodes */
1478 impl->prepareForPositionAssembly(cr, x);
1480 /* Initialize molecule blocks to make them whole later...*/
1483 impl->prepareMoleculesInImdGroup(top_global);
1490 bool ImdSession::Impl::run(int64_t step,
1497 if (!sessionPossible)
1502 wallcycle_start(wcycle, ewcIMD);
1504 /* read command from client and check if new incoming connection */
1507 /* If not already connected, check for new connections */
1520 /* Let's see if we have new IMD messages for us */
1527 /* is this an IMD communication step? */
1528 bool imdstep = do_per_step(step, nstimd);
1530 /* OK so this is an IMD step ... */
1533 /* First we sync all nodes to let everybody know whether we are connected to VMD */
1537 /* If a client is connected, we collect the positions
1538 * and put molecules back into the box before transfer */
1539 if ((imdstep && bConnected)
1540 || bNS) /* independent of imdstep, we communicate positions at each NS step */
1542 /* Transfer the IMD positions to the master node. Every node contributes
1543 * its local positions x and stores them in the assembled xa array. */
1544 communicate_group_positions(cr, xa, xa_shifts, xa_eshifts,
1545 true, x, nat, nat_loc,
1546 ind_loc, xa_ind, xa_old, box);
1548 /* If connected and master -> remove shifts */
1549 if ((imdstep && bConnected) && MASTER(cr))
1551 removeMolecularShifts(box);
1555 wallcycle_stop(wcycle, ewcIMD);
1560 bool ImdSession::run(int64_t step,
1566 return impl_->run(step, bNS, box, x, t);
1569 void ImdSession::fillEnergyRecord(int64_t step,
1570 bool bHaveNewEnergies)
1572 if (!impl_->sessionPossible || !impl_->clientsocket)
1577 IMDEnergyBlock *ene = impl_->energies;
1581 /* In MPI-parallel simulations the energies are not accessible a at every time step.
1582 * We update them if we have new values, otherwise, the energy values from the
1583 * last global communication step are still on display in the viewer. */
1584 if (bHaveNewEnergies)
1586 ene->T_abs = static_cast<float>(impl_->enerd->term[F_TEMP ]);
1587 ene->E_pot = static_cast<float>(impl_->enerd->term[F_EPOT ]);
1588 ene->E_tot = static_cast<float>(impl_->enerd->term[F_ETOT ]);
1589 ene->E_bond = static_cast<float>(impl_->enerd->term[F_BONDS ]);
1590 ene->E_angle = static_cast<float>(impl_->enerd->term[F_ANGLES ]);
1591 ene->E_dihe = static_cast<float>(impl_->enerd->term[F_PDIHS ]);
1592 ene->E_impr = static_cast<float>(impl_->enerd->term[F_IDIHS ]);
1593 ene->E_vdw = static_cast<float>(impl_->enerd->term[F_LJ ]);
1594 ene->E_coul = static_cast<float>(impl_->enerd->term[F_COUL_SR]);
1600 ImdSession::sendPositionsAndEnergies()
1602 if (!impl_->sessionPossible || !impl_->clientsocket)
1607 if (imd_send_energies(impl_->clientsocket, impl_->energies, impl_->energysendbuf))
1609 impl_->issueFatalError("Error sending updated energies. Disconnecting client.");
1612 if (imd_send_rvecs(impl_->clientsocket, impl_->nat, impl_->xa, impl_->coordsendbuf))
1614 impl_->issueFatalError("Error sending updated positions. Disconnecting client.");
1620 ImdSession::updateEnergyRecordAndSendPositionsAndEnergies(bool bIMDstep,
1622 bool bHaveNewEnergies)
1624 if (!impl_->sessionPossible)
1629 wallcycle_start(impl_->wcycle, ewcIMD);
1631 /* Update time step for IMD and prepare IMD energy record if we have new energies. */
1632 fillEnergyRecord(step, bHaveNewEnergies);
1636 /* Send positions and energies to VMD client via IMD */
1637 sendPositionsAndEnergies();
1640 wallcycle_stop(impl_->wcycle, ewcIMD);
1643 void ImdSession::applyForces(rvec *f)
1645 if (!impl_->sessionPossible || !impl_->bForceActivated)
1650 wallcycle_start(impl_->wcycle, ewcIMD);
1652 for (int i = 0; i < impl_->nforces; i++)
1654 /* j are the indices in the "System group".*/
1655 int j = impl_->ind[impl_->f_ind[i]];
1657 /* check if this is a local atom and find out locndx */
1659 const t_commrec *cr = impl_->cr;
1660 if (PAR(cr) && (locndx = cr->dd->ga2la->findHome(j)))
1665 rvec_inc(f[j], impl_->f[i]);
1668 wallcycle_stop(impl_->wcycle, ewcIMD);
1671 ImdSession::ImdSession(const MDLogger &mdlog)
1672 : impl_(new Impl(mdlog))
1675 ImdSession::~ImdSession() = default;