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
,
200 const gmx_output_env_t
*oenv
,
201 StartingBehavior startingBehavior
);
202 /*! \brief Creates the molecule start-end position array of molecules in the IMD group. */
203 void prepareMoleculesInImdGroup(const gmx_mtop_t
*top_global
);
204 /*! \brief Removes shifts of molecules diffused outside of the box. */
205 void removeMolecularShifts(const matrix box
);
206 /*! \brief Initialize arrays used to assemble the positions from the other nodes. */
207 void prepareForPositionAssembly(const t_commrec
*cr
, const rvec x
[]);
208 /*! \brief Interact with any connected VMD session */
209 bool run(int64_t step
,
215 // TODO rename all the data members to have underscore suffixes
217 //! True if tpr and mdrun input combine to permit IMD sessions
218 bool sessionPossible
= false;
219 //! Output file for IMD data, mainly forces.
220 FILE *outf
= nullptr;
222 //! Number of atoms that can be pulled via IMD.
224 //! Part of the atoms that are local.
226 //! Global indices of the IMD atoms.
228 //! Local indices of the IMD atoms.
229 int *ind_loc
= nullptr;
230 //! Allocation size for ind_loc.
232 //! Positions for all IMD atoms assembled on the master node.
234 //! Shifts for all IMD atoms, to make molecule(s) whole.
235 ivec
*xa_shifts
= nullptr;
236 //! Extra shifts since last DD step.
237 ivec
*xa_eshifts
= nullptr;
238 //! Old positions for all IMD atoms on master.
239 rvec
*xa_old
= nullptr;
240 //! Position of each local atom in the collective array.
241 int *xa_ind
= nullptr;
243 //! Global IMD frequency, known to all ranks.
245 //! New frequency from IMD client, master only.
247 //! Default IMD frequency when disconnected.
248 int defaultNstImd
= -1;
250 //! Port to use for network socket.
252 //! The IMD socket on the master node.
253 IMDSocket
*socket
= nullptr;
254 //! The IMD socket on the client.
255 IMDSocket
*clientsocket
= nullptr;
256 //! Length we got with last header.
259 //! Shall we block and wait for connection?
260 bool bWConnect
= false;
261 //! Set if MD is terminated.
262 bool bTerminated
= false;
263 //! Set if MD can be terminated.
264 bool bTerminatable
= false;
265 //! Set if connection is present.
266 bool bConnected
= false;
267 //! Set if we received new forces.
268 bool bNewForces
= false;
269 //! Set if pulling from VMD is allowed.
270 bool bForceActivated
= false;
272 //! Pointer to energies we send back.
273 IMDEnergyBlock
*energies
= nullptr;
275 //! Number of VMD forces.
276 int32_t vmd_nforces
= 0;
277 //! VMD forces indices.
278 int32_t *vmd_f_ind
= nullptr;
279 //! The VMD forces flat in memory.
280 float *vmd_forces
= nullptr;
281 //! Number of actual MD forces; this gets communicated to the clients.
284 int *f_ind
= nullptr;
285 //! The IMD pulling forces.
288 //! Buffer for force sending.
289 char *forcesendbuf
= nullptr;
290 //! Buffer for coordinate sending.
291 char *coordsendbuf
= nullptr;
292 //! Send buffer for energies.
293 char *energysendbuf
= nullptr;
294 //! Buffer to make molecules whole before sending.
295 rvec
*sendxbuf
= nullptr;
297 //! Molecules block in IMD group.
300 /* The next block is used on the master node only to reduce the output
301 * without sacrificing information. If any of these values changes,
302 * we need to write output */
303 //! Old value for nforces.
305 //! Old values for force indices.
306 int *old_f_ind
= nullptr;
307 //! Old values for IMD pulling forces.
308 rvec
*old_forces
= nullptr;
311 const MDLogger
&mdlog
;
312 //! Commmunication object
313 const t_commrec
*cr
= nullptr;
314 //! Wallcycle counting manager.
315 gmx_wallcycle
*wcycle
= nullptr;
316 //! Energy output handler
317 gmx_enerdata_t
*enerd
= nullptr;
321 * \brief Implement interactive molecular dynamics.
323 * \todo Some aspects of this module provides forces (when the user
324 * pulls on things in VMD), so in future it should have a class that
325 * models IForceProvider and is contributed to the collection of such
328 class InteractiveMolecularDynamics final
: public IMDModule
331 IMdpOptionProvider
*mdpOptionProvider() override
{ return nullptr; }
332 IMDOutputProvider
*outputProvider() override
{ return nullptr; }
333 void initForceProviders(ForceProviders
* /* forceProviders */) override
{}
336 std::unique_ptr
<IMDModule
> createInteractiveMolecularDynamicsModule()
338 return std::make_unique
<InteractiveMolecularDynamics
>();
341 /*! \brief Enum for types of IMD messages.
343 * We use the same records as the NAMD/VMD IMD implementation.
345 typedef enum IMDType_t
347 IMD_DISCONNECT
, /**< client disconnect */
348 IMD_ENERGIES
, /**< energy data */
349 IMD_FCOORDS
, /**< atomic coordinates */
350 IMD_GO
, /**< start command for the simulation */
351 IMD_HANDSHAKE
, /**< handshake to determine little/big endianness */
352 IMD_KILL
, /**< terminates the simulation */
353 IMD_MDCOMM
, /**< force data */
354 IMD_PAUSE
, /**< pauses the simulation */
355 IMD_TRATE
, /**< sets the IMD transmission and processing rate */
356 IMD_IOERROR
, /**< I/O error */
357 IMD_NR
/**< number of entries */
361 /*! \brief Names of the IMDType for error messages.
363 static const char *eIMDType_names
[IMD_NR
+ 1] = {
379 /*! \brief Fills the header with message and the length argument. */
380 static void fill_header(IMDHeader
*header
, IMDMessageType type
, int32_t length
)
382 /* We (ab-)use htonl network function for the correct endianness */
383 header
->type
= imd_htonl(static_cast<int32_t>(type
));
384 header
->length
= imd_htonl(length
);
388 /*! \brief Swaps the endianess of the header. */
389 static void swap_header(IMDHeader
*header
)
391 /* and vice versa... */
392 header
->type
= imd_ntohl(header
->type
);
393 header
->length
= imd_ntohl(header
->length
);
397 /*! \brief Reads multiple bytes from socket. */
398 static int32_t imd_read_multiple(IMDSocket
*socket
, char *datptr
, int32_t toread
)
400 int32_t leftcount
, countread
;
404 /* Try to read while we haven't reached toread */
405 while (leftcount
!= 0)
407 if ((countread
= imdsock_read(socket
, datptr
, leftcount
)) < 0)
409 /* interrupted function call, try again... */
414 /* this is an unexpected error, return what we got */
417 return toread
- leftcount
;
420 /* nothing read, finished */
422 else if (countread
== 0)
426 leftcount
-= countread
;
430 /* return nr of bytes read */
431 return toread
- leftcount
;
435 /*! \brief Writes multiple bytes to socket in analogy to imd_read_multiple. */
436 static int32_t imd_write_multiple(IMDSocket
*socket
, const char *datptr
, int32_t towrite
)
438 int32_t leftcount
, countwritten
;
442 while (leftcount
!= 0)
444 if ((countwritten
= imdsock_write(socket
, datptr
, leftcount
)) <= 0)
452 return towrite
- leftcount
;
455 leftcount
-= countwritten
;
456 datptr
+= countwritten
;
459 return towrite
- leftcount
;
463 /*! \brief Handshake with IMD client. */
464 static int imd_handshake(IMDSocket
*socket
)
469 fill_header(&header
, IMD_HANDSHAKE
, 1);
470 header
.length
= c_protocolVersion
; /* client wants unswapped version */
472 return static_cast<int>(imd_write_multiple(socket
, reinterpret_cast<char *>(&header
), c_headerSize
) != c_headerSize
);
476 /*! \brief Send energies using the energy block and the send buffer. */
477 static int imd_send_energies(IMDSocket
*socket
, const IMDEnergyBlock
*energies
, char *buffer
)
482 recsize
= c_headerSize
+ sizeof(IMDEnergyBlock
);
483 fill_header(reinterpret_cast<IMDHeader
*>(buffer
), IMD_ENERGIES
, 1);
484 memcpy(buffer
+ c_headerSize
, energies
, sizeof(IMDEnergyBlock
));
486 return static_cast<int>(imd_write_multiple(socket
, buffer
, recsize
) != recsize
);
490 /*! \brief Receive IMD header from socket, sets the length and returns the IMD message. */
491 static IMDMessageType
imd_recv_header(IMDSocket
*socket
, int32_t *length
)
496 if (imd_read_multiple(socket
, reinterpret_cast<char *>(&header
), c_headerSize
) != c_headerSize
)
500 swap_header(&header
);
501 *length
= header
.length
;
503 return static_cast<IMDMessageType
>(header
.type
);
507 /*! \brief Receive force indices and forces.
509 * The number of forces was previously communicated via the header.
511 static bool imd_recv_mdcomm(IMDSocket
*socket
, int32_t nforces
, int32_t *forcendx
, float *forces
)
513 /* reading indices */
514 int retsize
= sizeof(int32_t) * nforces
;
515 int retbytes
= imd_read_multiple(socket
, reinterpret_cast<char *>(forcendx
), retsize
);
516 if (retbytes
!= retsize
)
521 /* reading forces as float array */
522 retsize
= 3 * sizeof(float) * nforces
;
523 retbytes
= imd_read_multiple(socket
, reinterpret_cast<char *>(forces
), retsize
);
524 return (retbytes
== retsize
);
527 /* GROMACS specific functions for the IMD implementation */
528 void write_IMDgroup_to_file(bool bIMD
, t_inputrec
*ir
, const t_state
*state
,
529 const gmx_mtop_t
*sys
, int nfile
, const t_filenm fnm
[])
536 IMDatoms
= gmx_mtop_global_atoms(sys
);
537 write_sto_conf_indexed(opt2fn("-imd", nfile
, fnm
), "IMDgroup", &IMDatoms
,
538 state
->x
.rvec_array(), state
->v
.rvec_array(), ir
->ePBC
, state
->box
, ir
->imd
->nat
, ir
->imd
->ind
);
544 ImdSession::dd_make_local_IMD_atoms(const gmx_domdec_t
*dd
)
546 if (!impl_
->sessionPossible
)
551 dd_make_local_group_indices(dd
->ga2la
, impl_
->nat
, impl_
->ind
, &impl_
->nat_loc
,
552 &impl_
->ind_loc
, &impl_
->nalloc_loc
, impl_
->xa_ind
);
556 /*! \brief Send positions from rvec.
558 * We need a separate send buffer and conversion to Angstrom.
560 static int imd_send_rvecs(IMDSocket
*socket
, int nat
, rvec
*x
, char *buffer
)
565 int tuplesize
= 3 * sizeof(float);
568 /* Required size for the send buffer */
569 size
= c_headerSize
+ 3 * sizeof(float) * nat
;
572 fill_header(reinterpret_cast<IMDHeader
*>(buffer
), IMD_FCOORDS
, static_cast<int32_t>(nat
));
573 for (i
= 0; i
< nat
; i
++)
575 sendx
[0] = static_cast<float>(x
[i
][0]) * NM2A
;
576 sendx
[1] = static_cast<float>(x
[i
][1]) * NM2A
;
577 sendx
[2] = static_cast<float>(x
[i
][2]) * NM2A
;
578 memcpy(buffer
+ c_headerSize
+ i
* tuplesize
, sendx
, tuplesize
);
581 return static_cast<int>(imd_write_multiple(socket
, buffer
, size
) != size
);
586 ImdSession::Impl::prepareMasterSocket()
588 if (imdsock_winsockinit() == -1)
590 gmx_fatal(FARGS
, "%s Failed to initialize winsock.\n", IMDstr
);
593 /* The rest is identical, first create and bind a socket and set to listen then. */
594 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Setting up incoming socket.", IMDstr
);
595 socket
= imdsock_create();
598 gmx_fatal(FARGS
, "%s Failed to create socket.", IMDstr
);
602 int ret
= imdsock_bind(socket
, port
);
605 gmx_fatal(FARGS
, "%s binding socket to port %d failed with error %d.\n", IMDstr
, port
, ret
);
608 if (imd_sock_listen(socket
))
610 gmx_fatal(FARGS
, "%s socket listen failed with error %d.\n", IMDstr
, ret
);
613 if (imdsock_getport(socket
, &port
))
615 gmx_fatal(FARGS
, "%s Could not determine port number.\n", IMDstr
);
618 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Listening for IMD connection on port %d.", IMDstr
, port
);
623 ImdSession::Impl::disconnectClient()
625 /* Write out any buffered pulling data */
628 /* we first try to shut down the clientsocket */
629 imdsock_shutdown(clientsocket
);
630 if (!imdsock_destroy(clientsocket
))
632 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Failed to destroy socket.", IMDstr
);
635 /* then we reset the IMD step to its default, and reset the connection boolean */
636 nstimd_new
= defaultNstImd
;
637 clientsocket
= nullptr;
643 ImdSession::Impl::issueFatalError(const char *msg
)
645 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s %s", IMDstr
, msg
);
647 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s disconnected.", IMDstr
);
652 ImdSession::Impl::tryConnect()
654 if (imdsock_tryread(socket
, 0, 0) > 0)
656 /* yes, we got something, accept on clientsocket */
657 clientsocket
= imdsock_accept(socket
);
660 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Accepting the connection on the socket failed.", IMDstr
);
664 /* handshake with client */
665 if (imd_handshake(clientsocket
))
667 issueFatalError("Connection failed.");
671 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Connection established, checking if I got IMD_GO orders.", IMDstr
);
673 /* Check if we get the proper "GO" command from client. */
674 if (imdsock_tryread(clientsocket
, c_connectWait
, 0) != 1 || imd_recv_header(clientsocket
, &(length
)) != IMD_GO
)
676 issueFatalError("No IMD_GO order received. IMD connection failed.");
690 ImdSession::Impl::blockConnect()
692 /* do not wait for connection, when e.g. ctrl+c is pressed and we will terminate anyways. */
693 if (!(static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none
))
698 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Will wait until I have a connection and IMD_GO orders.", IMDstr
);
700 /* while we have no clientsocket... 2nd part: we should still react on ctrl+c */
701 while ((!clientsocket
) && (static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none
))
704 imd_sleep(c_loopWait
);
710 ImdSession::Impl::prepareVmdForces()
712 srenew((vmd_f_ind
), vmd_nforces
);
713 srenew((vmd_forces
), 3*vmd_nforces
);
718 ImdSession::Impl::readVmdForces()
720 /* the length of the previously received header tells us the nr of forces we will receive */
721 vmd_nforces
= length
;
722 /* prepare the arrays */
724 /* Now we read the forces... */
725 if (!(imd_recv_mdcomm(clientsocket
, vmd_nforces
, vmd_f_ind
, vmd_forces
)))
727 issueFatalError("Error while reading forces from remote. Disconnecting");
733 ImdSession::Impl::prepareMDForces()
735 srenew((f_ind
), nforces
);
736 srenew((f
), nforces
);
741 ImdSession::Impl::copyToMDForces()
744 real conversion
= CAL2JOULE
* NM2A
;
747 for (i
= 0; i
< nforces
; i
++)
749 /* Copy the indices, a copy is important because we may update the incoming forces
750 * whenever we receive new forces while the MD forces are only communicated upon
751 * IMD communication.*/
752 f_ind
[i
] = vmd_f_ind
[i
];
754 /* Convert to rvecs and do a proper unit conversion */
755 f
[i
][0] = vmd_forces
[3*i
] * conversion
;
756 f
[i
][1] = vmd_forces
[3*i
+ 1] * conversion
;
757 f
[i
][2] = vmd_forces
[3*i
+ 2] * conversion
;
762 /*! \brief Returns true if any component of the two rvecs differs. */
763 static inline bool rvecs_differ(const rvec v1
, const rvec v2
)
765 for (int i
= 0; i
< DIM
; i
++)
777 ImdSession::Impl::bForcesChanged() const
779 /* First, check whether the number of pulled atoms changed */
780 if (nforces
!= old_nforces
)
785 /* Second, check whether any of the involved atoms changed */
786 for (int i
= 0; i
< nforces
; i
++)
788 if (f_ind
[i
] != old_f_ind
[i
])
794 /* Third, check whether all forces are the same */
795 for (int i
= 0; i
< nforces
; i
++)
797 if (rvecs_differ(f
[i
], old_forces
[i
]))
803 /* All old and new forces are identical! */
809 ImdSession::Impl::keepOldValues()
811 old_nforces
= nforces
;
813 for (int i
= 0; i
< nforces
; i
++)
815 old_f_ind
[i
] = f_ind
[i
];
816 copy_rvec(f
[i
], old_forces
[i
]);
822 ImdSession::Impl::outputForces(double time
)
824 if (!bForcesChanged())
829 /* Write time and total number of applied IMD forces */
830 fprintf(outf
, "%14.6e%6d", time
, nforces
);
832 /* Write out the global atom indices of the pulled atoms and the forces itself,
833 * write out a force only if it has changed since the last output */
834 for (int i
= 0; i
< nforces
; i
++)
836 if (rvecs_differ(f
[i
], old_forces
[i
]))
838 fprintf(outf
, "%9d", ind
[f_ind
[i
]] + 1);
839 fprintf(outf
, "%12.4e%12.4e%12.4e", f
[i
][0], f
[i
][1], f
[i
][2]);
849 ImdSession::Impl::syncNodes(const t_commrec
*cr
, double t
)
851 /* Notify the other nodes whether we are still connected. */
854 block_bc(cr
, bConnected
);
857 /* ...if not connected, the job is done here. */
863 /* Let the other nodes know whether we got a new IMD synchronization frequency. */
866 block_bc(cr
, nstimd_new
);
869 /* Now we all set the (new) nstimd communication time step */
872 /* We're done if we don't allow pulling at all */
873 if (!(bForceActivated
))
878 /* OK, let's check if we have received forces which we need to communicate
879 * to the other nodes */
885 new_nforces
= vmd_nforces
;
887 /* make the "new_forces" negative, when we did not receive new ones */
890 new_nforces
= vmd_nforces
* -1;
894 /* make new_forces known to the clients */
897 block_bc(cr
, new_nforces
);
900 /* When new_natoms < 0 then we know that these are still the same forces
901 * so we don't communicate them, otherwise... */
907 /* set local VMD and nforces */
908 vmd_nforces
= new_nforces
;
909 nforces
= vmd_nforces
;
911 /* now everybody knows the number of forces in f_ind, so we can prepare
912 * the target arrays for indices and forces */
915 /* we first update the MD forces on the master by converting the VMD forces */
919 /* We also write out forces on every update, so that we know which
920 * forces are applied for every step */
927 /* In parallel mode we communicate the to-be-applied forces to the other nodes */
930 nblock_bc(cr
, nforces
, f_ind
);
931 nblock_bc(cr
, nforces
, f
);
934 /* done communicating the forces, reset bNewForces */
940 ImdSession::Impl::readCommand()
942 bool IMDpaused
= false;
944 while (clientsocket
&& (imdsock_tryread(clientsocket
, 0, 0) > 0 || IMDpaused
))
946 IMDMessageType itype
= imd_recv_header(clientsocket
, &(length
));
947 /* let's see what we got: */
950 /* IMD asks us to terminate the simulation, check if the user allowed this */
954 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Terminating connection and running simulation (if supported by integrator).", IMDstr
);
957 gmx_set_stop_condition(gmx_stop_cond_next
);
961 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Set -imdterm command line switch to allow mdrun termination from within IMD.", IMDstr
);
966 /* the client doen't want to talk to us anymore */
968 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Disconnecting client.", IMDstr
);
972 /* we got new forces, read them and set bNewForces flag */
978 /* the client asks us to (un)pause the simulation. So we toggle the IMDpaused state */
982 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Un-pause command received.", IMDstr
);
987 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Pause command received.", IMDstr
);
993 /* the client sets a new transfer rate, if we get 0, we reset the rate
994 * to the default. VMD filters 0 however */
996 nstimd_new
= (length
> 0) ? length
: defaultNstImd
;
997 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Update frequency will be set to %d.", IMDstr
, nstimd_new
);
1000 /* Catch all rule for the remaining IMD types which we don't expect */
1002 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Received unexpected %s.", IMDstr
, enum_name(static_cast<int>(itype
), IMD_NR
, eIMDType_names
));
1003 issueFatalError("Terminating connection");
1011 ImdSession::Impl::openOutputFile(const char *fn
,
1013 const gmx_output_env_t
*oenv
,
1014 const gmx::StartingBehavior startingBehavior
)
1016 /* Open log file of applied IMD forces if requested */
1019 fprintf(stdout
, "%s For a log of the IMD pull forces explicitly specify '-if' on the command line.\n"
1020 "%s (Not possible with energy minimization.)\n", IMDstr
, IMDstr
);
1024 /* If we append to an existing file, all the header information is already there */
1025 if (startingBehavior
== StartingBehavior::RestartWithAppending
)
1027 outf
= gmx_fio_fopen(fn
, "a+");
1031 outf
= gmx_fio_fopen(fn
, "w+");
1032 if (nat
== nat_total
)
1034 fprintf(outf
, "# Note that you can select an IMD index group in the .mdp file if a subset of the atoms suffices.\n");
1037 xvgr_header(outf
, "IMD Pull Forces", "Time (ps)", "# of Forces / Atom IDs / Forces (kJ/mol)", exvggtNONE
, oenv
);
1039 fprintf(outf
, "# Can display and manipulate %d (of a total of %d) atoms via IMD.\n", nat
, nat_total
);
1040 fprintf(outf
, "# column 1 : time (ps)\n");
1041 fprintf(outf
, "# column 2 : total number of atoms feeling an IMD pulling force at that time\n");
1042 fprintf(outf
, "# cols. 3.-6 : global atom number of pulled atom, x-force, y-force, z-force (kJ/mol)\n");
1043 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");
1044 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");
1048 /* To reduce the output file size we remember the old values and output only
1049 * when something changed */
1050 snew(old_f_ind
, nat
); /* One can never pull on more atoms */
1051 snew(old_forces
, nat
);
1055 ImdSession::Impl::Impl(const MDLogger
&mdlog
)
1061 ImdSession::Impl::~Impl()
1065 gmx_fio_fclose(outf
);
1072 ImdSession::Impl::prepareMoleculesInImdGroup(const gmx_mtop_t
*top_global
)
1074 /* check whether index is sorted */
1075 for (int i
= 0; i
< nat
-1; i
++)
1077 if (ind
[i
] > ind
[i
+1])
1079 gmx_fatal(FARGS
, "%s IMD index is not sorted. This is currently not supported.\n", IMDstr
);
1083 RangePartitioning gmols
= gmx_mtop_molecules(*top_global
);
1086 snew(lmols
.index
, gmols
.numBlocks() + 1);
1089 for (int i
= 0; i
< gmols
.numBlocks(); i
++)
1091 auto mol
= gmols
.block(i
);
1093 for (int ii
= 0; ii
< nat
; ii
++)
1095 if (mol
.inRange(ind
[ii
]))
1102 lmols
.index
[lmols
.nr
+1] = lmols
.index
[lmols
.nr
]+count
;
1106 srenew(lmols
.index
, lmols
.nr
+1);
1107 lmols
.nalloc_index
= lmols
.nr
+1;
1112 /*! \brief Copied and modified from groupcoord.c shift_positions_group(). */
1113 static void shift_positions(
1115 rvec x
[], /* The positions [0..nr] */
1116 const ivec is
, /* The shift [0..nr] */
1117 int nr
) /* The number of positions */
1121 /* Loop over the group's atoms */
1124 for (i
= 0; i
< nr
; i
++)
1130 x
[i
][XX
] = x
[i
][XX
]-tx
*box
[XX
][XX
]-ty
*box
[YY
][XX
]-tz
*box
[ZZ
][XX
];
1131 x
[i
][YY
] = x
[i
][YY
]-ty
*box
[YY
][YY
]-tz
*box
[ZZ
][YY
];
1132 x
[i
][ZZ
] = x
[i
][ZZ
]-tz
*box
[ZZ
][ZZ
];
1137 for (i
= 0; i
< nr
; i
++)
1143 x
[i
][XX
] = x
[i
][XX
]-tx
*box
[XX
][XX
];
1144 x
[i
][YY
] = x
[i
][YY
]-ty
*box
[YY
][YY
];
1145 x
[i
][ZZ
] = x
[i
][ZZ
]-tz
*box
[ZZ
][ZZ
];
1151 void ImdSession::Impl::removeMolecularShifts(const matrix box
)
1153 /* for each molecule also present in IMD group */
1154 for (int i
= 0; i
< mols
.nr
; i
++)
1156 /* first we determine the minimum and maximum shifts for each molecule */
1158 ivec largest
, smallest
, shift
;
1159 clear_ivec(largest
);
1160 clear_ivec(smallest
);
1163 copy_ivec(xa_shifts
[mols
.index
[i
]], largest
);
1164 copy_ivec(xa_shifts
[mols
.index
[i
]], smallest
);
1166 for (int ii
= mols
.index
[i
]+1; ii
< mols
.index
[i
+1]; ii
++)
1168 if (xa_shifts
[ii
][XX
] > largest
[XX
])
1170 largest
[XX
] = xa_shifts
[ii
][XX
];
1172 if (xa_shifts
[ii
][XX
] < smallest
[XX
])
1174 smallest
[XX
] = xa_shifts
[ii
][XX
];
1177 if (xa_shifts
[ii
][YY
] > largest
[YY
])
1179 largest
[YY
] = xa_shifts
[ii
][YY
];
1181 if (xa_shifts
[ii
][YY
] < smallest
[YY
])
1183 smallest
[YY
] = xa_shifts
[ii
][YY
];
1186 if (xa_shifts
[ii
][ZZ
] > largest
[ZZ
])
1188 largest
[ZZ
] = xa_shifts
[ii
][ZZ
];
1190 if (xa_shifts
[ii
][ZZ
] < smallest
[ZZ
])
1192 smallest
[ZZ
] = xa_shifts
[ii
][ZZ
];
1197 /* check if we what we can subtract/add to the positions
1198 * to put them back into the central box */
1199 if (smallest
[XX
] > 0)
1201 shift
[XX
] = smallest
[XX
];
1203 if (smallest
[YY
] > 0)
1205 shift
[YY
] = smallest
[YY
];
1207 if (smallest
[ZZ
] > 0)
1209 shift
[ZZ
] = smallest
[ZZ
];
1212 if (largest
[XX
] < 0)
1214 shift
[XX
] = largest
[XX
];
1216 if (largest
[YY
] < 0)
1218 shift
[YY
] = largest
[YY
];
1220 if (largest
[ZZ
] < 0)
1222 shift
[ZZ
] = largest
[ZZ
];
1225 /* is there a shift at all? */
1226 if ((shift
[XX
]) || (shift
[YY
]) || (shift
[ZZ
]))
1228 int molsize
= mols
.index
[i
+1]-mols
.index
[i
];
1229 /* shift the positions */
1230 shift_positions(box
, &(xa
[mols
.index
[i
]]), shift
, molsize
);
1238 ImdSession::Impl::prepareForPositionAssembly(const t_commrec
*cr
, const rvec x
[])
1242 snew(xa_shifts
, nat
);
1243 snew(xa_eshifts
, nat
);
1246 /* Save the original (whole) set of positions such that later the
1247 * molecule can always be made whole again */
1250 for (int i
= 0; i
< nat
; i
++)
1253 copy_rvec(x
[ii
], xa_old
[i
]);
1262 /* xa_ind[i] needs to be set to i for serial runs */
1263 for (int i
= 0; i
< nat
; i
++)
1269 /* Communicate initial coordinates xa_old to all processes */
1272 gmx_bcast(nat
* sizeof(xa_old
[0]), xa_old
, cr
);
1277 /*! \brief Check for non-working integrator / parallel options. */
1278 static void imd_check_integrator_parallel(const t_inputrec
*ir
, const t_commrec
*cr
)
1282 if (((ir
->eI
) == eiSteep
) || ((ir
->eI
) == eiCG
) || ((ir
->eI
) == eiLBFGS
) || ((ir
->eI
) == eiNM
))
1284 gmx_fatal(FARGS
, "%s Energy minimization via steep, CG, lbfgs and nm in parallel is currently not supported by IMD.\n", IMDstr
);
1289 std::unique_ptr
<ImdSession
>
1290 makeImdSession(const t_inputrec
*ir
,
1291 const t_commrec
*cr
,
1292 gmx_wallcycle
*wcycle
,
1293 gmx_enerdata_t
*enerd
,
1294 const gmx_multisim_t
*ms
,
1295 const gmx_mtop_t
*top_global
,
1296 const MDLogger
&mdlog
,
1299 const t_filenm fnm
[],
1300 const gmx_output_env_t
*oenv
,
1301 const ImdOptions
&options
,
1302 const gmx::StartingBehavior startingBehavior
)
1304 std::unique_ptr
<ImdSession
> session(new ImdSession(mdlog
));
1305 auto impl
= session
->impl_
.get();
1307 /* We will allow IMD sessions only if supported by the binary and
1308 explicitly enabled in the .tpr file */
1309 if (!GMX_IMD
|| !ir
->bIMD
)
1314 // TODO As IMD is intended for interactivity, and the .tpr file
1315 // opted in for that, it is acceptable to write more terminal
1316 // output than in a typical simulation. However, all the GMX_LOG
1317 // statements below should go to both the log file and to the
1318 // terminal. This is probably be implemented by adding a logging
1319 // stream named like ImdInfo, to separate it from warning and to
1320 // send it to both destinations.
1322 if (EI_DYNAMICS(ir
->eI
))
1324 impl
->defaultNstImd
= ir
->nstcalcenergy
;
1326 else if (EI_ENERGY_MINIMIZATION(ir
->eI
))
1328 impl
->defaultNstImd
= 1;
1332 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Integrator '%s' is not supported for Interactive Molecular Dynamics, running normally instead", IMDstr
, ei_names
[ir
->eI
]);
1337 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Cannot use IMD for multiple simulations or replica exchange, running normally instead", IMDstr
);
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
, startingBehavior
);
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;