Refactor SD update
[gromacs.git] / src / gromacs / imd / imd.cpp
blobbc242a4f34e601841d9beaf036dd0198e2ee2486
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018, 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.
36 /*! \internal \file
38 * \brief
39 * Implements functions of imd.h.
41 * Re-implementation of basic IMD functions from NAMD/VMD from scratch,
42 * see imdsocket.h for references to the IMD API.
44 * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
46 * \ingroup module_imd
48 #include "gmxpre.h"
50 #include "imd.h"
52 #include "config.h"
54 #include <errno.h>
55 #include <string.h>
57 #if GMX_NATIVE_WINDOWS
58 #include <windows.h>
59 #else
60 #include <unistd.h>
61 #endif
63 #include "gromacs/commandline/filenm.h"
64 #include "gromacs/domdec/domdec_struct.h"
65 #include "gromacs/domdec/ga2la.h"
66 #include "gromacs/fileio/confio.h"
67 #include "gromacs/fileio/gmxfio.h"
68 #include "gromacs/fileio/xvgr.h"
69 #include "gromacs/gmxlib/network.h"
70 #include "gromacs/imd/imdsocket.h"
71 #include "gromacs/math/units.h"
72 #include "gromacs/math/vec.h"
73 #include "gromacs/mdlib/broadcaststructs.h"
74 #include "gromacs/mdlib/groupcoord.h"
75 #include "gromacs/mdlib/mdrun.h"
76 #include "gromacs/mdlib/sighandler.h"
77 #include "gromacs/mdlib/sim_util.h"
78 #include "gromacs/mdtypes/inputrec.h"
79 #include "gromacs/mdtypes/md_enums.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/smalloc.h"
88 /*! \brief How long shall we wait in seconds until we check for a connection again? */
89 #define IMDLOOPWAIT 1
91 /*! \brief How long shall we check for the IMD_GO? */
92 #define IMDCONNECTWAIT 2
94 /*! \brief IMD Header Size. */
95 #define HEADERSIZE 8
96 /*! \brief IMD Protocol Version. */
97 #define IMDVERSION 2
100 /*! \internal
101 * \brief
102 * IMD (interactive molecular dynamics) energy record.
104 * As in the original IMD implementation. Energies in kcal/mol.
105 * NOTE: We return the energies in GROMACS / SI units,
106 * so they also show up as SI in VMD.
109 typedef struct
111 gmx_int32_t tstep; /**< time step */
112 float T_abs; /**< absolute temperature */
113 float E_tot; /**< total energy */
114 float E_pot; /**< potential energy */
115 float E_vdw; /**< van der Waals energy */
116 float E_coul; /**< Coulomb interaction energy */
117 float E_bond; /**< bonds energy */
118 float E_angle; /**< angles energy */
119 float E_dihe; /**< dihedrals energy */
120 float E_impr; /**< improper dihedrals energy */
121 } IMDEnergyBlock;
124 /*! \internal
125 * \brief IMD (interactive molecular dynamics) communication structure.
127 * This structure defines the IMD communication message header & protocol version.
129 typedef struct
131 gmx_int32_t type; /**< Type of IMD message, see IMDType_t above */
132 gmx_int32_t length; /**< Length */
133 } IMDHeader;
136 /*! \internal
137 * \brief IMD (interactive molecular dynamics) main data structure.
139 * Contains private IMD data
141 typedef struct t_gmx_IMD
143 FILE *outf; /**< Output file for IMD data, mainly forces. */
145 int nat; /**< Number of atoms that can be pulled via IMD. */
146 int nat_loc; /**< Part of the atoms that are local. */
147 int *ind; /**< Global indices of the IMD atoms. */
148 int *ind_loc; /**< Local indices of the IMD atoms. */
149 int nalloc_loc; /**< Allocation size for ind_loc. */
150 rvec *xa; /**< Positions for all IMD atoms assembled on
151 the master node. */
152 ivec *xa_shifts; /**< Shifts for all IMD atoms, to make
153 molecule(s) whole. */
154 ivec *xa_eshifts; /**< Extra shifts since last DD step. */
155 rvec *xa_old; /**< Old positions for all IMD atoms on master. */
156 int *xa_ind; /**< Position of each local atom in the
157 collective array. */
159 int nstimd; /**< Global IMD frequency, known to all nodes. */
160 int nstimd_new; /**< New frequency from IMD client, master only. */
161 int nstimd_def; /**< Default IMD frequency when disconnected. */
163 int port; /**< Port to use for network socket. */
164 IMDSocket *socket; /**< The IMD socket on the master node. */
165 IMDSocket *clientsocket; /**< The IMD socket on the client. */
166 int length; /**< Length we got with last header. */
168 gmx_bool bWConnect; /**< Shall we block and wait for connection? */
169 gmx_bool bTerminated; /**< Set if MD is terminated. */
170 gmx_bool bTerminatable; /**< Set if MD can be terminated. */
171 gmx_bool bConnected; /**< Set if connection is present. */
172 gmx_bool bNewForces; /**< Set if we received new forces. */
173 gmx_bool bForceActivated; /**< Set if pulling from VMD is allowed. */
175 IMDEnergyBlock *energies; /**< Pointer to energies we send back. */
177 gmx_int32_t vmd_nforces; /**< Number of VMD forces. */
178 gmx_int32_t *vmd_f_ind; /**< VMD forces indices. */
179 float *vmd_forces; /**< The VMD forces flat in memory. */
180 int nforces; /**< Number of actual MD forces;
181 this gets communicated to the clients. */
182 int *f_ind; /**< Force indices. */
183 rvec *f; /**< The IMD pulling forces. */
185 char *forcesendbuf; /**< Buffer for force sending. */
186 char *coordsendbuf; /**< Buffer for coordinate sending. */
187 char *energysendbuf; /**< Send buffer for energies. */
188 rvec *sendxbuf; /**< Buffer to make molecules whole before
189 sending. */
191 t_block mols; /**< Molecules block in IMD group. */
193 /* The next block is used on the master node only to reduce the output
194 * without sacrificing information. If any of these values changes,
195 * we need to write output */
196 int old_nforces; /**< Old value for nforces. */
197 int *old_f_ind; /**< Old values for force indices. */
198 rvec *old_forces; /**< Old values for IMD pulling forces. */
200 } t_gmx_IMD_setup;
203 /*! \internal
204 * \brief Enum for types of IMD messages.
206 * We use the same records as the NAMD/VMD IMD implementation.
208 typedef enum IMDType_t
210 IMD_DISCONNECT, /**< client disconnect */
211 IMD_ENERGIES, /**< energy data */
212 IMD_FCOORDS, /**< atomic coordinates */
213 IMD_GO, /**< start command for the simulation */
214 IMD_HANDSHAKE, /**< handshake to determine little/big endianness */
215 IMD_KILL, /**< terminates the simulation */
216 IMD_MDCOMM, /**< force data */
217 IMD_PAUSE, /**< pauses the simulation */
218 IMD_TRATE, /**< sets the IMD transmission and processing rate */
219 IMD_IOERROR, /**< I/O error */
220 IMD_NR /**< number of entries */
221 } IMDMessageType;
224 /*! \internal
225 * \brief Names of the IMDType for error messages.
227 const char *eIMDType_names[IMD_NR + 1] = {
228 "IMD_DISCONNECT",
229 "IMD_ENERGIES",
230 "IMD_FCOORDS",
231 "IMD_GO",
232 "IMD_HANDSHAKE",
233 "IMD_KILL",
234 "IMD_MDCOMM",
235 "IMD_PAUSE",
236 "IMD_TRATE",
237 "IMD_IOERROR",
238 nullptr
242 #ifdef GMX_IMD
245 /*! \brief Fills the header with message and the length argument. */
246 static void fill_header(IMDHeader *header, IMDMessageType type, gmx_int32_t length)
248 /* We (ab-)use htonl network function for the correct endianness */
249 header->type = htonl((gmx_int32_t) type);
250 header->length = htonl(length);
254 /*! \brief Swaps the endianess of the header. */
255 static void swap_header(IMDHeader *header)
257 /* and vice versa... */
258 header->type = ntohl(header->type);
259 header->length = ntohl(header->length);
263 /*! \brief Reads multiple bytes from socket. */
264 static gmx_int32_t imd_read_multiple(IMDSocket *socket, char *datptr, gmx_int32_t toread)
266 gmx_int32_t leftcount, countread;
269 leftcount = toread;
270 /* Try to read while we haven't reached toread */
271 while (leftcount != 0)
273 if ((countread = imdsock_read(socket, datptr, leftcount)) < 0)
275 /* interrupted function call, try again... */
276 if (errno == EINTR)
278 countread = 0;
280 /* this is an unexpected error, return what we got */
281 else
283 return toread - leftcount;
286 /* nothing read, finished */
288 else if (countread == 0)
290 break;
292 leftcount -= countread;
293 datptr += countread;
294 } /* end while */
296 /* return nr of bytes read */
297 return toread - leftcount;
301 /*! \brief Writes multiple bytes to socket in analogy to imd_read_multiple. */
302 static gmx_int32_t imd_write_multiple(IMDSocket *socket, const char *datptr, gmx_int32_t towrite)
304 gmx_int32_t leftcount, countwritten;
307 leftcount = towrite;
308 while (leftcount != 0)
310 if ((countwritten = imdsock_write(socket, datptr, leftcount)) <= 0)
312 if (errno == EINTR)
314 countwritten = 0;
316 else
318 return towrite - leftcount;
321 leftcount -= countwritten;
322 datptr += countwritten;
323 } /* end while */
325 return towrite - leftcount;
329 /*! \brief Handshake with IMD client. */
330 static int imd_handshake(IMDSocket *socket)
332 IMDHeader header;
335 fill_header(&header, IMD_HANDSHAKE, 1);
336 header.length = IMDVERSION; /* client wants unswapped version */
338 return (imd_write_multiple(socket, (char *) &header, HEADERSIZE) != HEADERSIZE);
342 /*! \brief Send energies using the energy block and the send buffer. */
343 static int imd_send_energies(IMDSocket *socket, const IMDEnergyBlock *energies, char *buffer)
345 gmx_int32_t recsize;
348 recsize = HEADERSIZE + sizeof(IMDEnergyBlock);
349 fill_header((IMDHeader *) buffer, IMD_ENERGIES, 1);
350 memcpy(buffer + HEADERSIZE, energies, sizeof(IMDEnergyBlock));
352 return (imd_write_multiple(socket, buffer, recsize) != recsize);
356 /*! \brief Receive IMD header from socket, sets the length and returns the IMD message. */
357 static IMDMessageType imd_recv_header(IMDSocket *socket, gmx_int32_t *length)
359 IMDHeader header;
362 if (imd_read_multiple(socket, (char *) &header, HEADERSIZE) != HEADERSIZE)
364 return IMD_IOERROR;
366 swap_header(&header);
367 *length = header.length;
369 return (IMDMessageType) header.type;
373 /*! \brief Receive force indices and forces.
375 * The number of forces was previously communicated via the header.
377 static int imd_recv_mdcomm(IMDSocket *socket, gmx_int32_t nforces, gmx_int32_t *forcendx, float *forces)
379 int retsize, retbytes;
382 /* reading indices */
383 retsize = sizeof(gmx_int32_t) * nforces;
384 retbytes = imd_read_multiple(socket, (char *) forcendx, retsize);
385 if (retbytes != retsize)
387 return FALSE;
390 /* reading forces as float array */
391 retsize = 3 * sizeof(float) * nforces;
392 retbytes = imd_read_multiple(socket, (char *) forces, retsize);
393 if (retbytes != retsize)
395 return FALSE;
398 return TRUE;
401 #endif
403 /* GROMACS specific functions for the IMD implementation */
404 void write_IMDgroup_to_file(gmx_bool bIMD, t_inputrec *ir, t_state *state,
405 gmx_mtop_t *sys, int nfile, const t_filenm fnm[])
407 t_atoms IMDatoms;
410 if (bIMD)
412 IMDatoms = gmx_mtop_global_atoms(sys);
413 write_sto_conf_indexed(opt2fn("-imd", nfile, fnm), "IMDgroup", &IMDatoms,
414 as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), ir->ePBC, state->box, ir->imd->nat, ir->imd->ind);
419 void dd_make_local_IMD_atoms(gmx_bool bIMD, gmx_domdec_t *dd, t_IMD *imd)
421 gmx_ga2la_t *ga2la;
422 t_gmx_IMD_setup *IMDsetup;
424 if (bIMD)
426 IMDsetup = imd->setup;
427 ga2la = dd->ga2la;
429 dd_make_local_group_indices(
430 ga2la, IMDsetup->nat, IMDsetup->ind, &IMDsetup->nat_loc,
431 &IMDsetup->ind_loc, &IMDsetup->nalloc_loc, IMDsetup->xa_ind);
436 #ifdef GMX_IMD
437 /*! \brief Send positions from rvec.
439 * We need a separate send buffer and conversion to Angstrom.
441 static int imd_send_rvecs(IMDSocket *socket, int nat, rvec *x, char *buffer)
443 gmx_int32_t size;
444 int i;
445 float sendx[3];
446 int tuplesize = 3 * sizeof(float);
449 /* Required size for the send buffer */
450 size = HEADERSIZE + 3 * sizeof(float) * nat;
452 /* Prepare header */
453 fill_header((IMDHeader *) buffer, IMD_FCOORDS, (gmx_int32_t) nat);
454 for (i = 0; i < nat; i++)
456 sendx[0] = (float) x[i][0] * NM2A;
457 sendx[1] = (float) x[i][1] * NM2A;
458 sendx[2] = (float) x[i][2] * NM2A;
459 memcpy(buffer + HEADERSIZE + i * tuplesize, sendx, tuplesize);
462 return (imd_write_multiple(socket, buffer, size) != size);
466 /*! \brief Initializes the IMD private data. */
467 static t_gmx_IMD_setup* imd_create(int imdatoms, int nstimddef, int imdport)
469 t_gmx_IMD_setup *IMDsetup = nullptr;
472 snew(IMDsetup, 1);
473 IMDsetup->nat = imdatoms;
474 IMDsetup->bTerminated = FALSE;
475 IMDsetup->bTerminatable = FALSE;
476 IMDsetup->bWConnect = FALSE;
477 IMDsetup->bConnected = FALSE;
478 IMDsetup->bForceActivated = FALSE;
479 IMDsetup->bNewForces = FALSE;
480 IMDsetup->bForceActivated = FALSE;
481 IMDsetup->nstimd = 1;
482 IMDsetup->nstimd_new = 1;
483 IMDsetup->nstimd_def = nstimddef;
484 if (imdport < 1)
486 IMDsetup->port = 0;
487 fprintf(stderr, "%s You chose a port number < 1. Will automatically assign a free port.\n", IMDstr);
489 else
491 IMDsetup->port = imdport;
494 return IMDsetup;
498 /*! \brief Prepare the socket on the MASTER. */
499 static void imd_prepare_master_socket(t_gmx_IMD_setup *IMDsetup)
501 int ret;
504 #if GMX_NATIVE_WINDOWS
505 /* Winsock requires separate initialization */
506 fprintf(stderr, "%s Initializing winsock.\n", IMDstr);
507 #ifdef GMX_HAVE_WINSOCK
508 if (imdsock_winsockinit())
510 gmx_fatal(FARGS, "%s Failed to initialize winsock.\n", IMDstr);
512 #endif
513 #endif
515 /* The rest is identical, first create and bind a socket and set to listen then. */
516 fprintf(stderr, "%s Setting up incoming socket.\n", IMDstr);
517 IMDsetup->socket = imdsock_create();
518 if (!IMDsetup->socket)
520 gmx_fatal(FARGS, "%s Failed to create socket.", IMDstr);
523 /* Bind to port */
524 ret = imdsock_bind(IMDsetup->socket, IMDsetup->port);
525 if (ret)
527 gmx_fatal(FARGS, "%s binding socket to port %d failed with error %d.\n", IMDstr, IMDsetup->port, ret);
530 if (imd_sock_listen(IMDsetup->socket))
532 gmx_fatal(FARGS, "%s socket listen failed with error %d.\n", IMDstr, ret);
535 if (imdsock_getport(IMDsetup->socket, &IMDsetup->port))
537 gmx_fatal(FARGS, "%s Could not determine port number.\n", IMDstr, ret);
540 fprintf(stderr, "%s Listening for IMD connection on port %d.\n", IMDstr, IMDsetup->port);
544 /*! \brief Disconnect the client. */
545 static void imd_disconnect(t_gmx_IMD_setup *IMDsetup)
547 /* Write out any buffered pulling data */
548 fflush(IMDsetup->outf);
550 /* we first try to shut down the clientsocket */
551 imdsock_shutdown(IMDsetup->clientsocket);
552 if (!imdsock_destroy(IMDsetup->clientsocket))
554 fprintf(stderr, "%s Failed to destroy socket.\n", IMDstr);
557 /* then we reset the IMD step to its default, and reset the connection boolean */
558 IMDsetup->nstimd_new = IMDsetup->nstimd_def;
559 IMDsetup->clientsocket = nullptr;
560 IMDsetup->bConnected = FALSE;
564 /*! \brief Prints an error message and disconnects the client.
566 * Does not terminate mdrun!
568 static void imd_fatal(t_gmx_IMD_setup *IMDsetup, const char *msg)
570 fprintf(stderr, "%s %s", IMDstr, msg);
571 imd_disconnect(IMDsetup);
572 fprintf(stderr, "%s disconnected.\n", IMDstr);
576 /*! \brief Check whether we got an incoming connection. */
577 static gmx_bool imd_tryconnect(t_gmx_IMD_setup *IMDsetup)
579 if (imdsock_tryread(IMDsetup->socket, 0, 0) > 0)
581 /* yes, we got something, accept on clientsocket */
582 IMDsetup->clientsocket = imdsock_accept(IMDsetup->socket);
583 if (!IMDsetup->clientsocket)
585 fprintf(stderr, "%s Accepting the connection on the socket failed.\n", IMDstr);
586 return FALSE;
589 /* handshake with client */
590 if (imd_handshake(IMDsetup->clientsocket))
592 imd_fatal(IMDsetup, "Connection failed.\n");
593 return FALSE;
596 fprintf(stderr, "%s Connection established, checking if I got IMD_GO orders.\n", IMDstr);
598 /* Check if we get the proper "GO" command from client. */
599 if (imdsock_tryread(IMDsetup->clientsocket, IMDCONNECTWAIT, 0) != 1 || imd_recv_header(IMDsetup->clientsocket, &(IMDsetup->length)) != IMD_GO)
601 imd_fatal(IMDsetup, "No IMD_GO order received. IMD connection failed.\n");
604 /* IMD connected */
605 IMDsetup->bConnected = TRUE;
607 return TRUE;
610 return FALSE;
614 /*! \brief Wrap imd_tryconnect in order to make it blocking.
616 * Used when the simulation should wait for an incoming connection.
618 static void imd_blockconnect(t_gmx_IMD_setup *IMDsetup)
620 /* do not wait for connection, when e.g. ctrl+c is pressed and we will terminate anyways. */
621 if (!((int) gmx_get_stop_condition() == gmx_stop_cond_none))
623 return;
626 fprintf(stderr, "%s Will wait until I have a connection and IMD_GO orders.\n", IMDstr);
628 /* while we have no clientsocket... 2nd part: we should still react on ctrl+c */
629 while ((!IMDsetup->clientsocket) && ((int) gmx_get_stop_condition() == gmx_stop_cond_none))
631 imd_tryconnect(IMDsetup);
632 #if GMX_NATIVE_WINDOWS
633 /* for whatever reason, it is called Sleep on windows */
634 Sleep(IMDLOOPWAIT);
635 #else
636 sleep(IMDLOOPWAIT);
637 #endif
642 /*! \brief Make sure that our array holding the forces received via IMD is large enough. */
643 static void imd_prepare_vmd_Forces(t_gmx_IMD_setup *IMDsetup)
645 srenew((IMDsetup->vmd_f_ind), IMDsetup->vmd_nforces);
646 srenew((IMDsetup->vmd_forces), 3*IMDsetup->vmd_nforces);
650 /*! \brief Reads forces received via IMD. */
651 static void imd_read_vmd_Forces(t_gmx_IMD_setup *IMDsetup)
653 /* the length of the previously received header tells us the nr of forces we will receive */
654 IMDsetup->vmd_nforces = IMDsetup->length;
655 /* prepare the arrays */
656 imd_prepare_vmd_Forces(IMDsetup);
657 /* Now we read the forces... */
658 if (!(imd_recv_mdcomm(IMDsetup->clientsocket, IMDsetup->vmd_nforces, IMDsetup->vmd_f_ind, IMDsetup->vmd_forces)))
660 imd_fatal(IMDsetup, "Error while reading forces from remote. Disconnecting\n");
665 /*! \brief Prepares the MD force arrays. */
666 static void imd_prepare_MD_Forces(t_gmx_IMD_setup *IMDsetup)
668 srenew((IMDsetup->f_ind), IMDsetup->nforces);
669 srenew((IMDsetup->f ), IMDsetup->nforces);
673 /*! \brief Copy IMD forces to MD forces.
675 * Do conversion from Cal->Joule and from
676 * Angstrom -> nm and from a pointer array to arrays to 3*N array.
678 static void imd_copyto_MD_Forces(t_gmx_IMD_setup *IMDsetup)
680 int i;
681 real conversion = CAL2JOULE * NM2A;
684 for (i = 0; i < IMDsetup->nforces; i++)
686 /* Copy the indices, a copy is important because we may update the incoming forces
687 * whenever we receive new forces while the MD forces are only communicated upon
688 * IMD communication.*/
689 IMDsetup->f_ind[i] = IMDsetup->vmd_f_ind[i];
691 /* Convert to rvecs and do a proper unit conversion */
692 IMDsetup->f[i][0] = IMDsetup->vmd_forces[3*i ] * conversion;
693 IMDsetup->f[i][1] = IMDsetup->vmd_forces[3*i + 1] * conversion;
694 IMDsetup->f[i][2] = IMDsetup->vmd_forces[3*i + 2] * conversion;
699 /*! \brief Return TRUE if any of the forces or indices changed. */
700 static gmx_bool bForcesChanged(t_gmx_IMD_setup *IMDsetup)
702 int i;
705 /* First, check whether the number of pulled atoms changed */
706 if (IMDsetup->nforces != IMDsetup->old_nforces)
708 return TRUE;
711 /* Second, check whether any of the involved atoms changed */
712 for (i = 0; i < IMDsetup->nforces; i++)
714 if (IMDsetup->f_ind[i] != IMDsetup->old_f_ind[i])
716 return TRUE;
720 /* Third, check whether all forces are the same */
721 for (i = 0; i < IMDsetup->nforces; i++)
723 if (IMDsetup->f[i][XX] != IMDsetup->old_forces[i][XX])
725 return TRUE;
727 if (IMDsetup->f[i][YY] != IMDsetup->old_forces[i][YY])
729 return TRUE;
731 if (IMDsetup->f[i][ZZ] != IMDsetup->old_forces[i][ZZ])
733 return TRUE;
737 /* All old and new forces are identical! */
738 return FALSE;
742 /*! \brief Fill the old_f_ind and old_forces arrays with the new, old values. */
743 static void keep_old_values(t_gmx_IMD_setup *IMDsetup)
745 int i;
748 IMDsetup->old_nforces = IMDsetup->nforces;
750 for (i = 0; i < IMDsetup->nforces; i++)
752 IMDsetup->old_f_ind[i] = IMDsetup->f_ind[i];
753 copy_rvec(IMDsetup->f[i], IMDsetup->old_forces[i]);
758 /*! \brief Returns TRUE if any component of the two rvecs differs. */
759 static inline gmx_bool rvecs_differ(const rvec v1, const rvec v2)
761 int i;
764 for (i = 0; i < DIM; i++)
766 if (v1[i] != v2[i])
768 return TRUE;
772 return FALSE;
776 /*! \brief Write the applied pull forces to logfile.
778 * Call on master only!
780 static void output_imd_forces(t_inputrec *ir, double time)
782 t_gmx_IMD_setup *IMDsetup;
783 int i;
786 IMDsetup = ir->imd->setup;
788 if (bForcesChanged(IMDsetup))
790 /* Write time and total number of applied IMD forces */
791 fprintf(IMDsetup->outf, "%14.6e%6d", time, IMDsetup->nforces);
793 /* Write out the global atom indices of the pulled atoms and the forces itself,
794 * write out a force only if it has changed since the last output */
795 for (i = 0; i < IMDsetup->nforces; i++)
797 if (rvecs_differ(IMDsetup->f[i], IMDsetup->old_forces[i]))
799 fprintf(IMDsetup->outf, "%9d", IMDsetup->ind[IMDsetup->f_ind[i]] + 1);
800 fprintf(IMDsetup->outf, "%12.4e%12.4e%12.4e", IMDsetup->f[i][0], IMDsetup->f[i][1], IMDsetup->f[i][2]);
803 fprintf(IMDsetup->outf, "\n");
805 keep_old_values(IMDsetup);
810 /*! \brief Synchronize the nodes. */
811 static void imd_sync_nodes(t_inputrec *ir, const t_commrec *cr, double t)
813 int new_nforces = 0;
814 t_gmx_IMD_setup *IMDsetup;
817 IMDsetup = ir->imd->setup;
819 /* Notify the other nodes whether we are still connected. */
820 if (PAR(cr))
822 block_bc(cr, IMDsetup->bConnected);
825 /* ...if not connected, the job is done here. */
826 if (!IMDsetup->bConnected)
828 return;
831 /* Let the other nodes know whether we got a new IMD synchronization frequency. */
832 if (PAR(cr))
834 block_bc(cr, IMDsetup->nstimd_new);
837 /* Now we all set the (new) nstimd communication time step */
838 IMDsetup->nstimd = IMDsetup->nstimd_new;
840 /* We're done if we don't allow pulling at all */
841 if (!(IMDsetup->bForceActivated))
843 return;
846 /* OK, let's check if we have received forces which we need to communicate
847 * to the other nodes */
848 if (MASTER(cr))
850 if (IMDsetup->bNewForces)
852 new_nforces = IMDsetup->vmd_nforces;
854 /* make the "new_forces" negative, when we did not receive new ones */
855 else
857 new_nforces = IMDsetup->vmd_nforces * -1;
861 /* make new_forces known to the clients */
862 if (PAR(cr))
864 block_bc(cr, new_nforces);
867 /* When new_natoms < 0 then we know that these are still the same forces
868 * so we don't communicate them, otherwise... */
869 if (new_nforces >= 0)
871 /* set local VMD and nforces */
872 IMDsetup->vmd_nforces = new_nforces;
873 IMDsetup->nforces = IMDsetup->vmd_nforces;
875 /* now everybody knows the number of forces in f_ind, so we can prepare
876 * the target arrays for indices and forces */
877 imd_prepare_MD_Forces(IMDsetup);
879 /* we first update the MD forces on the master by converting the VMD forces */
880 if (MASTER(cr))
882 imd_copyto_MD_Forces(IMDsetup);
883 /* We also write out forces on every update, so that we know which
884 * forces are applied for every step */
885 if (IMDsetup->outf)
887 output_imd_forces(ir, t);
891 /* In parallel mode we communicate the to-be-applied forces to the other nodes */
892 if (PAR(cr))
894 nblock_bc(cr, IMDsetup->nforces, IMDsetup->f_ind);
895 nblock_bc(cr, IMDsetup->nforces, IMDsetup->f );
898 /* done communicating the forces, reset bNewForces */
899 IMDsetup->bNewForces = FALSE;
904 /*! \brief Reads header from the client and decides what to do. */
905 static void imd_readcommand(t_gmx_IMD_setup *IMDsetup)
907 gmx_bool IMDpaused = FALSE;
908 IMDMessageType itype;
911 while (IMDsetup->clientsocket && (imdsock_tryread(IMDsetup->clientsocket, 0, 0) > 0 || IMDpaused))
913 itype = imd_recv_header(IMDsetup->clientsocket, &(IMDsetup->length));
914 /* let's see what we got: */
915 switch (itype)
917 /* IMD asks us to terminate the simulation, check if the user allowed this */
918 case IMD_KILL:
919 if (IMDsetup->bTerminatable)
921 fprintf(stderr, " %s Terminating connection and running simulation (if supported by integrator).\n", IMDstr);
922 IMDsetup->bTerminated = TRUE;
923 IMDsetup->bWConnect = FALSE;
924 gmx_set_stop_condition(gmx_stop_cond_next);
926 else
928 fprintf(stderr, " %s Set -imdterm command line switch to allow mdrun termination from within IMD.\n", IMDstr);
931 break;
933 /* the client doen't want to talk to us anymore */
934 case IMD_DISCONNECT:
935 fprintf(stderr, " %s Disconnecting client.\n", IMDstr);
936 imd_disconnect(IMDsetup);
937 break;
939 /* we got new forces, read them and set bNewForces flag */
940 case IMD_MDCOMM:
941 imd_read_vmd_Forces(IMDsetup);
942 IMDsetup->bNewForces = TRUE;
943 break;
945 /* the client asks us to (un)pause the simulation. So we toggle the IMDpaused state */
946 case IMD_PAUSE:
947 if (IMDpaused)
949 fprintf(stderr, " %s Un-pause command received.\n", IMDstr);
950 IMDpaused = FALSE;
952 else
954 fprintf(stderr, " %s Pause command received.\n", IMDstr);
955 IMDpaused = TRUE;
958 break;
960 /* the client sets a new transfer rate, if we get 0, we reset the rate
961 * to the default. VMD filters 0 however */
962 case IMD_TRATE:
963 IMDsetup->nstimd_new = (IMDsetup->length > 0) ? IMDsetup->length : IMDsetup->nstimd_def;
964 fprintf(stderr, " %s Update frequency will be set to %d.\n", IMDstr, IMDsetup->nstimd_new);
965 break;
967 /* Catch all rule for the remaining IMD types which we don't expect */
968 default:
969 fprintf(stderr, " %s Received unexpected %s.\n", IMDstr, enum_name((int)itype, IMD_NR, eIMDType_names));
970 imd_fatal(IMDsetup, "Terminating connection\n");
971 break;
972 } /* end switch */
973 } /* end while */
977 /*! \brief Open IMD output file and write header information.
979 * Call on master only.
981 static FILE *open_imd_out(const char *fn,
982 t_gmx_IMD_setup *IMDsetup,
983 int nat_total,
984 const gmx_output_env_t *oenv,
985 const ContinuationOptions &continuationOptions)
987 FILE *fp;
990 /* Open log file of applied IMD forces if requested */
991 if (fn && oenv)
993 /* If we append to an existing file, all the header information is already there */
994 if (continuationOptions.appendFiles)
996 fp = gmx_fio_fopen(fn, "a+");
998 else
1000 fp = gmx_fio_fopen(fn, "w+");
1001 if (IMDsetup->nat == nat_total)
1003 fprintf(fp, "# Note that you can select an IMD index group in the .mdp file if a subset of the atoms suffices.\n");
1006 xvgr_header(fp, "IMD Pull Forces", "Time (ps)", "# of Forces / Atom IDs / Forces (kJ/mol)", exvggtNONE, oenv);
1008 fprintf(fp, "# Can display and manipulate %d (of a total of %d) atoms via IMD.\n", IMDsetup->nat, nat_total);
1009 fprintf(fp, "# column 1 : time (ps)\n");
1010 fprintf(fp, "# column 2 : total number of atoms feeling an IMD pulling force at that time\n");
1011 fprintf(fp, "# cols. 3.-6 : global atom number of pulled atom, x-force, y-force, z-force (kJ/mol)\n");
1012 fprintf(fp, "# then follow : atom-ID, f[x], f[y], f[z] for more atoms in case the force on multiple atoms is changed simultaneously.\n");
1013 fprintf(fp, "# Note that the force on any atom is always equal to the last value for that atom-ID found in the data.\n");
1014 fflush(fp);
1017 /* To reduce the output file size we remember the old values and output only
1018 * when something changed */
1019 snew(IMDsetup->old_f_ind, IMDsetup->nat); /* One can never pull on more atoms */
1020 snew(IMDsetup->old_forces, IMDsetup->nat);
1022 return fp;
1025 fprintf(stdout, "%s For a log of the IMD pull forces explicitly specify '-if' on the command line.\n"
1026 "%s (Not possible with energy minimization.)\n", IMDstr, IMDstr);
1028 return nullptr;
1030 #endif
1033 void IMD_finalize(gmx_bool bIMD, t_IMD *imd)
1035 if (bIMD)
1037 if (imd->setup->outf)
1039 gmx_fio_fclose(imd->setup->outf);
1045 #ifdef GMX_IMD
1046 /*! \brief Creates the molecule start-end position array of molecules in the IMD group. */
1047 static void init_imd_prepare_mols_in_imdgroup(t_gmx_IMD_setup *IMDsetup, gmx_mtop_t *top_global)
1049 int i, ii;
1050 int gstart, gend, count;
1051 t_block lmols;
1052 int nat;
1053 int *ind;
1055 nat = IMDsetup->nat;
1056 ind = IMDsetup->ind;
1058 lmols.nr = 0;
1060 /* check whether index is sorted */
1061 for (i = 0; i < nat-1; i++)
1063 if (ind[i] > ind[i+1])
1065 gmx_fatal(FARGS, "%s IMD index is not sorted. This is currently not supported.\n", IMDstr);
1069 gmx::BlockRanges gmols = gmx_mtop_molecules(*top_global);
1070 snew(lmols.index, gmols.numBlocks() + 1);
1071 lmols.index[0] = 0;
1073 for (i = 0; i < gmols.numBlocks(); i++)
1075 gstart = gmols.index[i];
1076 gend = gmols.index[i+1];
1077 count = 0;
1078 for (ii = 0; ii < nat; ii++)
1080 if ((ind[ii] >= gstart) && (ind[ii] < gend))
1082 count += 1;
1085 if (count > 0)
1087 lmols.index[lmols.nr+1] = lmols.index[lmols.nr]+count;
1088 lmols.nr += 1;
1091 srenew(lmols.index, lmols.nr+1);
1092 lmols.nalloc_index = lmols.nr+1;
1093 IMDsetup->mols = lmols;
1097 /*! \brief Copied and modified from groupcoord.c shift_positions_group(). */
1098 static void shift_positions(
1099 matrix box,
1100 rvec x[], /* The positions [0..nr] */
1101 ivec is, /* The shift [0..nr] */
1102 int nr) /* The number of positions */
1104 int i, tx, ty, tz;
1106 /* Loop over the group's atoms */
1107 if (TRICLINIC(box))
1109 for (i = 0; i < nr; i++)
1111 tx = is[XX];
1112 ty = is[YY];
1113 tz = is[ZZ];
1115 x[i][XX] = x[i][XX]-tx*box[XX][XX]-ty*box[YY][XX]-tz*box[ZZ][XX];
1116 x[i][YY] = x[i][YY]-ty*box[YY][YY]-tz*box[ZZ][YY];
1117 x[i][ZZ] = x[i][ZZ]-tz*box[ZZ][ZZ];
1120 else
1122 for (i = 0; i < nr; i++)
1124 tx = is[XX];
1125 ty = is[YY];
1126 tz = is[ZZ];
1128 x[i][XX] = x[i][XX]-tx*box[XX][XX];
1129 x[i][YY] = x[i][YY]-ty*box[YY][YY];
1130 x[i][ZZ] = x[i][ZZ]-tz*box[ZZ][ZZ];
1136 /*! \brief Removes shifts of molecules diffused outside of the box. */
1137 static void imd_remove_molshifts(t_gmx_IMD_setup *IMDsetup, matrix box)
1139 int i, ii, molsize;
1140 ivec largest, smallest, shift;
1141 t_block mols;
1144 mols = IMDsetup->mols;
1146 /* for each molecule also present in IMD group */
1147 for (i = 0; i < mols.nr; i++)
1149 /* first we determine the minimum and maximum shifts for each molecule */
1151 clear_ivec(largest);
1152 clear_ivec(smallest);
1153 clear_ivec(shift);
1155 copy_ivec(IMDsetup->xa_shifts[mols.index[i]], largest);
1156 copy_ivec(IMDsetup->xa_shifts[mols.index[i]], smallest);
1158 for (ii = mols.index[i]+1; ii < mols.index[i+1]; ii++)
1160 if (IMDsetup->xa_shifts[ii][XX] > largest[XX])
1162 largest[XX] = IMDsetup->xa_shifts[ii][XX];
1164 if (IMDsetup->xa_shifts[ii][XX] < smallest[XX])
1166 smallest[XX] = IMDsetup->xa_shifts[ii][XX];
1169 if (IMDsetup->xa_shifts[ii][YY] > largest[YY])
1171 largest[YY] = IMDsetup->xa_shifts[ii][YY];
1173 if (IMDsetup->xa_shifts[ii][YY] < smallest[YY])
1175 smallest[YY] = IMDsetup->xa_shifts[ii][YY];
1178 if (IMDsetup->xa_shifts[ii][ZZ] > largest[ZZ])
1180 largest[ZZ] = IMDsetup->xa_shifts[ii][ZZ];
1182 if (IMDsetup->xa_shifts[ii][ZZ] < smallest[ZZ])
1184 smallest[ZZ] = IMDsetup->xa_shifts[ii][ZZ];
1189 /* check if we what we can subtract/add to the positions
1190 * to put them back into the central box */
1191 if (smallest[XX] > 0)
1193 shift[XX] = smallest[XX];
1195 if (smallest[YY] > 0)
1197 shift[YY] = smallest[YY];
1199 if (smallest[ZZ] > 0)
1201 shift[ZZ] = smallest[ZZ];
1204 if (largest[XX] < 0)
1206 shift[XX] = largest[XX];
1208 if (largest[YY] < 0)
1210 shift[YY] = largest[YY];
1212 if (largest[ZZ] < 0)
1214 shift[ZZ] = largest[ZZ];
1217 /* is there a shift at all? */
1218 if ((shift[XX]) || (shift[YY]) || (shift[ZZ]))
1220 molsize = mols.index[i+1]-mols.index[i];
1221 /* shift the positions */
1222 shift_positions(box, &(IMDsetup->xa[mols.index[i]]), shift, molsize);
1229 /*! \brief Initialize arrays used to assemble the positions from the other nodes. */
1230 static void init_imd_prepare_for_x_assembly(const t_commrec *cr, rvec x[], t_gmx_IMD_setup *IMDsetup)
1232 int i, ii;
1235 snew(IMDsetup->xa, IMDsetup->nat);
1236 snew(IMDsetup->xa_ind, IMDsetup->nat);
1237 snew(IMDsetup->xa_shifts, IMDsetup->nat);
1238 snew(IMDsetup->xa_eshifts, IMDsetup->nat);
1239 snew(IMDsetup->xa_old, IMDsetup->nat);
1241 /* Save the original (whole) set of positions such that later the
1242 * molecule can always be made whole again */
1243 if (MASTER(cr))
1245 for (i = 0; i < IMDsetup->nat; i++)
1247 ii = IMDsetup->ind[i];
1248 copy_rvec(x[ii], IMDsetup->xa_old[i]);
1252 if (!PAR(cr))
1254 IMDsetup->nat_loc = IMDsetup->nat;
1255 IMDsetup->ind_loc = IMDsetup->ind;
1257 /* xa_ind[i] needs to be set to i for serial runs */
1258 for (i = 0; i < IMDsetup->nat; i++)
1260 IMDsetup->xa_ind[i] = i;
1264 /* Communicate initial coordinates xa_old to all processes */
1265 #if GMX_MPI
1266 if (PAR(cr))
1268 gmx_bcast(IMDsetup->nat * sizeof(IMDsetup->xa_old[0]), IMDsetup->xa_old, cr);
1270 #endif
1272 #endif
1275 /*! \brief Check for non-working integrator / parallel options. */
1276 static void imd_check_integrator_parallel(t_inputrec *ir, const t_commrec *cr)
1278 if (PAR(cr))
1280 if (((ir->eI) == eiSteep) || ((ir->eI) == eiCG) || ((ir->eI) == eiLBFGS) || ((ir->eI) == eiNM))
1282 gmx_fatal(FARGS, "%s Energy minimization via steep, CG, lbfgs and nm in parallel is currently not supported by IMD.\n", IMDstr);
1283 return;
1288 void init_IMD(t_inputrec *ir,
1289 const t_commrec *cr,
1290 const gmx_multisim_t *ms,
1291 gmx_mtop_t *top_global,
1292 FILE *fplog,
1293 int defnstimd,
1294 rvec x[],
1295 int nfile,
1296 const t_filenm fnm[],
1297 const gmx_output_env_t *oenv,
1298 const MdrunOptions &mdrunOptions)
1300 int i;
1301 int nat_total;
1302 t_gmx_IMD_setup *IMDsetup;
1303 gmx_int32_t bufxsize;
1304 gmx_bool bIMD = FALSE;
1307 /* We will allow IMD sessions only if explicitly enabled in the .tpr file */
1308 if (FALSE == ir->bIMD)
1310 return;
1312 // TODO many of these error conditions were we can't do what the
1313 // user asked for should be handled with a fatal error, not just a
1314 // warning.
1316 const ImdOptions &options = mdrunOptions.imdOptions;
1318 /* It seems we have a .tpr file that defines an IMD group and thus allows IMD sessions.
1319 * Check whether we can actually provide the IMD functionality for this setting: */
1320 if (MASTER(cr))
1322 /* Check whether IMD was enabled by one of the command line switches: */
1323 if (options.wait || options.terminatable || options.pull)
1325 /* Multiple simulations or replica exchange */
1326 if (isMultiSim(ms))
1328 fprintf(stderr, "%s Cannot use IMD for multiple simulations or replica exchange.\n", IMDstr);
1330 /* OK, IMD seems to be allowed and turned on... */
1331 else
1333 fprintf(stderr, "%s Enabled. This simulation will accept incoming IMD connections.\n", IMDstr);
1334 bIMD = TRUE;
1337 else
1339 fprintf(stderr, "%s None of the -imd switches was used.\n"
1340 "%s This run will not accept incoming IMD connections\n", IMDstr, IMDstr);
1342 } /* end master only */
1344 /* Disable IMD if not all the needed functionality is there! */
1345 #if GMX_NATIVE_WINDOWS && !defined(GMX_HAVE_WINSOCK)
1346 bIMD = FALSE;
1347 fprintf(stderr, "Disabling IMD because the winsock library was not found at compile time.\n");
1348 #endif
1350 /* Let the other nodes know whether we want IMD */
1351 if (PAR(cr))
1353 block_bc(cr, bIMD);
1355 /* ... and update our local inputrec accordingly. */
1356 ir->bIMD = bIMD;
1358 /*... if not we are done.*/
1359 if (!ir->bIMD)
1361 return;
1365 /* check if we're using a sane integrator / parallel combination */
1366 imd_check_integrator_parallel(ir, cr);
1370 *****************************************************
1371 * From here on we assume that IMD is turned on *
1372 *****************************************************
1375 #ifdef GMX_IMD
1376 nat_total = top_global->natoms;
1378 /* Initialize IMD setup structure. If we read in a pre-IMD .tpr file, imd->nat
1379 * will be zero. For those cases we transfer _all_ atomic positions */
1380 ir->imd->setup = imd_create(ir->imd->nat > 0 ? ir->imd->nat : nat_total,
1381 defnstimd, options.port);
1382 IMDsetup = ir->imd->setup;
1384 /* We might need to open an output file for IMD forces data */
1385 if (MASTER(cr))
1387 IMDsetup->outf = open_imd_out(opt2fn("-if", nfile, fnm), ir->imd->setup, nat_total, oenv, mdrunOptions.continuationOptions);
1390 /* Make sure that we operate with a valid atom index array for the IMD atoms */
1391 if (ir->imd->nat > 0)
1393 /* Point to the user-supplied array of atom numbers */
1394 IMDsetup->ind = ir->imd->ind;
1396 else
1398 /* Make a dummy (ind[i] = i) array of all atoms */
1399 snew(IMDsetup->ind, nat_total);
1400 for (i = 0; i < nat_total; i++)
1402 IMDsetup->ind[i] = i;
1406 /* read environment on master and prepare socket for incoming connections */
1407 if (MASTER(cr))
1409 /* we allocate memory for our IMD energy structure */
1410 gmx_int32_t recsize = HEADERSIZE + sizeof(IMDEnergyBlock);
1411 snew(IMDsetup->energysendbuf, recsize);
1413 /* Shall we wait for a connection? */
1414 if (options.wait)
1416 IMDsetup->bWConnect = TRUE;
1417 fprintf(stderr, "%s Pausing simulation while no IMD connection present (-imdwait).\n", IMDstr);
1420 /* Will the IMD clients be able to terminate the simulation? */
1421 if (options.terminatable)
1423 IMDsetup->bTerminatable = TRUE;
1424 fprintf(stderr, "%s Allow termination of the simulation from IMD client (-imdterm).\n", IMDstr);
1427 /* Is pulling from IMD client allowed? */
1428 if (options.pull)
1430 IMDsetup->bForceActivated = TRUE;
1431 fprintf(stderr, "%s Pulling from IMD remote is enabled (-imdpull).\n", IMDstr);
1434 /* Initialize send buffers with constant size */
1435 snew(IMDsetup->sendxbuf, IMDsetup->nat);
1436 snew(IMDsetup->energies, 1);
1437 bufxsize = HEADERSIZE + 3 * sizeof(float) * IMDsetup->nat;
1438 snew(IMDsetup->coordsendbuf, bufxsize);
1441 /* do we allow interactive pulling? If so let the other nodes know. */
1442 if (PAR(cr))
1444 block_bc(cr, IMDsetup->bForceActivated);
1447 /* setup the listening socket on master process */
1448 if (MASTER(cr))
1450 fprintf(fplog, "%s Setting port for connection requests to %d.\n", IMDstr, IMDsetup->port);
1451 fprintf(stderr, "%s Turning on IMD - port for incoming requests is %d.\n", IMDstr, IMDsetup->port);
1452 imd_prepare_master_socket(IMDsetup);
1453 /* Wait until we have a connection if specified before */
1454 if (IMDsetup->bWConnect)
1456 imd_blockconnect(IMDsetup);
1458 else
1460 fprintf(stderr, "%s -imdwait not set, starting simulation.\n", IMDstr);
1463 /* Let the other nodes know whether we are connected */
1464 imd_sync_nodes(ir, cr, 0);
1466 /* Initialize arrays used to assemble the positions from the other nodes */
1467 init_imd_prepare_for_x_assembly(cr, x, IMDsetup);
1469 /* Initialize molecule blocks to make them whole later...*/
1470 if (MASTER(cr))
1472 init_imd_prepare_mols_in_imdgroup(IMDsetup, top_global);
1474 #else
1475 gmx_incons("init_IMD: this GROMACS version was not compiled with IMD support!");
1476 #endif
1480 gmx_bool do_IMD(gmx_bool bIMD,
1481 gmx_int64_t step,
1482 const t_commrec *cr,
1483 gmx_bool bNS,
1484 matrix box,
1485 rvec x[],
1486 t_inputrec *ir,
1487 double t,
1488 gmx_wallcycle *wcycle)
1490 gmx_bool imdstep = FALSE;
1491 t_gmx_IMD_setup *IMDsetup;
1494 /* IMD at all? */
1495 if (!bIMD)
1497 return FALSE;
1500 #ifdef GMX_IMD
1501 wallcycle_start(wcycle, ewcIMD);
1503 IMDsetup = ir->imd->setup;
1505 /* read command from client and check if new incoming connection */
1506 if (MASTER(cr))
1508 /* If not already connected, check for new connections */
1509 if (!IMDsetup->clientsocket)
1511 if (IMDsetup->bWConnect)
1513 imd_blockconnect(IMDsetup);
1515 else
1517 imd_tryconnect(IMDsetup);
1521 /* Let's see if we have new IMD messages for us */
1522 if (IMDsetup->clientsocket)
1524 imd_readcommand(IMDsetup);
1528 /* is this an IMD communication step? */
1529 imdstep = do_per_step(step, IMDsetup->nstimd);
1531 /* OK so this is an IMD step ... */
1532 if (imdstep)
1534 /* First we sync all nodes to let everybody know whether we are connected to VMD */
1535 imd_sync_nodes(ir, cr, t);
1538 /* If a client is connected, we collect the positions
1539 * and put molecules back into the box before transfer */
1540 if ((imdstep && IMDsetup->bConnected)
1541 || bNS) /* independent of imdstep, we communicate positions at each NS step */
1543 /* Transfer the IMD positions to the master node. Every node contributes
1544 * its local positions x and stores them in the assembled xa array. */
1545 communicate_group_positions(cr, IMDsetup->xa, IMDsetup->xa_shifts, IMDsetup->xa_eshifts,
1546 TRUE, x, IMDsetup->nat, IMDsetup->nat_loc,
1547 IMDsetup->ind_loc, IMDsetup->xa_ind, IMDsetup->xa_old, box);
1549 /* If connected and master -> remove shifts */
1550 if ((imdstep && IMDsetup->bConnected) && MASTER(cr))
1552 imd_remove_molshifts(IMDsetup, box);
1556 wallcycle_stop(wcycle, ewcIMD);
1557 #else
1558 gmx_incons("do_IMD called without IMD support!");
1559 #endif
1561 return imdstep;
1565 void IMD_fill_energy_record(gmx_bool bIMD, t_IMD *imd, gmx_enerdata_t *enerd,
1566 gmx_int64_t step, gmx_bool bHaveNewEnergies)
1568 IMDEnergyBlock *ene;
1569 t_gmx_IMD *IMDsetup;
1572 if (bIMD)
1574 #ifdef GMX_IMD
1575 IMDsetup = imd->setup;
1577 if (IMDsetup->clientsocket)
1579 ene = IMDsetup->energies;
1581 ene->tstep = step;
1583 /* In MPI-parallel simulations the energies are not accessible a at every time step.
1584 * We update them if we have new values, otherwise, the energy values from the
1585 * last global communication step are still on display in the viewer. */
1586 if (bHaveNewEnergies)
1588 ene->T_abs = (float) enerd->term[F_TEMP ];
1589 ene->E_pot = (float) enerd->term[F_EPOT ];
1590 ene->E_tot = (float) enerd->term[F_ETOT ];
1591 ene->E_bond = (float) enerd->term[F_BONDS ];
1592 ene->E_angle = (float) enerd->term[F_ANGLES ];
1593 ene->E_dihe = (float) enerd->term[F_PDIHS ];
1594 ene->E_impr = (float) enerd->term[F_IDIHS ];
1595 ene->E_vdw = (float) enerd->term[F_LJ ];
1596 ene->E_coul = (float) enerd->term[F_COUL_SR];
1599 #else
1600 gmx_incons("IMD_fill_energy_record called without IMD support.");
1601 #endif
1606 void IMD_send_positions(t_IMD *imd)
1608 #ifdef GMX_IMD
1609 t_gmx_IMD *IMDsetup;
1612 IMDsetup = imd->setup;
1614 if (IMDsetup->clientsocket)
1617 if (imd_send_energies(IMDsetup->clientsocket, IMDsetup->energies, IMDsetup->energysendbuf))
1619 imd_fatal(IMDsetup, "Error sending updated energies. Disconnecting client.\n");
1622 if (imd_send_rvecs(IMDsetup->clientsocket, IMDsetup->nat, IMDsetup->xa, IMDsetup->coordsendbuf))
1624 imd_fatal(IMDsetup, "Error sending updated positions. Disconnecting client.\n");
1627 #else
1628 gmx_incons("IMD_send_positions called without IMD support.");
1629 #endif
1633 void IMD_prep_energies_send_positions(gmx_bool bIMD, gmx_bool bIMDstep,
1634 t_IMD *imd, gmx_enerdata_t *enerd,
1635 gmx_int64_t step, gmx_bool bHaveNewEnergies,
1636 gmx_wallcycle *wcycle)
1638 if (bIMD)
1640 #ifdef GMX_IMD
1641 wallcycle_start(wcycle, ewcIMD);
1643 /* Update time step for IMD and prepare IMD energy record if we have new energies. */
1644 IMD_fill_energy_record(TRUE, imd, enerd, step, bHaveNewEnergies);
1646 if (bIMDstep)
1648 /* Send positions and energies to VMD client via IMD */
1649 IMD_send_positions(imd);
1652 wallcycle_stop(wcycle, ewcIMD);
1653 #else
1654 gmx_incons("IMD_prep_energies_send_positions called without IMD support.");
1655 #endif
1659 int IMD_get_step(t_gmx_IMD *IMDsetup)
1661 return IMDsetup->nstimd;
1665 void IMD_apply_forces(gmx_bool bIMD, t_IMD *imd, const t_commrec *cr, rvec *f,
1666 gmx_wallcycle *wcycle)
1668 int i, j;
1669 int locndx;
1670 t_gmx_IMD_setup *IMDsetup;
1673 if (bIMD)
1675 #ifdef GMX_IMD
1676 wallcycle_start(wcycle, ewcIMD);
1678 IMDsetup = imd->setup;
1680 /* Are forces allowed at all? If not we're done */
1681 if (!IMDsetup->bForceActivated)
1683 return;
1686 for (i = 0; i < IMDsetup->nforces; i++)
1688 /* j are the indices in the "System group".*/
1689 j = IMDsetup->ind[IMDsetup->f_ind[i]];
1691 /* check if this is a local atom and find out locndx */
1692 if (PAR(cr) && ga2la_get_home(cr->dd->ga2la, j, &locndx))
1694 j = locndx;
1697 rvec_inc(f[j], IMDsetup->f[i]);
1700 wallcycle_start(wcycle, ewcIMD);
1701 #else
1702 gmx_incons("IMD_apply_forces called without IMD support.");
1703 #endif