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/handlerestart.h"
72 #include "gromacs/mdrunutility/multisim.h"
73 #include "gromacs/mdtypes/enerdata.h"
74 #include "gromacs/mdtypes/imdmodule.h"
75 #include "gromacs/mdtypes/inputrec.h"
76 #include "gromacs/mdtypes/md_enums.h"
77 #include "gromacs/mdtypes/mdrunoptions.h"
78 #include "gromacs/mdtypes/state.h"
79 #include "gromacs/pbcutil/pbc.h"
80 #include "gromacs/timing/wallcycle.h"
81 #include "gromacs/topology/mtop_util.h"
82 #include "gromacs/topology/topology.h"
83 #include "gromacs/utility/fatalerror.h"
84 #include "gromacs/utility/logger.h"
85 #include "gromacs/utility/smalloc.h"
86 #include "gromacs/utility/stringutil.h"
91 /*! \brief How long shall we wait in seconds until we check for a connection again? */
92 constexpr int c_loopWait = 1;
94 /*! \brief How long shall we check for the IMD_GO? */
95 constexpr int c_connectWait = 1;
97 /*! \brief IMD Header Size. */
98 constexpr int c_headerSize = 8;
100 /*! \brief IMD Protocol Version. */
101 constexpr int c_protocolVersion = 2;
105 * IMD (interactive molecular dynamics) energy record.
107 * As in the original IMD implementation. Energies in kcal/mol.
108 * NOTE: We return the energies in GROMACS / SI units,
109 * so they also show up as SI in VMD.
114 int32_t tstep; /**< time step */
115 float T_abs; /**< absolute temperature */
116 float E_tot; /**< total energy */
117 float E_pot; /**< potential energy */
118 float E_vdw; /**< van der Waals energy */
119 float E_coul; /**< Coulomb interaction energy */
120 float E_bond; /**< bonds energy */
121 float E_angle; /**< angles energy */
122 float E_dihe; /**< dihedrals energy */
123 float E_impr; /**< improper dihedrals energy */
128 * \brief IMD (interactive molecular dynamics) communication structure.
130 * This structure defines the IMD communication message header & protocol version.
134 int32_t type; /**< Type of IMD message, see IMDType_t above */
135 int32_t length; /**< Length */
140 * \brief Implementation type for the IMD session
142 * \todo Make this class (or one extracted from it) model
144 * \todo Use RAII for files and allocations
146 class ImdSession::Impl
150 Impl(const MDLogger& mdlog);
153 /*! \brief Prepare the socket on the MASTER. */
154 void prepareMasterSocket();
155 /*! \brief Disconnect the client. */
156 void disconnectClient();
157 /*! \brief Prints an error message and disconnects the client.
159 * Does not terminate mdrun!
161 void issueFatalError(const char* msg);
162 /*! \brief Check whether we got an incoming connection. */
164 /*! \brief Wrap imd_tryconnect in order to make it blocking.
166 * Used when the simulation should wait for an incoming connection.
169 /*! \brief Make sure that our array holding the forces received via IMD is large enough. */
170 void prepareVmdForces();
171 /*! \brief Reads forces received via IMD. */
172 void readVmdForces();
173 /*! \brief Prepares the MD force arrays. */
174 void prepareMDForces();
175 /*! \brief Copy IMD forces to MD forces.
177 * Do conversion from Cal->Joule and from
178 * Angstrom -> nm and from a pointer array to arrays to 3*N array.
180 void copyToMDForces();
181 /*! \brief Return true if any of the forces or indices changed. */
182 bool bForcesChanged() const;
183 /*! \brief Update the old_f_ind and old_forces arrays to contain the current values. */
184 void keepOldValues();
185 /*! \brief Write the applied pull forces to logfile.
187 * Call on master only!
189 void outputForces(double time);
190 /*! \brief Synchronize the nodes. */
191 void syncNodes(const t_commrec* cr, double t);
192 /*! \brief Reads header from the client and decides what to do. */
194 /*! \brief Open IMD output file and write header information.
196 * Call on master only.
198 void openOutputFile(const char* fn, int nat_total, const gmx_output_env_t* oenv, StartingBehavior startingBehavior);
199 /*! \brief Creates the molecule start-end position array of molecules in the IMD group. */
200 void prepareMoleculesInImdGroup(const gmx_mtop_t* top_global);
201 /*! \brief Removes shifts of molecules diffused outside of the box. */
202 void removeMolecularShifts(const matrix box);
203 /*! \brief Initialize arrays used to assemble the positions from the other nodes. */
204 void prepareForPositionAssembly(const t_commrec* cr, const rvec x[]);
205 /*! \brief Interact with any connected VMD session */
206 bool run(int64_t step, bool bNS, const matrix box, const rvec x[], double t);
208 // TODO rename all the data members to have underscore suffixes
210 //! True if tpr and mdrun input combine to permit IMD sessions
211 bool sessionPossible = false;
212 //! Output file for IMD data, mainly forces.
213 FILE* outf = nullptr;
215 //! Number of atoms that can be pulled via IMD.
217 //! Part of the atoms that are local.
219 //! Global indices of the IMD atoms.
221 //! Local indices of the IMD atoms.
222 int* ind_loc = nullptr;
223 //! Allocation size for ind_loc.
225 //! Positions for all IMD atoms assembled on the master node.
227 //! Shifts for all IMD atoms, to make molecule(s) whole.
228 ivec* xa_shifts = nullptr;
229 //! Extra shifts since last DD step.
230 ivec* xa_eshifts = nullptr;
231 //! Old positions for all IMD atoms on master.
232 rvec* xa_old = nullptr;
233 //! Position of each local atom in the collective array.
234 int* xa_ind = nullptr;
236 //! Global IMD frequency, known to all ranks.
238 //! New frequency from IMD client, master only.
240 //! Default IMD frequency when disconnected.
241 int defaultNstImd = -1;
243 //! Port to use for network socket.
245 //! The IMD socket on the master node.
246 IMDSocket* socket = nullptr;
247 //! The IMD socket on the client.
248 IMDSocket* clientsocket = nullptr;
249 //! Length we got with last header.
252 //! Shall we block and wait for connection?
253 bool bWConnect = false;
254 //! Set if MD is terminated.
255 bool bTerminated = false;
256 //! Set if MD can be terminated.
257 bool bTerminatable = false;
258 //! Set if connection is present.
259 bool bConnected = false;
260 //! Set if we received new forces.
261 bool bNewForces = false;
262 //! Set if pulling from VMD is allowed.
263 bool bForceActivated = false;
265 //! Pointer to energies we send back.
266 IMDEnergyBlock* energies = nullptr;
268 //! Number of VMD forces.
269 int32_t vmd_nforces = 0;
270 //! VMD forces indices.
271 int32_t* vmd_f_ind = nullptr;
272 //! The VMD forces flat in memory.
273 float* vmd_forces = nullptr;
274 //! Number of actual MD forces; this gets communicated to the clients.
277 int* f_ind = nullptr;
278 //! The IMD pulling forces.
281 //! Buffer for force sending.
282 char* forcesendbuf = nullptr;
283 //! Buffer for coordinate sending.
284 char* coordsendbuf = nullptr;
285 //! Send buffer for energies.
286 char* energysendbuf = nullptr;
287 //! Buffer to make molecules whole before sending.
288 rvec* sendxbuf = nullptr;
290 //! Molecules block in IMD group.
293 /* The next block is used on the master node only to reduce the output
294 * without sacrificing information. If any of these values changes,
295 * we need to write output */
296 //! Old value for nforces.
298 //! Old values for force indices.
299 int* old_f_ind = nullptr;
300 //! Old values for IMD pulling forces.
301 rvec* old_forces = nullptr;
304 const MDLogger& mdlog;
305 //! Commmunication object
306 const t_commrec* cr = nullptr;
307 //! Wallcycle counting manager.
308 gmx_wallcycle* wcycle = nullptr;
309 //! Energy output handler
310 gmx_enerdata_t* enerd = nullptr;
314 * \brief Implement interactive molecular dynamics.
316 * \todo Some aspects of this module provides forces (when the user
317 * pulls on things in VMD), so in future it should have a class that
318 * models IForceProvider and is contributed to the collection of such
321 class InteractiveMolecularDynamics final : public IMDModule
324 IMdpOptionProvider* mdpOptionProvider() override { return nullptr; }
325 IMDOutputProvider* outputProvider() override { return nullptr; }
326 void initForceProviders(ForceProviders* /* forceProviders */) override {}
329 std::unique_ptr<IMDModule> createInteractiveMolecularDynamicsModule()
331 return std::make_unique<InteractiveMolecularDynamics>();
334 /*! \brief Enum for types of IMD messages.
336 * We use the same records as the NAMD/VMD IMD implementation.
338 typedef enum IMDType_t
340 IMD_DISCONNECT, /**< client disconnect */
341 IMD_ENERGIES, /**< energy data */
342 IMD_FCOORDS, /**< atomic coordinates */
343 IMD_GO, /**< start command for the simulation */
344 IMD_HANDSHAKE, /**< handshake to determine little/big endianness */
345 IMD_KILL, /**< terminates the simulation */
346 IMD_MDCOMM, /**< force data */
347 IMD_PAUSE, /**< pauses the simulation */
348 IMD_TRATE, /**< sets the IMD transmission and processing rate */
349 IMD_IOERROR, /**< I/O error */
350 IMD_NR /**< number of entries */
354 /*! \brief Names of the IMDType for error messages.
356 static const char* eIMDType_names[IMD_NR + 1] = { "IMD_DISCONNECT", "IMD_ENERGIES", "IMD_FCOORDS",
357 "IMD_GO", "IMD_HANDSHAKE", "IMD_KILL",
358 "IMD_MDCOMM", "IMD_PAUSE", "IMD_TRATE",
359 "IMD_IOERROR", nullptr };
362 /*! \brief Fills the header with message and the length argument. */
363 static void fill_header(IMDHeader* header, IMDMessageType type, int32_t length)
365 /* We (ab-)use htonl network function for the correct endianness */
366 header->type = imd_htonl(static_cast<int32_t>(type));
367 header->length = imd_htonl(length);
371 /*! \brief Swaps the endianess of the header. */
372 static void swap_header(IMDHeader* header)
374 /* and vice versa... */
375 header->type = imd_ntohl(header->type);
376 header->length = imd_ntohl(header->length);
380 /*! \brief Reads multiple bytes from socket. */
381 static int32_t imd_read_multiple(IMDSocket* socket, char* datptr, int32_t toread)
383 int32_t leftcount, countread;
387 /* Try to read while we haven't reached toread */
388 while (leftcount != 0)
390 if ((countread = imdsock_read(socket, datptr, leftcount)) < 0)
392 /* interrupted function call, try again... */
397 /* this is an unexpected error, return what we got */
400 return toread - leftcount;
403 /* nothing read, finished */
405 else if (countread == 0)
409 leftcount -= countread;
413 /* return nr of bytes read */
414 return toread - leftcount;
418 /*! \brief Writes multiple bytes to socket in analogy to imd_read_multiple. */
419 static int32_t imd_write_multiple(IMDSocket* socket, const char* datptr, int32_t towrite)
421 int32_t leftcount, countwritten;
425 while (leftcount != 0)
427 if ((countwritten = imdsock_write(socket, datptr, leftcount)) <= 0)
435 return towrite - leftcount;
438 leftcount -= countwritten;
439 datptr += countwritten;
442 return towrite - leftcount;
446 /*! \brief Handshake with IMD client. */
447 static int imd_handshake(IMDSocket* socket)
452 fill_header(&header, IMD_HANDSHAKE, 1);
453 header.length = c_protocolVersion; /* client wants unswapped version */
455 return static_cast<int>(imd_write_multiple(socket, reinterpret_cast<char*>(&header), c_headerSize)
460 /*! \brief Send energies using the energy block and the send buffer. */
461 static int imd_send_energies(IMDSocket* socket, const IMDEnergyBlock* energies, char* buffer)
466 recsize = c_headerSize + sizeof(IMDEnergyBlock);
467 fill_header(reinterpret_cast<IMDHeader*>(buffer), IMD_ENERGIES, 1);
468 memcpy(buffer + c_headerSize, energies, sizeof(IMDEnergyBlock));
470 return static_cast<int>(imd_write_multiple(socket, buffer, recsize) != recsize);
474 /*! \brief Receive IMD header from socket, sets the length and returns the IMD message. */
475 static IMDMessageType imd_recv_header(IMDSocket* socket, int32_t* length)
480 if (imd_read_multiple(socket, reinterpret_cast<char*>(&header), c_headerSize) != c_headerSize)
484 swap_header(&header);
485 *length = header.length;
487 return static_cast<IMDMessageType>(header.type);
491 /*! \brief Receive force indices and forces.
493 * The number of forces was previously communicated via the header.
495 static bool imd_recv_mdcomm(IMDSocket* socket, int32_t nforces, int32_t* forcendx, float* forces)
497 /* reading indices */
498 int retsize = sizeof(int32_t) * nforces;
499 int retbytes = imd_read_multiple(socket, reinterpret_cast<char*>(forcendx), retsize);
500 if (retbytes != retsize)
505 /* reading forces as float array */
506 retsize = 3 * sizeof(float) * nforces;
507 retbytes = imd_read_multiple(socket, reinterpret_cast<char*>(forces), retsize);
508 return (retbytes == retsize);
511 /* GROMACS specific functions for the IMD implementation */
512 void write_IMDgroup_to_file(bool bIMD,
514 const t_state* state,
515 const gmx_mtop_t* sys,
517 const t_filenm fnm[])
524 IMDatoms = gmx_mtop_global_atoms(sys);
525 write_sto_conf_indexed(opt2fn("-imd", nfile, fnm), "IMDgroup", &IMDatoms, state->x.rvec_array(),
526 state->v.rvec_array(), ir->ePBC, state->box, ir->imd->nat, ir->imd->ind);
531 void ImdSession::dd_make_local_IMD_atoms(const gmx_domdec_t* dd)
533 if (!impl_->sessionPossible)
538 dd_make_local_group_indices(dd->ga2la, impl_->nat, impl_->ind, &impl_->nat_loc, &impl_->ind_loc,
539 &impl_->nalloc_loc, impl_->xa_ind);
543 /*! \brief Send positions from rvec.
545 * We need a separate send buffer and conversion to Angstrom.
547 static int imd_send_rvecs(IMDSocket* socket, int nat, rvec* x, char* buffer)
552 int tuplesize = 3 * sizeof(float);
555 /* Required size for the send buffer */
556 size = c_headerSize + 3 * sizeof(float) * nat;
559 fill_header(reinterpret_cast<IMDHeader*>(buffer), IMD_FCOORDS, static_cast<int32_t>(nat));
560 for (i = 0; i < nat; i++)
562 sendx[0] = static_cast<float>(x[i][0]) * NM2A;
563 sendx[1] = static_cast<float>(x[i][1]) * NM2A;
564 sendx[2] = static_cast<float>(x[i][2]) * NM2A;
565 memcpy(buffer + c_headerSize + i * tuplesize, sendx, tuplesize);
568 return static_cast<int>(imd_write_multiple(socket, buffer, size) != size);
572 void ImdSession::Impl::prepareMasterSocket()
574 if (imdsock_winsockinit() == -1)
576 gmx_fatal(FARGS, "%s Failed to initialize winsock.\n", IMDstr);
579 /* The rest is identical, first create and bind a socket and set to listen then. */
580 GMX_LOG(mdlog.warning).appendTextFormatted("%s Setting up incoming socket.", IMDstr);
581 socket = imdsock_create();
584 gmx_fatal(FARGS, "%s Failed to create socket.", IMDstr);
588 int ret = imdsock_bind(socket, port);
591 gmx_fatal(FARGS, "%s binding socket to port %d failed with error %d.\n", IMDstr, port, ret);
594 if (imd_sock_listen(socket))
596 gmx_fatal(FARGS, "%s socket listen failed with error %d.\n", IMDstr, ret);
599 if (imdsock_getport(socket, &port))
601 gmx_fatal(FARGS, "%s Could not determine port number.\n", IMDstr);
604 GMX_LOG(mdlog.warning).appendTextFormatted("%s Listening for IMD connection on port %d.", IMDstr, port);
608 void ImdSession::Impl::disconnectClient()
610 /* Write out any buffered pulling data */
613 /* we first try to shut down the clientsocket */
614 imdsock_shutdown(clientsocket);
615 if (!imdsock_destroy(clientsocket))
617 GMX_LOG(mdlog.warning).appendTextFormatted("%s Failed to destroy socket.", IMDstr);
620 /* then we reset the IMD step to its default, and reset the connection boolean */
621 nstimd_new = defaultNstImd;
622 clientsocket = nullptr;
627 void ImdSession::Impl::issueFatalError(const char* msg)
629 GMX_LOG(mdlog.warning).appendTextFormatted("%s %s", IMDstr, msg);
631 GMX_LOG(mdlog.warning).appendTextFormatted("%s disconnected.", IMDstr);
635 bool ImdSession::Impl::tryConnect()
637 if (imdsock_tryread(socket, 0, 0) > 0)
639 /* yes, we got something, accept on clientsocket */
640 clientsocket = imdsock_accept(socket);
643 GMX_LOG(mdlog.warning)
644 .appendTextFormatted("%s Accepting the connection on the socket failed.", IMDstr);
648 /* handshake with client */
649 if (imd_handshake(clientsocket))
651 issueFatalError("Connection failed.");
655 GMX_LOG(mdlog.warning)
656 .appendTextFormatted("%s Connection established, checking if I got IMD_GO orders.", IMDstr);
658 /* Check if we get the proper "GO" command from client. */
659 if (imdsock_tryread(clientsocket, c_connectWait, 0) != 1
660 || imd_recv_header(clientsocket, &(length)) != IMD_GO)
662 issueFatalError("No IMD_GO order received. IMD connection failed.");
675 void ImdSession::Impl::blockConnect()
677 /* do not wait for connection, when e.g. ctrl+c is pressed and we will terminate anyways. */
678 if (!(static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none))
683 GMX_LOG(mdlog.warning)
684 .appendTextFormatted("%s Will wait until I have a connection and IMD_GO orders.", IMDstr);
686 /* while we have no clientsocket... 2nd part: we should still react on ctrl+c */
687 while ((!clientsocket) && (static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none))
690 imd_sleep(c_loopWait);
695 void ImdSession::Impl::prepareVmdForces()
697 srenew((vmd_f_ind), vmd_nforces);
698 srenew((vmd_forces), 3 * vmd_nforces);
702 void ImdSession::Impl::readVmdForces()
704 /* the length of the previously received header tells us the nr of forces we will receive */
705 vmd_nforces = length;
706 /* prepare the arrays */
708 /* Now we read the forces... */
709 if (!(imd_recv_mdcomm(clientsocket, vmd_nforces, vmd_f_ind, vmd_forces)))
711 issueFatalError("Error while reading forces from remote. Disconnecting");
716 void ImdSession::Impl::prepareMDForces()
718 srenew((f_ind), nforces);
719 srenew((f), nforces);
723 void ImdSession::Impl::copyToMDForces()
726 real conversion = CAL2JOULE * NM2A;
729 for (i = 0; i < nforces; i++)
731 /* Copy the indices, a copy is important because we may update the incoming forces
732 * whenever we receive new forces while the MD forces are only communicated upon
733 * IMD communication.*/
734 f_ind[i] = vmd_f_ind[i];
736 /* Convert to rvecs and do a proper unit conversion */
737 f[i][0] = vmd_forces[3 * i] * conversion;
738 f[i][1] = vmd_forces[3 * i + 1] * conversion;
739 f[i][2] = vmd_forces[3 * i + 2] * conversion;
744 /*! \brief Returns true if any component of the two rvecs differs. */
745 static inline bool rvecs_differ(const rvec v1, const rvec v2)
747 for (int i = 0; i < DIM; i++)
758 bool ImdSession::Impl::bForcesChanged() const
760 /* First, check whether the number of pulled atoms changed */
761 if (nforces != old_nforces)
766 /* Second, check whether any of the involved atoms changed */
767 for (int i = 0; i < nforces; i++)
769 if (f_ind[i] != old_f_ind[i])
775 /* Third, check whether all forces are the same */
776 for (int i = 0; i < nforces; i++)
778 if (rvecs_differ(f[i], old_forces[i]))
784 /* All old and new forces are identical! */
789 void ImdSession::Impl::keepOldValues()
791 old_nforces = nforces;
793 for (int i = 0; i < nforces; i++)
795 old_f_ind[i] = f_ind[i];
796 copy_rvec(f[i], old_forces[i]);
801 void ImdSession::Impl::outputForces(double time)
803 if (!bForcesChanged())
808 /* Write time and total number of applied IMD forces */
809 fprintf(outf, "%14.6e%6d", time, nforces);
811 /* Write out the global atom indices of the pulled atoms and the forces itself,
812 * write out a force only if it has changed since the last output */
813 for (int i = 0; i < nforces; i++)
815 if (rvecs_differ(f[i], old_forces[i]))
817 fprintf(outf, "%9d", ind[f_ind[i]] + 1);
818 fprintf(outf, "%12.4e%12.4e%12.4e", f[i][0], f[i][1], f[i][2]);
827 void ImdSession::Impl::syncNodes(const t_commrec* cr, double t)
829 /* Notify the other nodes whether we are still connected. */
832 block_bc(cr, bConnected);
835 /* ...if not connected, the job is done here. */
841 /* Let the other nodes know whether we got a new IMD synchronization frequency. */
844 block_bc(cr, nstimd_new);
847 /* Now we all set the (new) nstimd communication time step */
850 /* We're done if we don't allow pulling at all */
851 if (!(bForceActivated))
856 /* OK, let's check if we have received forces which we need to communicate
857 * to the other nodes */
863 new_nforces = vmd_nforces;
865 /* make the "new_forces" negative, when we did not receive new ones */
868 new_nforces = vmd_nforces * -1;
872 /* make new_forces known to the clients */
875 block_bc(cr, new_nforces);
878 /* When new_natoms < 0 then we know that these are still the same forces
879 * so we don't communicate them, otherwise... */
885 /* set local VMD and nforces */
886 vmd_nforces = new_nforces;
887 nforces = vmd_nforces;
889 /* now everybody knows the number of forces in f_ind, so we can prepare
890 * the target arrays for indices and forces */
893 /* we first update the MD forces on the master by converting the VMD forces */
897 /* We also write out forces on every update, so that we know which
898 * forces are applied for every step */
905 /* In parallel mode we communicate the to-be-applied forces to the other nodes */
908 nblock_bc(cr, nforces, f_ind);
909 nblock_bc(cr, nforces, f);
912 /* done communicating the forces, reset bNewForces */
917 void ImdSession::Impl::readCommand()
919 bool IMDpaused = false;
921 while (clientsocket && (imdsock_tryread(clientsocket, 0, 0) > 0 || IMDpaused))
923 IMDMessageType itype = imd_recv_header(clientsocket, &(length));
924 /* let's see what we got: */
927 /* IMD asks us to terminate the simulation, check if the user allowed this */
931 GMX_LOG(mdlog.warning)
932 .appendTextFormatted(
933 " %s Terminating connection and running simulation (if "
934 "supported by integrator).",
938 gmx_set_stop_condition(gmx_stop_cond_next);
942 GMX_LOG(mdlog.warning)
943 .appendTextFormatted(
944 " %s Set -imdterm command line switch to allow mdrun "
945 "termination from within IMD.",
951 /* the client doen't want to talk to us anymore */
953 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Disconnecting client.", IMDstr);
957 /* we got new forces, read them and set bNewForces flag */
963 /* the client asks us to (un)pause the simulation. So we toggle the IMDpaused state */
967 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Un-pause command received.", IMDstr);
972 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Pause command received.", IMDstr);
978 /* the client sets a new transfer rate, if we get 0, we reset the rate
979 * to the default. VMD filters 0 however */
981 nstimd_new = (length > 0) ? length : defaultNstImd;
982 GMX_LOG(mdlog.warning)
983 .appendTextFormatted(" %s Update frequency will be set to %d.", IMDstr, nstimd_new);
986 /* Catch all rule for the remaining IMD types which we don't expect */
988 GMX_LOG(mdlog.warning)
989 .appendTextFormatted(" %s Received unexpected %s.", IMDstr,
990 enum_name(static_cast<int>(itype), IMD_NR, eIMDType_names));
991 issueFatalError("Terminating connection");
998 void ImdSession::Impl::openOutputFile(const char* fn,
1000 const gmx_output_env_t* oenv,
1001 const gmx::StartingBehavior startingBehavior)
1003 /* Open log file of applied IMD forces if requested */
1007 "%s For a log of the IMD pull forces explicitly specify '-if' on the command "
1009 "%s (Not possible with energy minimization.)\n",
1014 /* If we append to an existing file, all the header information is already there */
1015 if (startingBehavior == StartingBehavior::RestartWithAppending)
1017 outf = gmx_fio_fopen(fn, "a+");
1021 outf = gmx_fio_fopen(fn, "w+");
1022 if (nat == nat_total)
1025 "# Note that you can select an IMD index group in the .mdp file if a subset of "
1026 "the atoms suffices.\n");
1029 xvgr_header(outf, "IMD Pull Forces", "Time (ps)",
1030 "# of Forces / Atom IDs / Forces (kJ/mol)", exvggtNONE, oenv);
1032 fprintf(outf, "# Can display and manipulate %d (of a total of %d) atoms via IMD.\n", nat, nat_total);
1033 fprintf(outf, "# column 1 : time (ps)\n");
1035 "# column 2 : total number of atoms feeling an IMD pulling force at that "
1038 "# cols. 3.-6 : global atom number of pulled atom, x-force, y-force, z-force "
1041 "# then follow : atom-ID, f[x], f[y], f[z] for more atoms in case the force on "
1042 "multiple atoms is changed simultaneously.\n");
1044 "# Note that the force on any atom is always equal to the last value for that "
1045 "atom-ID found in the data.\n");
1049 /* To reduce the output file size we remember the old values and output only
1050 * when something changed */
1051 snew(old_f_ind, nat); /* One can never pull on more atoms */
1052 snew(old_forces, nat);
1056 ImdSession::Impl::Impl(const MDLogger& mdlog) : mdlog(mdlog)
1061 ImdSession::Impl::~Impl()
1065 gmx_fio_fclose(outf);
1071 void 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.isInRange(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(const matrix box,
1113 rvec x[], /* The positions [0..nr] */
1114 const ivec is, /* The shift [0..nr] */
1115 int nr) /* The number of positions */
1119 /* Loop over the group's atoms */
1122 for (i = 0; i < nr; i++)
1128 x[i][XX] = x[i][XX] - tx * box[XX][XX] - ty * box[YY][XX] - tz * box[ZZ][XX];
1129 x[i][YY] = x[i][YY] - ty * box[YY][YY] - tz * box[ZZ][YY];
1130 x[i][ZZ] = x[i][ZZ] - tz * box[ZZ][ZZ];
1135 for (i = 0; i < nr; i++)
1141 x[i][XX] = x[i][XX] - tx * box[XX][XX];
1142 x[i][YY] = x[i][YY] - ty * box[YY][YY];
1143 x[i][ZZ] = x[i][ZZ] - tz * box[ZZ][ZZ];
1149 void ImdSession::Impl::removeMolecularShifts(const matrix box)
1151 /* for each molecule also present in IMD group */
1152 for (int i = 0; i < mols.nr; i++)
1154 /* first we determine the minimum and maximum shifts for each molecule */
1156 ivec largest, smallest, shift;
1157 clear_ivec(largest);
1158 clear_ivec(smallest);
1161 copy_ivec(xa_shifts[mols.index[i]], largest);
1162 copy_ivec(xa_shifts[mols.index[i]], smallest);
1164 for (int ii = mols.index[i] + 1; ii < mols.index[i + 1]; ii++)
1166 if (xa_shifts[ii][XX] > largest[XX])
1168 largest[XX] = xa_shifts[ii][XX];
1170 if (xa_shifts[ii][XX] < smallest[XX])
1172 smallest[XX] = xa_shifts[ii][XX];
1175 if (xa_shifts[ii][YY] > largest[YY])
1177 largest[YY] = xa_shifts[ii][YY];
1179 if (xa_shifts[ii][YY] < smallest[YY])
1181 smallest[YY] = xa_shifts[ii][YY];
1184 if (xa_shifts[ii][ZZ] > largest[ZZ])
1186 largest[ZZ] = xa_shifts[ii][ZZ];
1188 if (xa_shifts[ii][ZZ] < smallest[ZZ])
1190 smallest[ZZ] = xa_shifts[ii][ZZ];
1194 /* check if we what we can subtract/add to the positions
1195 * to put them back into the central box */
1196 if (smallest[XX] > 0)
1198 shift[XX] = smallest[XX];
1200 if (smallest[YY] > 0)
1202 shift[YY] = smallest[YY];
1204 if (smallest[ZZ] > 0)
1206 shift[ZZ] = smallest[ZZ];
1209 if (largest[XX] < 0)
1211 shift[XX] = largest[XX];
1213 if (largest[YY] < 0)
1215 shift[YY] = largest[YY];
1217 if (largest[ZZ] < 0)
1219 shift[ZZ] = largest[ZZ];
1222 /* is there a shift at all? */
1223 if ((shift[XX]) || (shift[YY]) || (shift[ZZ]))
1225 int molsize = mols.index[i + 1] - mols.index[i];
1226 /* shift the positions */
1227 shift_positions(box, &(xa[mols.index[i]]), shift, molsize);
1233 void ImdSession::Impl::prepareForPositionAssembly(const t_commrec* cr, const rvec x[])
1237 snew(xa_shifts, nat);
1238 snew(xa_eshifts, nat);
1241 /* Save the original (whole) set of positions such that later the
1242 * molecule can always be made whole again */
1245 for (int i = 0; i < nat; i++)
1248 copy_rvec(x[ii], xa_old[i]);
1257 /* xa_ind[i] needs to be set to i for serial runs */
1258 for (int i = 0; i < nat; i++)
1264 /* Communicate initial coordinates xa_old to all processes */
1267 gmx_bcast(nat * sizeof(xa_old[0]), xa_old, cr);
1272 /*! \brief Check for non-working integrator / parallel options. */
1273 static void imd_check_integrator_parallel(const t_inputrec* ir, const t_commrec* cr)
1277 if (((ir->eI) == eiSteep) || ((ir->eI) == eiCG) || ((ir->eI) == eiLBFGS) || ((ir->eI) == eiNM))
1280 "%s Energy minimization via steep, CG, lbfgs and nm in parallel is currently "
1281 "not supported by IMD.\n",
1287 std::unique_ptr<ImdSession> makeImdSession(const t_inputrec* ir,
1288 const t_commrec* cr,
1289 gmx_wallcycle* wcycle,
1290 gmx_enerdata_t* enerd,
1291 const gmx_multisim_t* ms,
1292 const gmx_mtop_t* top_global,
1293 const MDLogger& mdlog,
1296 const t_filenm fnm[],
1297 const gmx_output_env_t* oenv,
1298 const ImdOptions& options,
1299 const gmx::StartingBehavior startingBehavior)
1301 std::unique_ptr<ImdSession> session(new ImdSession(mdlog));
1302 auto impl = session->impl_.get();
1304 /* We will allow IMD sessions only if supported by the binary and
1305 explicitly enabled in the .tpr file */
1306 if (!GMX_IMD || !ir->bIMD)
1311 // TODO As IMD is intended for interactivity, and the .tpr file
1312 // opted in for that, it is acceptable to write more terminal
1313 // output than in a typical simulation. However, all the GMX_LOG
1314 // statements below should go to both the log file and to the
1315 // terminal. This is probably be implemented by adding a logging
1316 // stream named like ImdInfo, to separate it from warning and to
1317 // send it to both destinations.
1319 if (EI_DYNAMICS(ir->eI))
1321 impl->defaultNstImd = ir->nstcalcenergy;
1323 else if (EI_ENERGY_MINIMIZATION(ir->eI))
1325 impl->defaultNstImd = 1;
1329 GMX_LOG(mdlog.warning)
1330 .appendTextFormatted(
1331 "%s Integrator '%s' is not supported for Interactive Molecular Dynamics, "
1332 "running normally instead",
1333 IMDstr, ei_names[ir->eI]);
1338 GMX_LOG(mdlog.warning)
1339 .appendTextFormatted(
1340 "%s Cannot use IMD for multiple simulations or replica exchange, running "
1346 bool createSession = false;
1347 /* It seems we have a .tpr file that defines an IMD group and thus allows IMD connections.
1348 * Check whether we can actually provide the IMD functionality for this setting: */
1351 /* Check whether IMD was enabled by one of the command line switches: */
1352 if (options.wait || options.terminatable || options.pull)
1354 GMX_LOG(mdlog.warning)
1355 .appendTextFormatted(
1356 "%s Enabled. This simulation will accept incoming IMD connections.", IMDstr);
1357 createSession = true;
1361 GMX_LOG(mdlog.warning)
1362 .appendTextFormatted(
1363 "%s None of the -imd switches was used.\n"
1364 "%s This run will not accept incoming IMD connections",
1367 } /* end master only */
1369 /* Let the other nodes know whether we want IMD */
1372 block_bc(cr, createSession);
1375 /*... if not we are done.*/
1382 /* check if we're using a sane integrator / parallel combination */
1383 imd_check_integrator_parallel(ir, cr);
1387 *****************************************************
1388 * From here on we assume that IMD is turned on *
1389 *****************************************************
1392 int nat_total = top_global->natoms;
1394 /* Initialize IMD session. If we read in a pre-IMD .tpr file, ir->imd->nat
1395 * will be zero. For those cases we transfer _all_ atomic positions */
1396 impl->sessionPossible = true;
1397 impl->nat = ir->imd->nat > 0 ? ir->imd->nat : nat_total;
1398 if (options.port >= 1)
1400 impl->port = options.port;
1403 impl->wcycle = wcycle;
1404 impl->enerd = enerd;
1406 /* We might need to open an output file for IMD forces data */
1409 impl->openOutputFile(opt2fn("-if", nfile, fnm), nat_total, oenv, startingBehavior);
1412 /* Make sure that we operate with a valid atom index array for the IMD atoms */
1413 if (ir->imd->nat > 0)
1415 /* Point to the user-supplied array of atom numbers */
1416 impl->ind = ir->imd->ind;
1420 /* Make a dummy (ind[i] = i) array of all atoms */
1421 snew(impl->ind, nat_total);
1422 for (int i = 0; i < nat_total; i++)
1428 /* read environment on master and prepare socket for incoming connections */
1431 /* we allocate memory for our IMD energy structure */
1432 int32_t recsize = c_headerSize + sizeof(IMDEnergyBlock);
1433 snew(impl->energysendbuf, recsize);
1435 /* Shall we wait for a connection? */
1438 impl->bWConnect = true;
1439 GMX_LOG(mdlog.warning)
1440 .appendTextFormatted(
1441 "%s Pausing simulation while no IMD connection present (-imdwait).", IMDstr);
1444 /* Will the IMD clients be able to terminate the simulation? */
1445 if (options.terminatable)
1447 impl->bTerminatable = true;
1448 GMX_LOG(mdlog.warning)
1449 .appendTextFormatted(
1450 "%s Allow termination of the simulation from IMD client (-imdterm).", IMDstr);
1453 /* Is pulling from IMD client allowed? */
1456 impl->bForceActivated = true;
1457 GMX_LOG(mdlog.warning)
1458 .appendTextFormatted("%s Pulling from IMD remote is enabled (-imdpull).", IMDstr);
1461 /* Initialize send buffers with constant size */
1462 snew(impl->sendxbuf, impl->nat);
1463 snew(impl->energies, 1);
1464 int32_t bufxsize = c_headerSize + 3 * sizeof(float) * impl->nat;
1465 snew(impl->coordsendbuf, bufxsize);
1468 /* do we allow interactive pulling? If so let the other nodes know. */
1471 block_bc(cr, impl->bForceActivated);
1474 /* setup the listening socket on master process */
1477 GMX_LOG(mdlog.warning).appendTextFormatted("%s Setting port for connection requests to %d.", IMDstr, impl->port);
1478 impl->prepareMasterSocket();
1479 /* Wait until we have a connection if specified before */
1480 if (impl->bWConnect)
1482 impl->blockConnect();
1486 GMX_LOG(mdlog.warning).appendTextFormatted("%s -imdwait not set, starting simulation.", IMDstr);
1489 /* Let the other nodes know whether we are connected */
1490 impl->syncNodes(cr, 0);
1492 /* Initialize arrays used to assemble the positions from the other nodes */
1493 impl->prepareForPositionAssembly(cr, x);
1495 /* Initialize molecule blocks to make them whole later...*/
1498 impl->prepareMoleculesInImdGroup(top_global);
1505 bool ImdSession::Impl::run(int64_t step, bool bNS, const matrix box, const rvec x[], double t)
1508 if (!sessionPossible)
1513 wallcycle_start(wcycle, ewcIMD);
1515 /* read command from client and check if new incoming connection */
1518 /* If not already connected, check for new connections */
1531 /* Let's see if we have new IMD messages for us */
1538 /* is this an IMD communication step? */
1539 bool imdstep = do_per_step(step, nstimd);
1541 /* OK so this is an IMD step ... */
1544 /* First we sync all nodes to let everybody know whether we are connected to VMD */
1548 /* If a client is connected, we collect the positions
1549 * and put molecules back into the box before transfer */
1550 if ((imdstep && bConnected) || bNS) /* independent of imdstep, we communicate positions at each NS step */
1552 /* Transfer the IMD positions to the master node. Every node contributes
1553 * its local positions x and stores them in the assembled xa array. */
1554 communicate_group_positions(cr, xa, xa_shifts, xa_eshifts, true, x, nat, nat_loc, ind_loc,
1555 xa_ind, xa_old, box);
1557 /* If connected and master -> remove shifts */
1558 if ((imdstep && bConnected) && MASTER(cr))
1560 removeMolecularShifts(box);
1564 wallcycle_stop(wcycle, ewcIMD);
1569 bool ImdSession::run(int64_t step, bool bNS, const matrix box, const rvec x[], double t)
1571 return impl_->run(step, bNS, box, x, t);
1574 void ImdSession::fillEnergyRecord(int64_t step, bool bHaveNewEnergies)
1576 if (!impl_->sessionPossible || !impl_->clientsocket)
1581 IMDEnergyBlock* ene = impl_->energies;
1585 /* In MPI-parallel simulations the energies are not accessible a at every time step.
1586 * We update them if we have new values, otherwise, the energy values from the
1587 * last global communication step are still on display in the viewer. */
1588 if (bHaveNewEnergies)
1590 ene->T_abs = static_cast<float>(impl_->enerd->term[F_TEMP]);
1591 ene->E_pot = static_cast<float>(impl_->enerd->term[F_EPOT]);
1592 ene->E_tot = static_cast<float>(impl_->enerd->term[F_ETOT]);
1593 ene->E_bond = static_cast<float>(impl_->enerd->term[F_BONDS]);
1594 ene->E_angle = static_cast<float>(impl_->enerd->term[F_ANGLES]);
1595 ene->E_dihe = static_cast<float>(impl_->enerd->term[F_PDIHS]);
1596 ene->E_impr = static_cast<float>(impl_->enerd->term[F_IDIHS]);
1597 ene->E_vdw = static_cast<float>(impl_->enerd->term[F_LJ]);
1598 ene->E_coul = static_cast<float>(impl_->enerd->term[F_COUL_SR]);
1603 void ImdSession::sendPositionsAndEnergies()
1605 if (!impl_->sessionPossible || !impl_->clientsocket)
1610 if (imd_send_energies(impl_->clientsocket, impl_->energies, impl_->energysendbuf))
1612 impl_->issueFatalError("Error sending updated energies. Disconnecting client.");
1615 if (imd_send_rvecs(impl_->clientsocket, impl_->nat, impl_->xa, impl_->coordsendbuf))
1617 impl_->issueFatalError("Error sending updated positions. Disconnecting client.");
1622 void ImdSession::updateEnergyRecordAndSendPositionsAndEnergies(bool bIMDstep, int64_t step, 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) : impl_(new Impl(mdlog)) {}
1673 ImdSession::~ImdSession() = default;