2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
40 * Implements Interactive Molecular Dynamics
42 * Re-implementation of basic IMD functions to work with VMD,
43 * see imdsocket.h for references to the IMD API.
45 * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
58 #include "gromacs/commandline/filenm.h"
59 #include "gromacs/domdec/domdec_struct.h"
60 #include "gromacs/domdec/ga2la.h"
61 #include "gromacs/fileio/confio.h"
62 #include "gromacs/fileio/gmxfio.h"
63 #include "gromacs/fileio/xvgr.h"
64 #include "gromacs/gmxlib/network.h"
65 #include "gromacs/imd/imdsocket.h"
66 #include "gromacs/math/units.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/mdlib/broadcaststructs.h"
69 #include "gromacs/mdlib/groupcoord.h"
70 #include "gromacs/mdlib/sighandler.h"
71 #include "gromacs/mdlib/stat.h"
72 #include "gromacs/mdrunutility/handlerestart.h"
73 #include "gromacs/mdrunutility/multisim.h"
74 #include "gromacs/mdtypes/commrec.h"
75 #include "gromacs/mdtypes/enerdata.h"
76 #include "gromacs/mdtypes/imdmodule.h"
77 #include "gromacs/mdtypes/inputrec.h"
78 #include "gromacs/mdtypes/md_enums.h"
79 #include "gromacs/mdtypes/mdrunoptions.h"
80 #include "gromacs/mdtypes/state.h"
81 #include "gromacs/pbcutil/pbc.h"
82 #include "gromacs/timing/wallcycle.h"
83 #include "gromacs/topology/mtop_util.h"
84 #include "gromacs/topology/topology.h"
85 #include "gromacs/utility/fatalerror.h"
86 #include "gromacs/utility/logger.h"
87 #include "gromacs/utility/smalloc.h"
88 #include "gromacs/utility/stringutil.h"
93 /*! \brief How long shall we wait in seconds until we check for a connection again? */
94 constexpr int c_loopWait
= 1;
96 /*! \brief How long shall we check for the IMD_GO? */
97 constexpr int c_connectWait
= 1;
99 /*! \brief IMD Header Size. */
100 constexpr int c_headerSize
= 8;
102 /*! \brief IMD Protocol Version. */
103 constexpr int c_protocolVersion
= 2;
107 * IMD (interactive molecular dynamics) energy record.
109 * As in the original IMD implementation. Energies in kcal/mol.
110 * NOTE: We return the energies in GROMACS / SI units,
111 * so they also show up as SI in VMD.
116 int32_t tstep
; /**< time step */
117 float T_abs
; /**< absolute temperature */
118 float E_tot
; /**< total energy */
119 float E_pot
; /**< potential energy */
120 float E_vdw
; /**< van der Waals energy */
121 float E_coul
; /**< Coulomb interaction energy */
122 float E_bond
; /**< bonds energy */
123 float E_angle
; /**< angles energy */
124 float E_dihe
; /**< dihedrals energy */
125 float E_impr
; /**< improper dihedrals energy */
130 * \brief IMD (interactive molecular dynamics) communication structure.
132 * This structure defines the IMD communication message header & protocol version.
136 int32_t type
; /**< Type of IMD message, see IMDType_t above */
137 int32_t length
; /**< Length */
142 * \brief Implementation type for the IMD session
144 * \todo Make this class (or one extracted from it) model
146 * \todo Use RAII for files and allocations
148 class ImdSession::Impl
152 Impl(const MDLogger
& mdlog
);
155 /*! \brief Prepare the socket on the MASTER. */
156 void prepareMasterSocket();
157 /*! \brief Disconnect the client. */
158 void disconnectClient();
159 /*! \brief Prints an error message and disconnects the client.
161 * Does not terminate mdrun!
163 void issueFatalError(const char* msg
);
164 /*! \brief Check whether we got an incoming connection. */
166 /*! \brief Wrap imd_tryconnect in order to make it blocking.
168 * Used when the simulation should wait for an incoming connection.
171 /*! \brief Make sure that our array holding the forces received via IMD is large enough. */
172 void prepareVmdForces();
173 /*! \brief Reads forces received via IMD. */
174 void readVmdForces();
175 /*! \brief Prepares the MD force arrays. */
176 void prepareMDForces();
177 /*! \brief Copy IMD forces to MD forces.
179 * Do conversion from Cal->Joule and from
180 * Angstrom -> nm and from a pointer array to arrays to 3*N array.
182 void copyToMDForces();
183 /*! \brief Return true if any of the forces or indices changed. */
184 bool bForcesChanged() const;
185 /*! \brief Update the old_f_ind and old_forces arrays to contain the current values. */
186 void keepOldValues();
187 /*! \brief Write the applied pull forces to logfile.
189 * Call on master only!
191 void outputForces(double time
);
192 /*! \brief Synchronize the nodes. */
193 void syncNodes(const t_commrec
* cr
, double t
);
194 /*! \brief Reads header from the client and decides what to do. */
196 /*! \brief Open IMD output file and write header information.
198 * Call on master only.
200 void openOutputFile(const char* fn
, int nat_total
, const gmx_output_env_t
* oenv
, StartingBehavior startingBehavior
);
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
, bool bNS
, const matrix box
, const rvec x
[], double t
);
210 // TODO rename all the data members to have underscore suffixes
212 //! True if tpr and mdrun input combine to permit IMD sessions
213 bool sessionPossible
= false;
214 //! Output file for IMD data, mainly forces.
215 FILE* outf
= nullptr;
217 //! Number of atoms that can be pulled via IMD.
219 //! Part of the atoms that are local.
221 //! Global indices of the IMD atoms.
223 //! Local indices of the IMD atoms.
224 int* ind_loc
= nullptr;
225 //! Allocation size for ind_loc.
227 //! Positions for all IMD atoms assembled on the master node.
229 //! Shifts for all IMD atoms, to make molecule(s) whole.
230 ivec
* xa_shifts
= nullptr;
231 //! Extra shifts since last DD step.
232 ivec
* xa_eshifts
= nullptr;
233 //! Old positions for all IMD atoms on master.
234 rvec
* xa_old
= nullptr;
235 //! Position of each local atom in the collective array.
236 int* xa_ind
= nullptr;
238 //! Global IMD frequency, known to all ranks.
240 //! New frequency from IMD client, master only.
242 //! Default IMD frequency when disconnected.
243 int defaultNstImd
= -1;
245 //! Port to use for network socket.
247 //! The IMD socket on the master node.
248 IMDSocket
* socket
= nullptr;
249 //! The IMD socket on the client.
250 IMDSocket
* clientsocket
= nullptr;
251 //! Length we got with last header.
254 //! Shall we block and wait for connection?
255 bool bWConnect
= false;
256 //! Set if MD is terminated.
257 bool bTerminated
= false;
258 //! Set if MD can be terminated.
259 bool bTerminatable
= false;
260 //! Set if connection is present.
261 bool bConnected
= false;
262 //! Set if we received new forces.
263 bool bNewForces
= false;
264 //! Set if pulling from VMD is allowed.
265 bool bForceActivated
= false;
267 //! Pointer to energies we send back.
268 IMDEnergyBlock
* energies
= nullptr;
270 //! Number of VMD forces.
271 int32_t vmd_nforces
= 0;
272 //! VMD forces indices.
273 int32_t* vmd_f_ind
= nullptr;
274 //! The VMD forces flat in memory.
275 float* vmd_forces
= nullptr;
276 //! Number of actual MD forces; this gets communicated to the clients.
279 int* f_ind
= nullptr;
280 //! The IMD pulling forces.
283 //! Buffer for force sending.
284 char* forcesendbuf
= nullptr;
285 //! Buffer for coordinate sending.
286 char* coordsendbuf
= nullptr;
287 //! Send buffer for energies.
288 char* energysendbuf
= nullptr;
289 //! Buffer to make molecules whole before sending.
290 rvec
* sendxbuf
= nullptr;
292 //! Molecules block in IMD group.
295 /* The next block is used on the master node only to reduce the output
296 * without sacrificing information. If any of these values changes,
297 * we need to write output */
298 //! Old value for nforces.
300 //! Old values for force indices.
301 int* old_f_ind
= nullptr;
302 //! Old values for IMD pulling forces.
303 rvec
* old_forces
= nullptr;
306 const MDLogger
& mdlog
;
307 //! Commmunication object
308 const t_commrec
* cr
= nullptr;
309 //! Wallcycle counting manager.
310 gmx_wallcycle
* wcycle
= nullptr;
311 //! Energy output handler
312 gmx_enerdata_t
* enerd
= nullptr;
316 * \brief Implement interactive molecular dynamics.
318 * \todo Some aspects of this module provides forces (when the user
319 * pulls on things in VMD), so in future it should have a class that
320 * models IForceProvider and is contributed to the collection of such
323 class InteractiveMolecularDynamics final
: public IMDModule
326 IMdpOptionProvider
* mdpOptionProvider() override
{ return nullptr; }
327 IMDOutputProvider
* outputProvider() override
{ return nullptr; }
328 void initForceProviders(ForceProviders
* /* forceProviders */) override
{}
329 void subscribeToSimulationSetupNotifications(MdModulesNotifier
* /* notifier */) override
{}
330 void subscribeToPreProcessingNotifications(MdModulesNotifier
* /* notifier */) override
{}
333 std::unique_ptr
<IMDModule
> createInteractiveMolecularDynamicsModule()
335 return std::make_unique
<InteractiveMolecularDynamics
>();
338 /*! \brief Enum for types of IMD messages.
340 * We use the same records as the NAMD/VMD IMD implementation.
342 typedef enum IMDType_t
344 IMD_DISCONNECT
, /**< client disconnect */
345 IMD_ENERGIES
, /**< energy data */
346 IMD_FCOORDS
, /**< atomic coordinates */
347 IMD_GO
, /**< start command for the simulation */
348 IMD_HANDSHAKE
, /**< handshake to determine little/big endianness */
349 IMD_KILL
, /**< terminates the simulation */
350 IMD_MDCOMM
, /**< force data */
351 IMD_PAUSE
, /**< pauses the simulation */
352 IMD_TRATE
, /**< sets the IMD transmission and processing rate */
353 IMD_IOERROR
, /**< I/O error */
354 IMD_NR
/**< number of entries */
358 /*! \brief Names of the IMDType for error messages.
360 static const char* eIMDType_names
[IMD_NR
+ 1] = { "IMD_DISCONNECT", "IMD_ENERGIES", "IMD_FCOORDS",
361 "IMD_GO", "IMD_HANDSHAKE", "IMD_KILL",
362 "IMD_MDCOMM", "IMD_PAUSE", "IMD_TRATE",
363 "IMD_IOERROR", nullptr };
366 /*! \brief Fills the header with message and the length argument. */
367 static void fill_header(IMDHeader
* header
, IMDMessageType type
, int32_t length
)
369 /* We (ab-)use htonl network function for the correct endianness */
370 header
->type
= imd_htonl(static_cast<int32_t>(type
));
371 header
->length
= imd_htonl(length
);
375 /*! \brief Swaps the endianess of the header. */
376 static void swap_header(IMDHeader
* header
)
378 /* and vice versa... */
379 header
->type
= imd_ntohl(header
->type
);
380 header
->length
= imd_ntohl(header
->length
);
384 /*! \brief Reads multiple bytes from socket. */
385 static int32_t imd_read_multiple(IMDSocket
* socket
, char* datptr
, int32_t toread
)
387 int32_t leftcount
, countread
;
391 /* Try to read while we haven't reached toread */
392 while (leftcount
!= 0)
394 if ((countread
= imdsock_read(socket
, datptr
, leftcount
)) < 0)
396 /* interrupted function call, try again... */
401 /* this is an unexpected error, return what we got */
404 return toread
- leftcount
;
407 /* nothing read, finished */
409 else if (countread
== 0)
413 leftcount
-= countread
;
417 /* return nr of bytes read */
418 return toread
- leftcount
;
422 /*! \brief Writes multiple bytes to socket in analogy to imd_read_multiple. */
423 static int32_t imd_write_multiple(IMDSocket
* socket
, const char* datptr
, int32_t towrite
)
425 int32_t leftcount
, countwritten
;
429 while (leftcount
!= 0)
431 if ((countwritten
= imdsock_write(socket
, datptr
, leftcount
)) <= 0)
439 return towrite
- leftcount
;
442 leftcount
-= countwritten
;
443 datptr
+= countwritten
;
446 return towrite
- leftcount
;
450 /*! \brief Handshake with IMD client. */
451 static int imd_handshake(IMDSocket
* socket
)
456 fill_header(&header
, IMD_HANDSHAKE
, 1);
457 header
.length
= c_protocolVersion
; /* client wants unswapped version */
459 return static_cast<int>(imd_write_multiple(socket
, reinterpret_cast<char*>(&header
), c_headerSize
)
464 /*! \brief Send energies using the energy block and the send buffer. */
465 static int imd_send_energies(IMDSocket
* socket
, const IMDEnergyBlock
* energies
, char* buffer
)
470 recsize
= c_headerSize
+ sizeof(IMDEnergyBlock
);
471 fill_header(reinterpret_cast<IMDHeader
*>(buffer
), IMD_ENERGIES
, 1);
472 memcpy(buffer
+ c_headerSize
, energies
, sizeof(IMDEnergyBlock
));
474 return static_cast<int>(imd_write_multiple(socket
, buffer
, recsize
) != recsize
);
478 /*! \brief Receive IMD header from socket, sets the length and returns the IMD message. */
479 static IMDMessageType
imd_recv_header(IMDSocket
* socket
, int32_t* length
)
484 if (imd_read_multiple(socket
, reinterpret_cast<char*>(&header
), c_headerSize
) != c_headerSize
)
488 swap_header(&header
);
489 *length
= header
.length
;
491 return static_cast<IMDMessageType
>(header
.type
);
495 /*! \brief Receive force indices and forces.
497 * The number of forces was previously communicated via the header.
499 static bool imd_recv_mdcomm(IMDSocket
* socket
, int32_t nforces
, int32_t* forcendx
, float* forces
)
501 /* reading indices */
502 int retsize
= sizeof(int32_t) * nforces
;
503 int retbytes
= imd_read_multiple(socket
, reinterpret_cast<char*>(forcendx
), retsize
);
504 if (retbytes
!= retsize
)
509 /* reading forces as float array */
510 retsize
= 3 * sizeof(float) * nforces
;
511 retbytes
= imd_read_multiple(socket
, reinterpret_cast<char*>(forces
), retsize
);
512 return (retbytes
== retsize
);
515 /* GROMACS specific functions for the IMD implementation */
516 void write_IMDgroup_to_file(bool bIMD
,
518 const t_state
* state
,
519 const gmx_mtop_t
* sys
,
521 const t_filenm fnm
[])
528 IMDatoms
= gmx_mtop_global_atoms(sys
);
529 write_sto_conf_indexed(opt2fn("-imd", nfile
, fnm
), "IMDgroup", &IMDatoms
,
530 state
->x
.rvec_array(), state
->v
.rvec_array(), ir
->pbcType
,
531 state
->box
, ir
->imd
->nat
, ir
->imd
->ind
);
536 void ImdSession::dd_make_local_IMD_atoms(const gmx_domdec_t
* dd
)
538 if (!impl_
->sessionPossible
)
543 dd_make_local_group_indices(dd
->ga2la
, impl_
->nat
, impl_
->ind
, &impl_
->nat_loc
, &impl_
->ind_loc
,
544 &impl_
->nalloc_loc
, impl_
->xa_ind
);
548 /*! \brief Send positions from rvec.
550 * We need a separate send buffer and conversion to Angstrom.
552 static int imd_send_rvecs(IMDSocket
* socket
, int nat
, rvec
* x
, char* buffer
)
557 int tuplesize
= 3 * sizeof(float);
560 /* Required size for the send buffer */
561 size
= c_headerSize
+ 3 * sizeof(float) * nat
;
564 fill_header(reinterpret_cast<IMDHeader
*>(buffer
), IMD_FCOORDS
, static_cast<int32_t>(nat
));
565 for (i
= 0; i
< nat
; i
++)
567 sendx
[0] = static_cast<float>(x
[i
][0]) * NM2A
;
568 sendx
[1] = static_cast<float>(x
[i
][1]) * NM2A
;
569 sendx
[2] = static_cast<float>(x
[i
][2]) * NM2A
;
570 memcpy(buffer
+ c_headerSize
+ i
* tuplesize
, sendx
, tuplesize
);
573 return static_cast<int>(imd_write_multiple(socket
, buffer
, size
) != size
);
577 void ImdSession::Impl::prepareMasterSocket()
579 if (imdsock_winsockinit() == -1)
581 gmx_fatal(FARGS
, "%s Failed to initialize winsock.\n", IMDstr
);
584 /* The rest is identical, first create and bind a socket and set to listen then. */
585 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Setting up incoming socket.", IMDstr
);
586 socket
= imdsock_create();
589 gmx_fatal(FARGS
, "%s Failed to create socket.", IMDstr
);
593 int ret
= imdsock_bind(socket
, port
);
596 gmx_fatal(FARGS
, "%s binding socket to port %d failed with error %d.\n", IMDstr
, port
, ret
);
599 if (imd_sock_listen(socket
))
601 gmx_fatal(FARGS
, "%s socket listen failed with error %d.\n", IMDstr
, ret
);
604 if (imdsock_getport(socket
, &port
))
606 gmx_fatal(FARGS
, "%s Could not determine port number.\n", IMDstr
);
609 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Listening for IMD connection on port %d.", IMDstr
, port
);
613 void ImdSession::Impl::disconnectClient()
615 /* Write out any buffered pulling data */
618 /* we first try to shut down the clientsocket */
619 imdsock_shutdown(clientsocket
);
620 if (!imdsock_destroy(clientsocket
))
622 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Failed to destroy socket.", IMDstr
);
625 /* then we reset the IMD step to its default, and reset the connection boolean */
626 nstimd_new
= defaultNstImd
;
627 clientsocket
= nullptr;
632 void ImdSession::Impl::issueFatalError(const char* msg
)
634 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s %s", IMDstr
, msg
);
636 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s disconnected.", IMDstr
);
640 bool ImdSession::Impl::tryConnect()
642 if (imdsock_tryread(socket
, 0, 0) > 0)
644 /* yes, we got something, accept on clientsocket */
645 clientsocket
= imdsock_accept(socket
);
648 GMX_LOG(mdlog
.warning
)
649 .appendTextFormatted("%s Accepting the connection on the socket failed.", IMDstr
);
653 /* handshake with client */
654 if (imd_handshake(clientsocket
))
656 issueFatalError("Connection failed.");
660 GMX_LOG(mdlog
.warning
)
661 .appendTextFormatted("%s Connection established, checking if I got IMD_GO orders.", IMDstr
);
663 /* Check if we get the proper "GO" command from client. */
664 if (imdsock_tryread(clientsocket
, c_connectWait
, 0) != 1
665 || imd_recv_header(clientsocket
, &(length
)) != IMD_GO
)
667 issueFatalError("No IMD_GO order received. IMD connection failed.");
680 void ImdSession::Impl::blockConnect()
682 /* do not wait for connection, when e.g. ctrl+c is pressed and we will terminate anyways. */
683 if (!(static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none
))
688 GMX_LOG(mdlog
.warning
)
689 .appendTextFormatted("%s Will wait until I have a connection and IMD_GO orders.", IMDstr
);
691 /* while we have no clientsocket... 2nd part: we should still react on ctrl+c */
692 while ((!clientsocket
) && (static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none
))
695 imd_sleep(c_loopWait
);
700 void ImdSession::Impl::prepareVmdForces()
702 srenew((vmd_f_ind
), vmd_nforces
);
703 srenew((vmd_forces
), 3 * vmd_nforces
);
707 void ImdSession::Impl::readVmdForces()
709 /* the length of the previously received header tells us the nr of forces we will receive */
710 vmd_nforces
= length
;
711 /* prepare the arrays */
713 /* Now we read the forces... */
714 if (!(imd_recv_mdcomm(clientsocket
, vmd_nforces
, vmd_f_ind
, vmd_forces
)))
716 issueFatalError("Error while reading forces from remote. Disconnecting");
721 void ImdSession::Impl::prepareMDForces()
723 srenew((f_ind
), nforces
);
724 srenew((f
), nforces
);
728 void ImdSession::Impl::copyToMDForces()
731 real conversion
= CAL2JOULE
* NM2A
;
734 for (i
= 0; i
< nforces
; i
++)
736 /* Copy the indices, a copy is important because we may update the incoming forces
737 * whenever we receive new forces while the MD forces are only communicated upon
738 * IMD communication.*/
739 f_ind
[i
] = vmd_f_ind
[i
];
741 /* Convert to rvecs and do a proper unit conversion */
742 f
[i
][0] = vmd_forces
[3 * i
] * conversion
;
743 f
[i
][1] = vmd_forces
[3 * i
+ 1] * conversion
;
744 f
[i
][2] = vmd_forces
[3 * i
+ 2] * conversion
;
749 /*! \brief Returns true if any component of the two rvecs differs. */
750 static inline bool rvecs_differ(const rvec v1
, const rvec v2
)
752 for (int i
= 0; i
< DIM
; i
++)
763 bool ImdSession::Impl::bForcesChanged() const
765 /* First, check whether the number of pulled atoms changed */
766 if (nforces
!= old_nforces
)
771 /* Second, check whether any of the involved atoms changed */
772 for (int i
= 0; i
< nforces
; i
++)
774 if (f_ind
[i
] != old_f_ind
[i
])
780 /* Third, check whether all forces are the same */
781 for (int i
= 0; i
< nforces
; i
++)
783 if (rvecs_differ(f
[i
], old_forces
[i
]))
789 /* All old and new forces are identical! */
794 void ImdSession::Impl::keepOldValues()
796 old_nforces
= nforces
;
798 for (int i
= 0; i
< nforces
; i
++)
800 old_f_ind
[i
] = f_ind
[i
];
801 copy_rvec(f
[i
], old_forces
[i
]);
806 void ImdSession::Impl::outputForces(double time
)
808 if (!bForcesChanged())
813 /* Write time and total number of applied IMD forces */
814 fprintf(outf
, "%14.6e%6d", time
, nforces
);
816 /* Write out the global atom indices of the pulled atoms and the forces itself,
817 * write out a force only if it has changed since the last output */
818 for (int i
= 0; i
< nforces
; i
++)
820 if (rvecs_differ(f
[i
], old_forces
[i
]))
822 fprintf(outf
, "%9d", ind
[f_ind
[i
]] + 1);
823 fprintf(outf
, "%12.4e%12.4e%12.4e", f
[i
][0], f
[i
][1], f
[i
][2]);
832 void ImdSession::Impl::syncNodes(const t_commrec
* cr
, double t
)
834 /* Notify the other nodes whether we are still connected. */
837 block_bc(cr
->mpi_comm_mygroup
, bConnected
);
840 /* ...if not connected, the job is done here. */
846 /* Let the other nodes know whether we got a new IMD synchronization frequency. */
849 block_bc(cr
->mpi_comm_mygroup
, nstimd_new
);
852 /* Now we all set the (new) nstimd communication time step */
855 /* We're done if we don't allow pulling at all */
856 if (!(bForceActivated
))
861 /* OK, let's check if we have received forces which we need to communicate
862 * to the other nodes */
868 new_nforces
= vmd_nforces
;
870 /* make the "new_forces" negative, when we did not receive new ones */
873 new_nforces
= vmd_nforces
* -1;
877 /* make new_forces known to the clients */
880 block_bc(cr
->mpi_comm_mygroup
, new_nforces
);
883 /* When new_natoms < 0 then we know that these are still the same forces
884 * so we don't communicate them, otherwise... */
890 /* set local VMD and nforces */
891 vmd_nforces
= new_nforces
;
892 nforces
= vmd_nforces
;
894 /* now everybody knows the number of forces in f_ind, so we can prepare
895 * the target arrays for indices and forces */
898 /* we first update the MD forces on the master by converting the VMD forces */
902 /* We also write out forces on every update, so that we know which
903 * forces are applied for every step */
910 /* In parallel mode we communicate the to-be-applied forces to the other nodes */
913 nblock_bc(cr
->mpi_comm_mygroup
, nforces
, f_ind
);
914 nblock_bc(cr
->mpi_comm_mygroup
, nforces
, f
);
917 /* done communicating the forces, reset bNewForces */
922 void ImdSession::Impl::readCommand()
924 bool IMDpaused
= false;
926 while (clientsocket
&& (imdsock_tryread(clientsocket
, 0, 0) > 0 || IMDpaused
))
928 IMDMessageType itype
= imd_recv_header(clientsocket
, &(length
));
929 /* let's see what we got: */
932 /* IMD asks us to terminate the simulation, check if the user allowed this */
936 GMX_LOG(mdlog
.warning
)
937 .appendTextFormatted(
938 " %s Terminating connection and running simulation (if "
939 "supported by integrator).",
943 gmx_set_stop_condition(gmx_stop_cond_next
);
947 GMX_LOG(mdlog
.warning
)
948 .appendTextFormatted(
949 " %s Set -imdterm command line switch to allow mdrun "
950 "termination from within IMD.",
956 /* the client doen't want to talk to us anymore */
958 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Disconnecting client.", IMDstr
);
962 /* we got new forces, read them and set bNewForces flag */
968 /* the client asks us to (un)pause the simulation. So we toggle the IMDpaused state */
972 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Un-pause command received.", IMDstr
);
977 GMX_LOG(mdlog
.warning
).appendTextFormatted(" %s Pause command received.", IMDstr
);
983 /* the client sets a new transfer rate, if we get 0, we reset the rate
984 * to the default. VMD filters 0 however */
986 nstimd_new
= (length
> 0) ? length
: defaultNstImd
;
987 GMX_LOG(mdlog
.warning
)
988 .appendTextFormatted(" %s Update frequency will be set to %d.", IMDstr
, nstimd_new
);
991 /* Catch all rule for the remaining IMD types which we don't expect */
993 GMX_LOG(mdlog
.warning
)
994 .appendTextFormatted(" %s Received unexpected %s.", IMDstr
,
995 enum_name(static_cast<int>(itype
), IMD_NR
, eIMDType_names
));
996 issueFatalError("Terminating connection");
1003 void ImdSession::Impl::openOutputFile(const char* fn
,
1005 const gmx_output_env_t
* oenv
,
1006 const gmx::StartingBehavior startingBehavior
)
1008 /* Open log file of applied IMD forces if requested */
1012 "%s For a log of the IMD pull forces explicitly specify '-if' on the command "
1014 "%s (Not possible with energy minimization.)\n",
1019 /* If we append to an existing file, all the header information is already there */
1020 if (startingBehavior
== StartingBehavior::RestartWithAppending
)
1022 outf
= gmx_fio_fopen(fn
, "a+");
1026 outf
= gmx_fio_fopen(fn
, "w+");
1027 if (nat
== nat_total
)
1030 "# Note that you can select an IMD index group in the .mdp file if a subset of "
1031 "the atoms suffices.\n");
1034 xvgr_header(outf
, "IMD Pull Forces", "Time (ps)",
1035 "# of Forces / Atom IDs / Forces (kJ/mol)", exvggtNONE
, oenv
);
1037 fprintf(outf
, "# Can display and manipulate %d (of a total of %d) atoms via IMD.\n", nat
, nat_total
);
1038 fprintf(outf
, "# column 1 : time (ps)\n");
1040 "# column 2 : total number of atoms feeling an IMD pulling force at that "
1043 "# cols. 3.-6 : global atom number of pulled atom, x-force, y-force, z-force "
1046 "# then follow : atom-ID, f[x], f[y], f[z] for more atoms in case the force on "
1047 "multiple atoms is changed simultaneously.\n");
1049 "# Note that the force on any atom is always equal to the last value for that "
1050 "atom-ID found in the data.\n");
1054 /* To reduce the output file size we remember the old values and output only
1055 * when something changed */
1056 snew(old_f_ind
, nat
); /* One can never pull on more atoms */
1057 snew(old_forces
, nat
);
1061 ImdSession::Impl::Impl(const MDLogger
& mdlog
) : mdlog(mdlog
)
1066 ImdSession::Impl::~Impl()
1070 gmx_fio_fclose(outf
);
1076 void ImdSession::Impl::prepareMoleculesInImdGroup(const gmx_mtop_t
* top_global
)
1078 /* check whether index is sorted */
1079 for (int i
= 0; i
< nat
- 1; i
++)
1081 if (ind
[i
] > ind
[i
+ 1])
1083 gmx_fatal(FARGS
, "%s IMD index is not sorted. This is currently not supported.\n", IMDstr
);
1087 RangePartitioning gmols
= gmx_mtop_molecules(*top_global
);
1090 snew(lmols
.index
, gmols
.numBlocks() + 1);
1093 for (int i
= 0; i
< gmols
.numBlocks(); i
++)
1095 auto mol
= gmols
.block(i
);
1097 for (int ii
= 0; ii
< nat
; ii
++)
1099 if (mol
.isInRange(ind
[ii
]))
1106 lmols
.index
[lmols
.nr
+ 1] = lmols
.index
[lmols
.nr
] + count
;
1110 srenew(lmols
.index
, lmols
.nr
+ 1);
1111 lmols
.nalloc_index
= lmols
.nr
+ 1;
1116 /*! \brief Copied and modified from groupcoord.c shift_positions_group(). */
1117 static void shift_positions(const matrix box
,
1118 rvec x
[], /* The positions [0..nr] */
1119 const ivec is
, /* The shift [0..nr] */
1120 int nr
) /* The number of positions */
1124 /* Loop over the group's atoms */
1127 for (i
= 0; i
< nr
; i
++)
1133 x
[i
][XX
] = x
[i
][XX
] - tx
* box
[XX
][XX
] - ty
* box
[YY
][XX
] - tz
* box
[ZZ
][XX
];
1134 x
[i
][YY
] = x
[i
][YY
] - ty
* box
[YY
][YY
] - tz
* box
[ZZ
][YY
];
1135 x
[i
][ZZ
] = x
[i
][ZZ
] - tz
* box
[ZZ
][ZZ
];
1140 for (i
= 0; i
< nr
; i
++)
1146 x
[i
][XX
] = x
[i
][XX
] - tx
* box
[XX
][XX
];
1147 x
[i
][YY
] = x
[i
][YY
] - ty
* box
[YY
][YY
];
1148 x
[i
][ZZ
] = x
[i
][ZZ
] - tz
* box
[ZZ
][ZZ
];
1154 void ImdSession::Impl::removeMolecularShifts(const matrix box
)
1156 /* for each molecule also present in IMD group */
1157 for (int i
= 0; i
< mols
.nr
; i
++)
1159 /* first we determine the minimum and maximum shifts for each molecule */
1161 ivec largest
, smallest
, shift
;
1162 clear_ivec(largest
);
1163 clear_ivec(smallest
);
1166 copy_ivec(xa_shifts
[mols
.index
[i
]], largest
);
1167 copy_ivec(xa_shifts
[mols
.index
[i
]], smallest
);
1169 for (int ii
= mols
.index
[i
] + 1; ii
< mols
.index
[i
+ 1]; ii
++)
1171 if (xa_shifts
[ii
][XX
] > largest
[XX
])
1173 largest
[XX
] = xa_shifts
[ii
][XX
];
1175 if (xa_shifts
[ii
][XX
] < smallest
[XX
])
1177 smallest
[XX
] = xa_shifts
[ii
][XX
];
1180 if (xa_shifts
[ii
][YY
] > largest
[YY
])
1182 largest
[YY
] = xa_shifts
[ii
][YY
];
1184 if (xa_shifts
[ii
][YY
] < smallest
[YY
])
1186 smallest
[YY
] = xa_shifts
[ii
][YY
];
1189 if (xa_shifts
[ii
][ZZ
] > largest
[ZZ
])
1191 largest
[ZZ
] = xa_shifts
[ii
][ZZ
];
1193 if (xa_shifts
[ii
][ZZ
] < smallest
[ZZ
])
1195 smallest
[ZZ
] = xa_shifts
[ii
][ZZ
];
1199 /* check if we what we can subtract/add to the positions
1200 * to put them back into the central box */
1201 if (smallest
[XX
] > 0)
1203 shift
[XX
] = smallest
[XX
];
1205 if (smallest
[YY
] > 0)
1207 shift
[YY
] = smallest
[YY
];
1209 if (smallest
[ZZ
] > 0)
1211 shift
[ZZ
] = smallest
[ZZ
];
1214 if (largest
[XX
] < 0)
1216 shift
[XX
] = largest
[XX
];
1218 if (largest
[YY
] < 0)
1220 shift
[YY
] = largest
[YY
];
1222 if (largest
[ZZ
] < 0)
1224 shift
[ZZ
] = largest
[ZZ
];
1227 /* is there a shift at all? */
1228 if ((shift
[XX
]) || (shift
[YY
]) || (shift
[ZZ
]))
1230 int molsize
= mols
.index
[i
+ 1] - mols
.index
[i
];
1231 /* shift the positions */
1232 shift_positions(box
, &(xa
[mols
.index
[i
]]), shift
, molsize
);
1238 void 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
->mpi_comm_mygroup
);
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
))
1285 "%s Energy minimization via steep, CG, lbfgs and nm in parallel is currently "
1286 "not supported by IMD.\n",
1292 std::unique_ptr
<ImdSession
> makeImdSession(const t_inputrec
* ir
,
1293 const t_commrec
* cr
,
1294 gmx_wallcycle
* wcycle
,
1295 gmx_enerdata_t
* enerd
,
1296 const gmx_multisim_t
* ms
,
1297 const gmx_mtop_t
* top_global
,
1298 const MDLogger
& mdlog
,
1301 const t_filenm fnm
[],
1302 const gmx_output_env_t
* oenv
,
1303 const ImdOptions
& options
,
1304 const gmx::StartingBehavior startingBehavior
)
1306 std::unique_ptr
<ImdSession
> session(new ImdSession(mdlog
));
1307 auto impl
= session
->impl_
.get();
1309 /* We will allow IMD sessions only if supported by the binary and
1310 explicitly enabled in the .tpr file */
1311 if (!GMX_IMD
|| !ir
->bIMD
)
1316 // TODO As IMD is intended for interactivity, and the .tpr file
1317 // opted in for that, it is acceptable to write more terminal
1318 // output than in a typical simulation. However, all the GMX_LOG
1319 // statements below should go to both the log file and to the
1320 // terminal. This is probably be implemented by adding a logging
1321 // stream named like ImdInfo, to separate it from warning and to
1322 // send it to both destinations.
1324 if (EI_DYNAMICS(ir
->eI
))
1326 impl
->defaultNstImd
= ir
->nstcalcenergy
;
1328 else if (EI_ENERGY_MINIMIZATION(ir
->eI
))
1330 impl
->defaultNstImd
= 1;
1334 GMX_LOG(mdlog
.warning
)
1335 .appendTextFormatted(
1336 "%s Integrator '%s' is not supported for Interactive Molecular Dynamics, "
1337 "running normally instead",
1338 IMDstr
, ei_names
[ir
->eI
]);
1343 GMX_LOG(mdlog
.warning
)
1344 .appendTextFormatted(
1345 "%s Cannot use IMD for multiple simulations or replica exchange, running "
1351 bool createSession
= false;
1352 /* It seems we have a .tpr file that defines an IMD group and thus allows IMD connections.
1353 * Check whether we can actually provide the IMD functionality for this setting: */
1356 /* Check whether IMD was enabled by one of the command line switches: */
1357 if (options
.wait
|| options
.terminatable
|| options
.pull
)
1359 GMX_LOG(mdlog
.warning
)
1360 .appendTextFormatted(
1361 "%s Enabled. This simulation will accept incoming IMD connections.", IMDstr
);
1362 createSession
= true;
1366 GMX_LOG(mdlog
.warning
)
1367 .appendTextFormatted(
1368 "%s None of the -imd switches was used.\n"
1369 "%s This run will not accept incoming IMD connections",
1372 } /* end master only */
1374 /* Let the other nodes know whether we want IMD */
1377 block_bc(cr
->mpi_comm_mygroup
, createSession
);
1380 /*... if not we are done.*/
1387 /* check if we're using a sane integrator / parallel combination */
1388 imd_check_integrator_parallel(ir
, cr
);
1392 *****************************************************
1393 * From here on we assume that IMD is turned on *
1394 *****************************************************
1397 int nat_total
= top_global
->natoms
;
1399 /* Initialize IMD session. If we read in a pre-IMD .tpr file, ir->imd->nat
1400 * will be zero. For those cases we transfer _all_ atomic positions */
1401 impl
->sessionPossible
= true;
1402 impl
->nat
= ir
->imd
->nat
> 0 ? ir
->imd
->nat
: nat_total
;
1403 if (options
.port
>= 1)
1405 impl
->port
= options
.port
;
1408 impl
->wcycle
= wcycle
;
1409 impl
->enerd
= enerd
;
1411 /* We might need to open an output file for IMD forces data */
1414 impl
->openOutputFile(opt2fn("-if", nfile
, fnm
), nat_total
, oenv
, startingBehavior
);
1417 /* Make sure that we operate with a valid atom index array for the IMD atoms */
1418 if (ir
->imd
->nat
> 0)
1420 /* Point to the user-supplied array of atom numbers */
1421 impl
->ind
= ir
->imd
->ind
;
1425 /* Make a dummy (ind[i] = i) array of all atoms */
1426 snew(impl
->ind
, nat_total
);
1427 for (int i
= 0; i
< nat_total
; i
++)
1433 /* read environment on master and prepare socket for incoming connections */
1436 /* we allocate memory for our IMD energy structure */
1437 int32_t recsize
= c_headerSize
+ sizeof(IMDEnergyBlock
);
1438 snew(impl
->energysendbuf
, recsize
);
1440 /* Shall we wait for a connection? */
1443 impl
->bWConnect
= true;
1444 GMX_LOG(mdlog
.warning
)
1445 .appendTextFormatted(
1446 "%s Pausing simulation while no IMD connection present (-imdwait).", IMDstr
);
1449 /* Will the IMD clients be able to terminate the simulation? */
1450 if (options
.terminatable
)
1452 impl
->bTerminatable
= true;
1453 GMX_LOG(mdlog
.warning
)
1454 .appendTextFormatted(
1455 "%s Allow termination of the simulation from IMD client (-imdterm).", IMDstr
);
1458 /* Is pulling from IMD client allowed? */
1461 impl
->bForceActivated
= true;
1462 GMX_LOG(mdlog
.warning
)
1463 .appendTextFormatted("%s Pulling from IMD remote is enabled (-imdpull).", IMDstr
);
1466 /* Initialize send buffers with constant size */
1467 snew(impl
->sendxbuf
, impl
->nat
);
1468 snew(impl
->energies
, 1);
1469 int32_t bufxsize
= c_headerSize
+ 3 * sizeof(float) * impl
->nat
;
1470 snew(impl
->coordsendbuf
, bufxsize
);
1473 /* do we allow interactive pulling? If so let the other nodes know. */
1476 block_bc(cr
->mpi_comm_mygroup
, impl
->bForceActivated
);
1479 /* setup the listening socket on master process */
1482 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s Setting port for connection requests to %d.", IMDstr
, impl
->port
);
1483 impl
->prepareMasterSocket();
1484 /* Wait until we have a connection if specified before */
1485 if (impl
->bWConnect
)
1487 impl
->blockConnect();
1491 GMX_LOG(mdlog
.warning
).appendTextFormatted("%s -imdwait not set, starting simulation.", IMDstr
);
1494 /* Let the other nodes know whether we are connected */
1495 impl
->syncNodes(cr
, 0);
1497 /* Initialize arrays used to assemble the positions from the other nodes */
1498 impl
->prepareForPositionAssembly(cr
, x
);
1500 /* Initialize molecule blocks to make them whole later...*/
1503 impl
->prepareMoleculesInImdGroup(top_global
);
1510 bool ImdSession::Impl::run(int64_t step
, bool bNS
, const matrix box
, const rvec x
[], double t
)
1513 if (!sessionPossible
)
1518 wallcycle_start(wcycle
, ewcIMD
);
1520 /* read command from client and check if new incoming connection */
1523 /* If not already connected, check for new connections */
1536 /* Let's see if we have new IMD messages for us */
1543 /* is this an IMD communication step? */
1544 bool imdstep
= do_per_step(step
, nstimd
);
1546 /* OK so this is an IMD step ... */
1549 /* First we sync all nodes to let everybody know whether we are connected to VMD */
1553 /* If a client is connected, we collect the positions
1554 * and put molecules back into the box before transfer */
1555 if ((imdstep
&& bConnected
) || bNS
) /* independent of imdstep, we communicate positions at each NS step */
1557 /* Transfer the IMD positions to the master node. Every node contributes
1558 * its local positions x and stores them in the assembled xa array. */
1559 communicate_group_positions(cr
, xa
, xa_shifts
, xa_eshifts
, true, x
, nat
, nat_loc
, ind_loc
,
1560 xa_ind
, xa_old
, box
);
1562 /* If connected and master -> remove shifts */
1563 if ((imdstep
&& bConnected
) && MASTER(cr
))
1565 removeMolecularShifts(box
);
1569 wallcycle_stop(wcycle
, ewcIMD
);
1574 bool ImdSession::run(int64_t step
, bool bNS
, const matrix box
, const rvec x
[], double t
)
1576 return impl_
->run(step
, bNS
, box
, x
, t
);
1579 void ImdSession::fillEnergyRecord(int64_t step
, bool bHaveNewEnergies
)
1581 if (!impl_
->sessionPossible
|| !impl_
->clientsocket
)
1586 IMDEnergyBlock
* ene
= impl_
->energies
;
1590 /* In MPI-parallel simulations the energies are not accessible a at every time step.
1591 * We update them if we have new values, otherwise, the energy values from the
1592 * last global communication step are still on display in the viewer. */
1593 if (bHaveNewEnergies
)
1595 ene
->T_abs
= static_cast<float>(impl_
->enerd
->term
[F_TEMP
]);
1596 ene
->E_pot
= static_cast<float>(impl_
->enerd
->term
[F_EPOT
]);
1597 ene
->E_tot
= static_cast<float>(impl_
->enerd
->term
[F_ETOT
]);
1598 ene
->E_bond
= static_cast<float>(impl_
->enerd
->term
[F_BONDS
]);
1599 ene
->E_angle
= static_cast<float>(impl_
->enerd
->term
[F_ANGLES
]);
1600 ene
->E_dihe
= static_cast<float>(impl_
->enerd
->term
[F_PDIHS
]);
1601 ene
->E_impr
= static_cast<float>(impl_
->enerd
->term
[F_IDIHS
]);
1602 ene
->E_vdw
= static_cast<float>(impl_
->enerd
->term
[F_LJ
]);
1603 ene
->E_coul
= static_cast<float>(impl_
->enerd
->term
[F_COUL_SR
]);
1608 void ImdSession::sendPositionsAndEnergies()
1610 if (!impl_
->sessionPossible
|| !impl_
->clientsocket
)
1615 if (imd_send_energies(impl_
->clientsocket
, impl_
->energies
, impl_
->energysendbuf
))
1617 impl_
->issueFatalError("Error sending updated energies. Disconnecting client.");
1620 if (imd_send_rvecs(impl_
->clientsocket
, impl_
->nat
, impl_
->xa
, impl_
->coordsendbuf
))
1622 impl_
->issueFatalError("Error sending updated positions. Disconnecting client.");
1627 void ImdSession::updateEnergyRecordAndSendPositionsAndEnergies(bool bIMDstep
, int64_t step
, bool bHaveNewEnergies
)
1629 if (!impl_
->sessionPossible
)
1634 wallcycle_start(impl_
->wcycle
, ewcIMD
);
1636 /* Update time step for IMD and prepare IMD energy record if we have new energies. */
1637 fillEnergyRecord(step
, bHaveNewEnergies
);
1641 /* Send positions and energies to VMD client via IMD */
1642 sendPositionsAndEnergies();
1645 wallcycle_stop(impl_
->wcycle
, ewcIMD
);
1648 void ImdSession::applyForces(rvec
* f
)
1650 if (!impl_
->sessionPossible
|| !impl_
->bForceActivated
)
1655 wallcycle_start(impl_
->wcycle
, ewcIMD
);
1657 for (int i
= 0; i
< impl_
->nforces
; i
++)
1659 /* j are the indices in the "System group".*/
1660 int j
= impl_
->ind
[impl_
->f_ind
[i
]];
1662 /* check if this is a local atom and find out locndx */
1664 const t_commrec
* cr
= impl_
->cr
;
1665 if (PAR(cr
) && (locndx
= cr
->dd
->ga2la
->findHome(j
)))
1670 rvec_inc(f
[j
], impl_
->f
[i
]);
1673 wallcycle_stop(impl_
->wcycle
, ewcIMD
);
1676 ImdSession::ImdSession(const MDLogger
& mdlog
) : impl_(new Impl(mdlog
)) {}
1678 ImdSession::~ImdSession() = default;