Make ImdSession into a Pimpl-ed class with factory function
[alexxy/gromacs.git] / src / gromacs / imd / imd.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief
39  * Implements Interactive Molecular Dynamics
40  *
41  * Re-implementation of basic IMD functions to work with VMD,
42  * see imdsocket.h for references to the IMD API.
43  *
44  * \author Martin Hoefling, Carsten Kutzner <ckutzne@gwdg.de>
45  *
46  * \ingroup module_imd
47  */
48 #include "gmxpre.h"
49
50 #include "imd.h"
51
52 #include "config.h"
53
54 #include <cerrno>
55 #include <cstring>
56
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/mdtypes/enerdata.h"
72 #include "gromacs/mdtypes/inputrec.h"
73 #include "gromacs/mdtypes/md_enums.h"
74 #include "gromacs/mdtypes/mdrunoptions.h"
75 #include "gromacs/mdtypes/state.h"
76 #include "gromacs/pbcutil/pbc.h"
77 #include "gromacs/timing/wallcycle.h"
78 #include "gromacs/topology/mtop_util.h"
79 #include "gromacs/topology/topology.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/logger.h"
82 #include "gromacs/utility/smalloc.h"
83 #include "gromacs/utility/stringutil.h"
84
85 namespace gmx
86 {
87
88 /*! \brief How long shall we wait in seconds until we check for a connection again? */
89 constexpr int c_loopWait = 1;
90
91 /*! \brief How long shall we check for the IMD_GO? */
92 constexpr int c_connectWait = 1;
93
94 /*! \brief IMD Header Size. */
95 constexpr int c_headerSize = 8;
96
97 /*! \brief IMD Protocol Version. */
98 constexpr int c_protocolVersion = 2;
99
100
101 /*! \internal
102  * \brief
103  * IMD (interactive molecular dynamics) energy record.
104  *
105  * As in the original IMD implementation. Energies in kcal/mol.
106  * NOTE: We return the energies in GROMACS / SI units,
107  * so they also show up as SI in VMD.
108  *
109  */
110 typedef struct
111 {
112     int32_t     tstep;   /**< time step                                     */
113     float       T_abs;   /**< absolute temperature                          */
114     float       E_tot;   /**< total energy                                  */
115     float       E_pot;   /**< potential energy                              */
116     float       E_vdw;   /**< van der Waals energy                          */
117     float       E_coul;  /**< Coulomb interaction energy                    */
118     float       E_bond;  /**< bonds energy                                  */
119     float       E_angle; /**< angles energy                                 */
120     float       E_dihe;  /**< dihedrals energy                              */
121     float       E_impr;  /**< improper dihedrals energy                     */
122 } IMDEnergyBlock;
123
124
125 /*! \internal
126  * \brief IMD (interactive molecular dynamics) communication structure.
127  *
128  * This structure defines the IMD communication message header & protocol version.
129  */
130 typedef struct
131 {
132     int32_t type;      /**< Type of IMD message, see IMDType_t above      */
133     int32_t length;    /**< Length                                        */
134 } IMDHeader;
135
136
137 /*! \internal
138  * \brief Implementation type for the IMD session
139  *
140  * \todo Make this class (or one extracted from it) model
141  * IForceProvider.
142  * \todo Use RAII for files and allocations
143  */
144 class ImdSession::Impl
145 {
146     public:
147         //! Constructor
148         Impl(const MDLogger &mdlog);
149         ~Impl();
150
151         /*! \brief Prepare the socket on the MASTER. */
152         void prepareMasterSocket();
153         /*! \brief Disconnect the client. */
154         void disconnectClient();
155         /*! \brief Prints an error message and disconnects the client.
156          *
157          *  Does not terminate mdrun!
158          */
159         void issueFatalError(const char *msg);
160         /*! \brief Check whether we got an incoming connection. */
161         bool tryConnect();
162         /*! \brief Wrap imd_tryconnect in order to make it blocking.
163          *
164          * Used when the simulation should wait for an incoming connection.
165          */
166         void blockConnect();
167         /*! \brief Make sure that our array holding the forces received via IMD is large enough. */
168         void prepareVmdForces();
169         /*! \brief Reads forces received via IMD. */
170         void readVmdForces();
171         /*! \brief Prepares the MD force arrays. */
172         void prepareMDForces();
173         /*! \brief Copy IMD forces to MD forces.
174          *
175          * Do conversion from Cal->Joule and from
176          * Angstrom -> nm and from a pointer array to arrays to 3*N array.
177          */
178         void copyToMDForces();
179         /*! \brief Return true if any of the forces or indices changed. */
180         bool bForcesChanged() const;
181         /*! \brief Update the old_f_ind and old_forces arrays to contain the current values. */
182         void keepOldValues();
183         /*! \brief Write the applied pull forces to logfile.
184          *
185          * Call on master only!
186          */
187         void outputForces(double time);
188         /*! \brief Synchronize the nodes. */
189         void syncNodes(const t_commrec *cr, double t);
190         /*! \brief Reads header from the client and decides what to do. */
191         void readCommand();
192         /*! \brief Open IMD output file and write header information.
193          *
194          * Call on master only.
195          */
196         void openOutputFile(const char                *fn,
197                             int                        nat_total,
198                             const gmx_output_env_t    *oenv,
199                             const ContinuationOptions &continuationOptions);
200         /*! \brief Creates the molecule start-end position array of molecules in the IMD group. */
201         void prepareMoleculesInImdGroup(const gmx_mtop_t *top_global);
202         /*! \brief Removes shifts of molecules diffused outside of the box. */
203         void removeMolecularShifts(const matrix box);
204         /*! \brief Initialize arrays used to assemble the positions from the other nodes. */
205         void prepareForPositionAssembly(const t_commrec *cr, const rvec x[]);
206         /*! \brief Interact with any connected VMD session */
207         bool run(int64_t      step,
208                  bool         bNS,
209                  const matrix box,
210                  const rvec   x[],
211                  double       t);
212
213         // TODO rename all the data members to have underscore suffixes
214
215         //! True if tpr and mdrun input combine to permit IMD sessions
216         bool       sessionPossible = false;
217         //! Output file for IMD data, mainly forces.
218         FILE      *outf = nullptr;
219
220         //! Number of atoms that can be pulled via IMD.
221         int        nat = 0;
222         //! Part of the atoms that are local.
223         int        nat_loc = 0;
224         //! Global indices of the IMD atoms.
225         int       *ind = nullptr;
226         //! Local indices of the IMD atoms.
227         int       *ind_loc = nullptr;
228         //! Allocation size for ind_loc.
229         int        nalloc_loc = 0;
230         //! Positions for all IMD atoms assembled on the master node.
231         rvec      *xa = nullptr;
232         //! Shifts for all IMD atoms, to make molecule(s) whole.
233         ivec      *xa_shifts = nullptr;
234         //! Extra shifts since last DD step.
235         ivec      *xa_eshifts = nullptr;
236         //! Old positions for all IMD atoms on master.
237         rvec      *xa_old = nullptr;
238         //! Position of each local atom in the collective array.
239         int       *xa_ind = nullptr;
240
241         //! Global IMD frequency, known to all ranks.
242         int             nstimd = 1;
243         //! New frequency from IMD client, master only.
244         int             nstimd_new = 1;
245         //! Default IMD frequency when disconnected.
246         int             defaultNstImd = -1;
247
248         //! Port to use for network socket.
249         int             port = 0;
250         //! The IMD socket on the master node.
251         IMDSocket      *socket = nullptr;
252         //! The IMD socket on the client.
253         IMDSocket      *clientsocket = nullptr;
254         //! Length we got with last header.
255         int             length = 0;
256
257         //! Shall we block and wait for connection?
258         bool        bWConnect = false;
259         //! Set if MD is terminated.
260         bool        bTerminated = false;
261         //! Set if MD can be terminated.
262         bool        bTerminatable = false;
263         //! Set if connection is present.
264         bool        bConnected = false;
265         //! Set if we received new forces.
266         bool        bNewForces = false;
267         //! Set if pulling from VMD is allowed.
268         bool        bForceActivated = false;
269
270         //! Pointer to energies we send back.
271         IMDEnergyBlock *energies = nullptr;
272
273         //! Number of VMD forces.
274         int32_t         vmd_nforces = 0;
275         //! VMD forces indices.
276         int32_t        *vmd_f_ind = nullptr;
277         //! The VMD forces flat in memory.
278         float          *vmd_forces = nullptr;
279         //! Number of actual MD forces; this gets communicated to the clients.
280         int             nforces = 0;
281         //! Force indices.
282         int            *f_ind = nullptr;
283         //! The IMD pulling forces.
284         rvec           *f = nullptr;
285
286         //! Buffer for force sending.
287         char           *forcesendbuf = nullptr;
288         //! Buffer for coordinate sending.
289         char           *coordsendbuf = nullptr;
290         //! Send buffer for energies.
291         char           *energysendbuf = nullptr;
292         //! Buffer to make molecules whole before sending.
293         rvec           *sendxbuf = nullptr;
294
295         //! Molecules block in IMD group.
296         t_block         mols;
297
298         /* The next block is used on the master node only to reduce the output
299          * without sacrificing information. If any of these values changes,
300          * we need to write output */
301         //! Old value for nforces.
302         int       old_nforces = 0;
303         //! Old values for force indices.
304         int      *old_f_ind = nullptr;
305         //! Old values for IMD pulling forces.
306         rvec     *old_forces = nullptr;
307
308         //! Logger
309         const MDLogger  &mdlog;
310         //! Commmunication object
311         const t_commrec *cr = nullptr;
312         //! Wallcycle counting manager.
313         gmx_wallcycle   *wcycle = nullptr;
314         //! Energy output handler
315         gmx_enerdata_t  *enerd = nullptr;
316 };
317
318 /*! \brief Enum for types of IMD messages.
319  *
320  * We use the same records as the NAMD/VMD IMD implementation.
321  */
322 typedef enum IMDType_t
323 {
324     IMD_DISCONNECT,      /**< client disconnect                               */
325     IMD_ENERGIES,        /**< energy data                                     */
326     IMD_FCOORDS,         /**< atomic coordinates                              */
327     IMD_GO,              /**< start command for the simulation                */
328     IMD_HANDSHAKE,       /**< handshake to determine little/big endianness    */
329     IMD_KILL,            /**< terminates the simulation                       */
330     IMD_MDCOMM,          /**< force data                                      */
331     IMD_PAUSE,           /**< pauses the simulation                           */
332     IMD_TRATE,           /**< sets the IMD transmission and processing rate   */
333     IMD_IOERROR,         /**< I/O error                                       */
334     IMD_NR               /**< number of entries                               */
335 } IMDMessageType;
336
337
338 /*! \brief Names of the IMDType for error messages.
339  */
340 static const char *eIMDType_names[IMD_NR + 1] = {
341     "IMD_DISCONNECT",
342     "IMD_ENERGIES",
343     "IMD_FCOORDS",
344     "IMD_GO",
345     "IMD_HANDSHAKE",
346     "IMD_KILL",
347     "IMD_MDCOMM",
348     "IMD_PAUSE",
349     "IMD_TRATE",
350     "IMD_IOERROR",
351     nullptr
352 };
353
354
355
356 /*! \brief Fills the header with message and the length argument. */
357 static void fill_header(IMDHeader *header, IMDMessageType type, int32_t length)
358 {
359     /* We (ab-)use htonl network function for the correct endianness */
360     header->type   = imd_htonl(static_cast<int32_t>(type));
361     header->length = imd_htonl(length);
362 }
363
364
365 /*! \brief Swaps the endianess of the header. */
366 static void swap_header(IMDHeader *header)
367 {
368     /* and vice versa... */
369     header->type   = imd_ntohl(header->type);
370     header->length = imd_ntohl(header->length);
371 }
372
373
374 /*! \brief Reads multiple bytes from socket. */
375 static int32_t imd_read_multiple(IMDSocket *socket, char *datptr, int32_t toread)
376 {
377     int32_t leftcount, countread;
378
379
380     leftcount = toread;
381     /* Try to read while we haven't reached toread */
382     while (leftcount != 0)
383     {
384         if ((countread = imdsock_read(socket, datptr, leftcount)) < 0)
385         {
386             /* interrupted function call, try again... */
387             if (errno == EINTR)
388             {
389                 countread = 0;
390             }
391             /* this is an unexpected error, return what we got */
392             else
393             {
394                 return toread - leftcount;
395             }
396
397             /* nothing read, finished */
398         }
399         else if (countread == 0)
400         {
401             break;
402         }
403         leftcount -= countread;
404         datptr    += countread;
405     }   /* end while */
406
407     /* return nr of bytes read */
408     return toread - leftcount;
409 }
410
411
412 /*! \brief Writes multiple bytes to socket in analogy to imd_read_multiple. */
413 static int32_t imd_write_multiple(IMDSocket *socket, const char *datptr, int32_t towrite)
414 {
415     int32_t leftcount, countwritten;
416
417
418     leftcount = towrite;
419     while (leftcount != 0)
420     {
421         if ((countwritten = imdsock_write(socket, datptr, leftcount)) <= 0)
422         {
423             if (errno == EINTR)
424             {
425                 countwritten = 0;
426             }
427             else
428             {
429                 return towrite - leftcount;
430             }
431         }
432         leftcount -= countwritten;
433         datptr    += countwritten;
434     }   /* end while */
435
436     return towrite - leftcount;
437 }
438
439
440 /*! \brief Handshake with IMD client. */
441 static int imd_handshake(IMDSocket *socket)
442 {
443     IMDHeader header;
444
445
446     fill_header(&header, IMD_HANDSHAKE, 1);
447     header.length = c_protocolVersion; /* client wants unswapped version */
448
449     return static_cast<int>(imd_write_multiple(socket, reinterpret_cast<char *>(&header), c_headerSize) != c_headerSize);
450 }
451
452
453 /*! \brief Send energies using the energy block and the send buffer. */
454 static int imd_send_energies(IMDSocket *socket, const IMDEnergyBlock *energies, char *buffer)
455 {
456     int32_t recsize;
457
458
459     recsize = c_headerSize + sizeof(IMDEnergyBlock);
460     fill_header(reinterpret_cast<IMDHeader *>(buffer), IMD_ENERGIES, 1);
461     memcpy(buffer + c_headerSize, energies, sizeof(IMDEnergyBlock));
462
463     return static_cast<int>(imd_write_multiple(socket, buffer, recsize) != recsize);
464 }
465
466
467 /*! \brief Receive IMD header from socket, sets the length and returns the IMD message. */
468 static IMDMessageType imd_recv_header(IMDSocket *socket, int32_t *length)
469 {
470     IMDHeader header;
471
472
473     if (imd_read_multiple(socket, reinterpret_cast<char *>(&header), c_headerSize) != c_headerSize)
474     {
475         return IMD_IOERROR;
476     }
477     swap_header(&header);
478     *length = header.length;
479
480     return static_cast<IMDMessageType>(header.type);
481 }
482
483
484 /*! \brief Receive force indices and forces.
485  *
486  * The number of forces was previously communicated via the header.
487  */
488 static bool imd_recv_mdcomm(IMDSocket *socket, int32_t nforces, int32_t *forcendx, float *forces)
489 {
490     /* reading indices */
491     int retsize  = sizeof(int32_t) * nforces;
492     int retbytes = imd_read_multiple(socket, reinterpret_cast<char *>(forcendx), retsize);
493     if (retbytes != retsize)
494     {
495         return false;
496     }
497
498     /* reading forces as float array */
499     retsize  = 3 * sizeof(float) * nforces;
500     retbytes = imd_read_multiple(socket, reinterpret_cast<char *>(forces), retsize);
501     return (retbytes == retsize);
502 }
503
504 /* GROMACS specific functions for the IMD implementation */
505 void write_IMDgroup_to_file(bool bIMD, t_inputrec *ir, const t_state *state,
506                             const gmx_mtop_t *sys, int nfile, const t_filenm fnm[])
507 {
508     t_atoms IMDatoms;
509
510
511     if (bIMD)
512     {
513         IMDatoms = gmx_mtop_global_atoms(sys);
514         write_sto_conf_indexed(opt2fn("-imd", nfile, fnm), "IMDgroup", &IMDatoms,
515                                state->x.rvec_array(), state->v.rvec_array(), ir->ePBC, state->box, ir->imd->nat, ir->imd->ind);
516     }
517 }
518
519
520 void
521 ImdSession::dd_make_local_IMD_atoms(const gmx_domdec_t *dd)
522 {
523     if (!impl_->sessionPossible)
524     {
525         return;
526     }
527
528     dd_make_local_group_indices(dd->ga2la, impl_->nat, impl_->ind, &impl_->nat_loc,
529                                 &impl_->ind_loc, &impl_->nalloc_loc, impl_->xa_ind);
530 }
531
532
533 /*! \brief Send positions from rvec.
534  *
535  * We need a separate send buffer and conversion to Angstrom.
536  */
537 static int imd_send_rvecs(IMDSocket *socket, int nat, rvec *x, char *buffer)
538 {
539     int32_t     size;
540     int         i;
541     float       sendx[3];
542     int         tuplesize = 3 * sizeof(float);
543
544
545     /* Required size for the send buffer */
546     size = c_headerSize + 3 * sizeof(float) * nat;
547
548     /* Prepare header */
549     fill_header(reinterpret_cast<IMDHeader *>(buffer), IMD_FCOORDS, static_cast<int32_t>(nat));
550     for (i = 0; i < nat; i++)
551     {
552         sendx[0] = static_cast<float>(x[i][0]) * NM2A;
553         sendx[1] = static_cast<float>(x[i][1]) * NM2A;
554         sendx[2] = static_cast<float>(x[i][2]) * NM2A;
555         memcpy(buffer + c_headerSize + i * tuplesize, sendx, tuplesize);
556     }
557
558     return static_cast<int>(imd_write_multiple(socket, buffer, size) != size);
559 }
560
561
562 void
563 ImdSession::Impl::prepareMasterSocket()
564 {
565     if (imdsock_winsockinit() == -1)
566     {
567         gmx_fatal(FARGS, "%s Failed to initialize winsock.\n", IMDstr);
568     }
569
570     /* The rest is identical, first create and bind a socket and set to listen then. */
571     GMX_LOG(mdlog.warning).appendTextFormatted("%s Setting up incoming socket.", IMDstr);
572     socket = imdsock_create();
573     if (!socket)
574     {
575         gmx_fatal(FARGS, "%s Failed to create socket.", IMDstr);
576     }
577
578     /* Bind to port */
579     int ret = imdsock_bind(socket, port);
580     if (ret)
581     {
582         gmx_fatal(FARGS, "%s binding socket to port %d failed with error %d.\n", IMDstr, port, ret);
583     }
584
585     if (imd_sock_listen(socket))
586     {
587         gmx_fatal(FARGS, "%s socket listen failed with error %d.\n", IMDstr, ret);
588     }
589
590     if (imdsock_getport(socket, &port))
591     {
592         gmx_fatal(FARGS, "%s Could not determine port number.\n", IMDstr);
593     }
594
595     GMX_LOG(mdlog.warning).appendTextFormatted("%s Listening for IMD connection on port %d.", IMDstr, port);
596 }
597
598
599 void
600 ImdSession::Impl::disconnectClient()
601 {
602     /* Write out any buffered pulling data */
603     fflush(outf);
604
605     /* we first try to shut down the clientsocket */
606     imdsock_shutdown(clientsocket);
607     if (!imdsock_destroy(clientsocket))
608     {
609         GMX_LOG(mdlog.warning).appendTextFormatted("%s Failed to destroy socket.", IMDstr);
610     }
611
612     /* then we reset the IMD step to its default, and reset the connection boolean */
613     nstimd_new   = defaultNstImd;
614     clientsocket = nullptr;
615     bConnected   = false;
616 }
617
618
619 void
620 ImdSession::Impl::issueFatalError(const char *msg)
621 {
622     GMX_LOG(mdlog.warning).appendTextFormatted("%s %s", IMDstr, msg);
623     disconnectClient();
624     GMX_LOG(mdlog.warning).appendTextFormatted("%s disconnected.", IMDstr);
625 }
626
627
628 bool
629 ImdSession::Impl::tryConnect()
630 {
631     if (imdsock_tryread(socket, 0, 0) > 0)
632     {
633         /* yes, we got something, accept on clientsocket */
634         clientsocket = imdsock_accept(socket);
635         if (!clientsocket)
636         {
637             GMX_LOG(mdlog.warning).appendTextFormatted("%s Accepting the connection on the socket failed.", IMDstr);
638             return false;
639         }
640
641         /* handshake with client */
642         if (imd_handshake(clientsocket))
643         {
644             issueFatalError("Connection failed.");
645             return false;
646         }
647
648         GMX_LOG(mdlog.warning).appendTextFormatted("%s Connection established, checking if I got IMD_GO orders.", IMDstr);
649
650         /* Check if we get the proper "GO" command from client. */
651         if (imdsock_tryread(clientsocket, c_connectWait, 0) != 1 || imd_recv_header(clientsocket, &(length)) != IMD_GO)
652         {
653             issueFatalError("No IMD_GO order received. IMD connection failed.");
654         }
655
656         /* IMD connected */
657         bConnected = true;
658
659         return true;
660     }
661
662     return false;
663 }
664
665
666 void
667 ImdSession::Impl::blockConnect()
668 {
669     /* do not wait for connection, when e.g. ctrl+c is pressed and we will terminate anyways. */
670     if (!(static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none))
671     {
672         return;
673     }
674
675     GMX_LOG(mdlog.warning).appendTextFormatted("%s Will wait until I have a connection and IMD_GO orders.", IMDstr);
676
677     /* while we have no clientsocket... 2nd part: we should still react on ctrl+c */
678     while ((!clientsocket) && (static_cast<int>(gmx_get_stop_condition()) == gmx_stop_cond_none))
679     {
680         tryConnect();
681         imd_sleep(c_loopWait);
682     }
683 }
684
685
686 void
687 ImdSession::Impl::prepareVmdForces()
688 {
689     srenew((vmd_f_ind), vmd_nforces);
690     srenew((vmd_forces), 3*vmd_nforces);
691 }
692
693
694 void
695 ImdSession::Impl::readVmdForces()
696 {
697     /* the length of the previously received header tells us the nr of forces we will receive */
698     vmd_nforces = length;
699     /* prepare the arrays */
700     prepareVmdForces();
701     /* Now we read the forces... */
702     if (!(imd_recv_mdcomm(clientsocket, vmd_nforces, vmd_f_ind, vmd_forces)))
703     {
704         issueFatalError("Error while reading forces from remote. Disconnecting");
705     }
706 }
707
708
709 void
710 ImdSession::Impl::prepareMDForces()
711 {
712     srenew((f_ind), nforces);
713     srenew((f    ), nforces);
714 }
715
716
717 void
718 ImdSession::Impl::copyToMDForces()
719 {
720     int  i;
721     real conversion = CAL2JOULE * NM2A;
722
723
724     for (i = 0; i < nforces; i++)
725     {
726         /* Copy the indices, a copy is important because we may update the incoming forces
727          * whenever we receive new forces while the MD forces are only communicated upon
728          * IMD communication.*/
729         f_ind[i] = vmd_f_ind[i];
730
731         /* Convert to rvecs and do a proper unit conversion */
732         f[i][0] = vmd_forces[3*i    ] * conversion;
733         f[i][1] = vmd_forces[3*i + 1] * conversion;
734         f[i][2] = vmd_forces[3*i + 2] * conversion;
735     }
736 }
737
738
739 bool
740 ImdSession::Impl::bForcesChanged() const
741 {
742     /* First, check whether the number of pulled atoms changed */
743     if (nforces != old_nforces)
744     {
745         return true;
746     }
747
748     /* Second, check whether any of the involved atoms changed */
749     for (int i = 0; i < nforces; i++)
750     {
751         if (f_ind[i] != old_f_ind[i])
752         {
753             return true;
754         }
755     }
756
757     // TODO replace this with rvecs_differ
758     /* Third, check whether all forces are the same */
759     for (int i = 0; i < nforces; i++)
760     {
761         if (f[i][XX] != old_forces[i][XX])
762         {
763             return true;
764         }
765         if (f[i][YY] != old_forces[i][YY])
766         {
767             return true;
768         }
769         if (f[i][ZZ] != old_forces[i][ZZ])
770         {
771             return true;
772         }
773     }
774
775     /* All old and new forces are identical! */
776     return false;
777 }
778
779
780 void
781 ImdSession::Impl::keepOldValues()
782 {
783     old_nforces = nforces;
784
785     for (int i = 0; i < nforces; i++)
786     {
787         old_f_ind[i] = f_ind[i];
788         copy_rvec(f[i], old_forces[i]);
789     }
790 }
791
792
793 /*! \brief Returns true if any component of the two rvecs differs. */
794 static inline bool rvecs_differ(const rvec v1, const rvec v2)
795 {
796     for (int i = 0; i < DIM; i++)
797     {
798         if (v1[i] != v2[i])
799         {
800             return true;
801         }
802     }
803
804     return false;
805 }
806
807
808 void
809 ImdSession::Impl::outputForces(double time)
810 {
811     if (!bForcesChanged())
812     {
813         return;
814     }
815
816     /* Write time and total number of applied IMD forces */
817     fprintf(outf, "%14.6e%6d", time, nforces);
818
819     /* Write out the global atom indices of the pulled atoms and the forces itself,
820      * write out a force only if it has changed since the last output */
821     for (int i = 0; i < nforces; i++)
822     {
823         if (rvecs_differ(f[i], old_forces[i]))
824         {
825             fprintf(outf, "%9d", ind[f_ind[i]] + 1);
826             fprintf(outf, "%12.4e%12.4e%12.4e", f[i][0], f[i][1], f[i][2]);
827         }
828     }
829     fprintf(outf, "\n");
830
831     keepOldValues();
832 }
833
834
835 void
836 ImdSession::Impl::syncNodes(const t_commrec *cr, double t)
837 {
838     /* Notify the other nodes whether we are still connected. */
839     if (PAR(cr))
840     {
841         block_bc(cr, bConnected);
842     }
843
844     /* ...if not connected, the job is done here. */
845     if (!bConnected)
846     {
847         return;
848     }
849
850     /* Let the other nodes know whether we got a new IMD synchronization frequency. */
851     if (PAR(cr))
852     {
853         block_bc(cr, nstimd_new);
854     }
855
856     /* Now we all set the (new) nstimd communication time step */
857     nstimd = nstimd_new;
858
859     /* We're done if we don't allow pulling at all */
860     if (!(bForceActivated))
861     {
862         return;
863     }
864
865     /* OK, let's check if we have received forces which we need to communicate
866      * to the other nodes */
867     int new_nforces = 0;
868     if (MASTER(cr))
869     {
870         if (bNewForces)
871         {
872             new_nforces = vmd_nforces;
873         }
874         /* make the "new_forces" negative, when we did not receive new ones */
875         else
876         {
877             new_nforces = vmd_nforces * -1;
878         }
879     }
880
881     /* make new_forces known to the clients */
882     if (PAR(cr))
883     {
884         block_bc(cr, new_nforces);
885     }
886
887     /* When new_natoms < 0 then we know that these are still the same forces
888      * so we don't communicate them, otherwise... */
889     if (new_nforces < 0)
890     {
891         return;
892     }
893
894     /* set local VMD and nforces */
895     vmd_nforces = new_nforces;
896     nforces     = vmd_nforces;
897
898     /* now everybody knows the number of forces in f_ind, so we can prepare
899      * the target arrays for indices and forces */
900     prepareMDForces();
901
902     /* we first update the MD forces on the master by converting the VMD forces */
903     if (MASTER(cr))
904     {
905         copyToMDForces();
906         /* We also write out forces on every update, so that we know which
907          * forces are applied for every step */
908         if (outf)
909         {
910             outputForces(t);
911         }
912     }
913
914     /* In parallel mode we communicate the to-be-applied forces to the other nodes */
915     if (PAR(cr))
916     {
917         nblock_bc(cr, nforces, f_ind);
918         nblock_bc(cr, nforces, f    );
919     }
920
921     /* done communicating the forces, reset bNewForces */
922     bNewForces = false;
923 }
924
925
926 void
927 ImdSession::Impl::readCommand()
928 {
929     bool       IMDpaused = false;
930
931     while (clientsocket && (imdsock_tryread(clientsocket, 0, 0) > 0 || IMDpaused))
932     {
933         IMDMessageType itype = imd_recv_header(clientsocket, &(length));
934         /* let's see what we got: */
935         switch (itype)
936         {
937             /* IMD asks us to terminate the simulation, check if the user allowed this */
938             case IMD_KILL:
939                 if (bTerminatable)
940                 {
941                     GMX_LOG(mdlog.warning).appendTextFormatted(" %s Terminating connection and running simulation (if supported by integrator).", IMDstr);
942                     bTerminated = true;
943                     bWConnect   = false;
944                     gmx_set_stop_condition(gmx_stop_cond_next);
945                 }
946                 else
947                 {
948                     GMX_LOG(mdlog.warning).appendTextFormatted(" %s Set -imdterm command line switch to allow mdrun termination from within IMD.", IMDstr);
949                 }
950
951                 break;
952
953             /* the client doen't want to talk to us anymore */
954             case IMD_DISCONNECT:
955                 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Disconnecting client.", IMDstr);
956                 disconnectClient();
957                 break;
958
959             /* we got new forces, read them and set bNewForces flag */
960             case IMD_MDCOMM:
961                 readVmdForces();
962                 bNewForces = true;
963                 break;
964
965             /* the client asks us to (un)pause the simulation. So we toggle the IMDpaused state */
966             case IMD_PAUSE:
967                 if (IMDpaused)
968                 {
969                     GMX_LOG(mdlog.warning).appendTextFormatted(" %s Un-pause command received.", IMDstr);
970                     IMDpaused = false;
971                 }
972                 else
973                 {
974                     GMX_LOG(mdlog.warning).appendTextFormatted(" %s Pause command received.", IMDstr);
975                     IMDpaused = true;
976                 }
977
978                 break;
979
980             /* the client sets a new transfer rate, if we get 0, we reset the rate
981              * to the default. VMD filters 0 however */
982             case IMD_TRATE:
983                 nstimd_new = (length > 0) ? length : defaultNstImd;
984                 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Update frequency will be set to %d.", IMDstr, nstimd_new);
985                 break;
986
987             /* Catch all rule for the remaining IMD types which we don't expect */
988             default:
989                 GMX_LOG(mdlog.warning).appendTextFormatted(" %s Received unexpected %s.", IMDstr, enum_name(static_cast<int>(itype), IMD_NR, eIMDType_names));
990                 issueFatalError("Terminating connection");
991                 break;
992         } /* end switch */
993     }     /* end while  */
994 }
995
996
997 void
998 ImdSession::Impl::openOutputFile(const char                *fn,
999                                  int                        nat_total,
1000                                  const gmx_output_env_t    *oenv,
1001                                  const ContinuationOptions &continuationOptions)
1002 {
1003     /* Open log file of applied IMD forces if requested */
1004     if (!fn || !oenv)
1005     {
1006         fprintf(stdout, "%s For a log of the IMD pull forces explicitly specify '-if' on the command line.\n"
1007                 "%s (Not possible with energy minimization.)\n", IMDstr, IMDstr);
1008         return;
1009     }
1010
1011     /* If we append to an existing file, all the header information is already there */
1012     if (continuationOptions.appendFiles)
1013     {
1014         outf = gmx_fio_fopen(fn, "a+");
1015     }
1016     else
1017     {
1018         outf = gmx_fio_fopen(fn, "w+");
1019         if (nat == nat_total)
1020         {
1021             fprintf(outf, "# Note that you can select an IMD index group in the .mdp file if a subset of the atoms suffices.\n");
1022         }
1023
1024         xvgr_header(outf, "IMD Pull Forces", "Time (ps)", "# of Forces / Atom IDs / Forces (kJ/mol)", exvggtNONE, oenv);
1025
1026         fprintf(outf, "# Can display and manipulate %d (of a total of %d) atoms via IMD.\n", nat, nat_total);
1027         fprintf(outf, "# column 1    : time (ps)\n");
1028         fprintf(outf, "# column 2    : total number of atoms feeling an IMD pulling force at that time\n");
1029         fprintf(outf, "# cols. 3.-6  : global atom number of pulled atom, x-force, y-force, z-force (kJ/mol)\n");
1030         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");
1031         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");
1032         fflush(outf);
1033     }
1034
1035     /* To reduce the output file size we remember the old values and output only
1036      * when something changed */
1037     snew(old_f_ind, nat);  /* One can never pull on more atoms */
1038     snew(old_forces, nat);
1039 }
1040
1041
1042 ImdSession::Impl::Impl(const MDLogger &mdlog)
1043     : mdlog(mdlog)
1044 {
1045     init_block(&mols);
1046 }
1047
1048 ImdSession::Impl::~Impl()
1049 {
1050     if (outf)
1051     {
1052         gmx_fio_fclose(outf);
1053     }
1054     done_block(&mols);
1055 }
1056
1057
1058 void
1059 ImdSession::Impl::prepareMoleculesInImdGroup(const gmx_mtop_t *top_global)
1060 {
1061     /* check whether index is sorted */
1062     for (int i = 0; i < nat-1; i++)
1063     {
1064         if (ind[i] > ind[i+1])
1065         {
1066             gmx_fatal(FARGS, "%s IMD index is not sorted. This is currently not supported.\n", IMDstr);
1067         }
1068     }
1069
1070     RangePartitioning gmols = gmx_mtop_molecules(*top_global);
1071     t_block           lmols;
1072     lmols.nr = 0;
1073     snew(lmols.index, gmols.numBlocks() + 1);
1074     lmols.index[0] = 0;
1075
1076     for (int i = 0; i < gmols.numBlocks(); i++)
1077     {
1078         auto mol   = gmols.block(i);
1079         int  count = 0;
1080         for (int ii = 0; ii < nat; ii++)
1081         {
1082             if (mol.inRange(ind[ii]))
1083             {
1084                 count += 1;
1085             }
1086         }
1087         if (count > 0)
1088         {
1089             lmols.index[lmols.nr+1] = lmols.index[lmols.nr]+count;
1090             lmols.nr               += 1;
1091         }
1092     }
1093     srenew(lmols.index, lmols.nr+1);
1094     lmols.nalloc_index = lmols.nr+1;
1095     mols               = lmols;
1096 }
1097
1098
1099 /*! \brief Copied and modified from groupcoord.c shift_positions_group(). */
1100 static void shift_positions(
1101         const matrix box,
1102         rvec         x[], /* The positions [0..nr] */
1103         const ivec   is,  /* The shift [0..nr] */
1104         int          nr)  /* The number of positions */
1105 {
1106     int      i, tx, ty, tz;
1107
1108     /* Loop over the group's atoms */
1109     if (TRICLINIC(box))
1110     {
1111         for (i = 0; i < nr; i++)
1112         {
1113             tx = is[XX];
1114             ty = is[YY];
1115             tz = is[ZZ];
1116
1117             x[i][XX] = x[i][XX]-tx*box[XX][XX]-ty*box[YY][XX]-tz*box[ZZ][XX];
1118             x[i][YY] = x[i][YY]-ty*box[YY][YY]-tz*box[ZZ][YY];
1119             x[i][ZZ] = x[i][ZZ]-tz*box[ZZ][ZZ];
1120         }
1121     }
1122     else
1123     {
1124         for (i = 0; i < nr; i++)
1125         {
1126             tx = is[XX];
1127             ty = is[YY];
1128             tz = is[ZZ];
1129
1130             x[i][XX] = x[i][XX]-tx*box[XX][XX];
1131             x[i][YY] = x[i][YY]-ty*box[YY][YY];
1132             x[i][ZZ] = x[i][ZZ]-tz*box[ZZ][ZZ];
1133         }
1134     }
1135 }
1136
1137
1138 void ImdSession::Impl::removeMolecularShifts(const matrix box)
1139 {
1140     /* for each molecule also present in IMD group */
1141     for (int i = 0; i < mols.nr; i++)
1142     {
1143         /* first we determine the minimum and maximum shifts for each molecule */
1144
1145         ivec largest, smallest, shift;
1146         clear_ivec(largest);
1147         clear_ivec(smallest);
1148         clear_ivec(shift);
1149
1150         copy_ivec(xa_shifts[mols.index[i]], largest);
1151         copy_ivec(xa_shifts[mols.index[i]], smallest);
1152
1153         for (int ii = mols.index[i]+1; ii < mols.index[i+1]; ii++)
1154         {
1155             if (xa_shifts[ii][XX] > largest[XX])
1156             {
1157                 largest[XX]  = xa_shifts[ii][XX];
1158             }
1159             if (xa_shifts[ii][XX] < smallest[XX])
1160             {
1161                 smallest[XX] = xa_shifts[ii][XX];
1162             }
1163
1164             if (xa_shifts[ii][YY] > largest[YY])
1165             {
1166                 largest[YY]  = xa_shifts[ii][YY];
1167             }
1168             if (xa_shifts[ii][YY] < smallest[YY])
1169             {
1170                 smallest[YY] = xa_shifts[ii][YY];
1171             }
1172
1173             if (xa_shifts[ii][ZZ] > largest[ZZ])
1174             {
1175                 largest[ZZ]  = xa_shifts[ii][ZZ];
1176             }
1177             if (xa_shifts[ii][ZZ] < smallest[ZZ])
1178             {
1179                 smallest[ZZ] = xa_shifts[ii][ZZ];
1180             }
1181
1182         }
1183
1184         /* check if we what we can subtract/add to the positions
1185          * to put them back into the central box */
1186         if (smallest[XX] > 0)
1187         {
1188             shift[XX] = smallest[XX];
1189         }
1190         if (smallest[YY] > 0)
1191         {
1192             shift[YY] = smallest[YY];
1193         }
1194         if (smallest[ZZ] > 0)
1195         {
1196             shift[ZZ] = smallest[ZZ];
1197         }
1198
1199         if (largest[XX] < 0)
1200         {
1201             shift[XX] = largest[XX];
1202         }
1203         if (largest[YY] < 0)
1204         {
1205             shift[YY] = largest[YY];
1206         }
1207         if (largest[ZZ] < 0)
1208         {
1209             shift[ZZ] = largest[ZZ];
1210         }
1211
1212         /* is there a shift at all? */
1213         if ((shift[XX]) || (shift[YY]) || (shift[ZZ]))
1214         {
1215             int molsize = mols.index[i+1]-mols.index[i];
1216             /* shift the positions */
1217             shift_positions(box, &(xa[mols.index[i]]), shift, molsize);
1218         }
1219
1220     }
1221 }
1222
1223
1224 void
1225 ImdSession::Impl::prepareForPositionAssembly(const t_commrec *cr, const rvec x[])
1226 {
1227     snew(xa,         nat);
1228     snew(xa_ind,     nat);
1229     snew(xa_shifts,  nat);
1230     snew(xa_eshifts, nat);
1231     snew(xa_old,     nat);
1232
1233     /* Save the original (whole) set of positions such that later the
1234      * molecule can always be made whole again */
1235     if (MASTER(cr))
1236     {
1237         for (int i = 0; i < nat; i++)
1238         {
1239             int ii = ind[i];
1240             copy_rvec(x[ii], xa_old[i]);
1241         }
1242     }
1243
1244     if (!PAR(cr))
1245     {
1246         nat_loc = nat;
1247         ind_loc = ind;
1248
1249         /* xa_ind[i] needs to be set to i for serial runs */
1250         for (int i = 0; i < nat; i++)
1251         {
1252             xa_ind[i] = i;
1253         }
1254     }
1255
1256     /* Communicate initial coordinates xa_old to all processes */
1257     if (PAR(cr))
1258     {
1259         gmx_bcast(nat * sizeof(xa_old[0]), xa_old, cr);
1260     }
1261 }
1262
1263
1264 /*! \brief Check for non-working integrator / parallel options. */
1265 static void imd_check_integrator_parallel(const t_inputrec *ir, const t_commrec *cr)
1266 {
1267     if (PAR(cr))
1268     {
1269         if (((ir->eI) == eiSteep) || ((ir->eI) == eiCG) || ((ir->eI) == eiLBFGS) || ((ir->eI) == eiNM))
1270         {
1271             gmx_fatal(FARGS, "%s Energy minimization via steep, CG, lbfgs and nm in parallel is currently not supported by IMD.\n", IMDstr);
1272         }
1273     }
1274 }
1275
1276 std::unique_ptr<ImdSession>
1277 makeImdSession(const t_inputrec        *ir,
1278                const t_commrec         *cr,
1279                gmx_wallcycle           *wcycle,
1280                gmx_enerdata_t          *enerd,
1281                const gmx_multisim_t    *ms,
1282                const gmx_mtop_t        *top_global,
1283                const MDLogger          &mdlog,
1284                const rvec               x[],
1285                int                      nfile,
1286                const t_filenm           fnm[],
1287                const gmx_output_env_t  *oenv,
1288                const MdrunOptions      &mdrunOptions)
1289 {
1290     std::unique_ptr<ImdSession> session(new ImdSession(mdlog));
1291     auto impl = session->impl_.get();
1292
1293     /* We will allow IMD sessions only if supported by the binary and
1294        explicitly enabled in the .tpr file */
1295     if (!GMX_IMD || !ir->bIMD)
1296     {
1297         return session;
1298     }
1299
1300     // TODO As IMD is intended for interactivity, and the .tpr file
1301     // opted in for that, it is acceptable to write more terminal
1302     // output than in a typical simulation. However, all the GMX_LOG
1303     // statements below should go to both the log file and to the
1304     // terminal. This is probably be implemented by adding a logging
1305     // stream named like ImdInfo, to separate it from warning and to
1306     // send it to both destinations.
1307
1308     if (EI_DYNAMICS(ir->eI))
1309     {
1310         impl->defaultNstImd = ir->nstcalcenergy;
1311     }
1312     else if (EI_ENERGY_MINIMIZATION(ir->eI))
1313     {
1314         impl->defaultNstImd = 1;
1315     }
1316     else
1317     {
1318         GMX_LOG(mdlog.warning).appendTextFormatted("%s Integrator '%s' is not supported for Interactive Molecular Dynamics, running normally instead", IMDstr, ei_names[ir->eI]);
1319         return session;
1320     }
1321     if (isMultiSim(ms))
1322     {
1323         GMX_LOG(mdlog.warning).appendTextFormatted("%s Cannot use IMD for multiple simulations or replica exchange, running normally instead", IMDstr);
1324         return session;
1325     }
1326
1327     const auto &options = mdrunOptions.imdOptions;
1328
1329     bool        createSession = false;
1330     /* It seems we have a .tpr file that defines an IMD group and thus allows IMD connections.
1331      * Check whether we can actually provide the IMD functionality for this setting: */
1332     if (MASTER(cr))
1333     {
1334         /* Check whether IMD was enabled by one of the command line switches: */
1335         if (options.wait || options.terminatable || options.pull)
1336         {
1337             GMX_LOG(mdlog.warning).appendTextFormatted("%s Enabled. This simulation will accept incoming IMD connections.", IMDstr);
1338             createSession = true;
1339         }
1340         else
1341         {
1342             GMX_LOG(mdlog.warning).appendTextFormatted("%s None of the -imd switches was used.\n"
1343                                                        "%s This run will not accept incoming IMD connections", IMDstr, IMDstr);
1344         }
1345     }   /* end master only */
1346
1347     /* Let the other nodes know whether we want IMD */
1348     if (PAR(cr))
1349     {
1350         block_bc(cr, createSession);
1351     }
1352
1353     /*... if not we are done.*/
1354     if (!createSession)
1355     {
1356         return session;
1357     }
1358
1359
1360     /* check if we're using a sane integrator / parallel combination */
1361     imd_check_integrator_parallel(ir, cr);
1362
1363
1364     /*
1365      *****************************************************
1366      * From here on we assume that IMD is turned on      *
1367      *****************************************************
1368      */
1369
1370     int nat_total = top_global->natoms;
1371
1372     /* Initialize IMD session. If we read in a pre-IMD .tpr file, ir->imd->nat
1373      * will be zero. For those cases we transfer _all_ atomic positions */
1374     impl->sessionPossible = true;
1375     impl->nat             = ir->imd->nat > 0 ? ir->imd->nat : nat_total;
1376     if (options.port >= 1)
1377     {
1378         impl->port        = options.port;
1379     }
1380     impl->cr     = cr;
1381     impl->wcycle = wcycle;
1382     impl->enerd  = enerd;
1383
1384     /* We might need to open an output file for IMD forces data */
1385     if (MASTER(cr))
1386     {
1387         impl->openOutputFile(opt2fn("-if", nfile, fnm), nat_total, oenv, mdrunOptions.continuationOptions);
1388     }
1389
1390     /* Make sure that we operate with a valid atom index array for the IMD atoms */
1391     if (ir->imd->nat > 0)
1392     {
1393         /* Point to the user-supplied array of atom numbers */
1394         impl->ind = ir->imd->ind;
1395     }
1396     else
1397     {
1398         /* Make a dummy (ind[i] = i) array of all atoms */
1399         snew(impl->ind, nat_total);
1400         for (int i = 0; i < nat_total; i++)
1401         {
1402             impl->ind[i] = i;
1403         }
1404     }
1405
1406     /* read environment on master and prepare socket for incoming connections */
1407     if (MASTER(cr))
1408     {
1409         /* we allocate memory for our IMD energy structure */
1410         int32_t recsize = c_headerSize + sizeof(IMDEnergyBlock);
1411         snew(impl->energysendbuf, recsize);
1412
1413         /* Shall we wait for a connection? */
1414         if (options.wait)
1415         {
1416             impl->bWConnect = true;
1417             GMX_LOG(mdlog.warning).appendTextFormatted("%s Pausing simulation while no IMD connection present (-imdwait).", IMDstr);
1418         }
1419
1420         /* Will the IMD clients be able to terminate the simulation? */
1421         if (options.terminatable)
1422         {
1423             impl->bTerminatable = true;
1424             GMX_LOG(mdlog.warning).appendTextFormatted("%s Allow termination of the simulation from IMD client (-imdterm).", IMDstr);
1425         }
1426
1427         /* Is pulling from IMD client allowed? */
1428         if (options.pull)
1429         {
1430             impl->bForceActivated = true;
1431             GMX_LOG(mdlog.warning).appendTextFormatted("%s Pulling from IMD remote is enabled (-imdpull).", IMDstr);
1432         }
1433
1434         /* Initialize send buffers with constant size */
1435         snew(impl->sendxbuf, impl->nat);
1436         snew(impl->energies, 1);
1437         int32_t bufxsize = c_headerSize + 3 * sizeof(float) * impl->nat;
1438         snew(impl->coordsendbuf, bufxsize);
1439     }
1440
1441     /* do we allow interactive pulling? If so let the other nodes know. */
1442     if (PAR(cr))
1443     {
1444         block_bc(cr, impl->bForceActivated);
1445     }
1446
1447     /* setup the listening socket on master process */
1448     if (MASTER(cr))
1449     {
1450         GMX_LOG(mdlog.warning).appendTextFormatted("%s Setting port for connection requests to %d.", IMDstr, impl->port);
1451         impl->prepareMasterSocket();
1452         /* Wait until we have a connection if specified before */
1453         if (impl->bWConnect)
1454         {
1455             impl->blockConnect();
1456         }
1457         else
1458         {
1459             GMX_LOG(mdlog.warning).appendTextFormatted("%s -imdwait not set, starting simulation.", IMDstr);
1460         }
1461     }
1462     /* Let the other nodes know whether we are connected */
1463     impl->syncNodes(cr, 0);
1464
1465     /* Initialize arrays used to assemble the positions from the other nodes */
1466     impl->prepareForPositionAssembly(cr, x);
1467
1468     /* Initialize molecule blocks to make them whole later...*/
1469     if (MASTER(cr))
1470     {
1471         impl->prepareMoleculesInImdGroup(top_global);
1472     }
1473
1474     return session;
1475 }
1476
1477
1478 bool ImdSession::Impl::run(int64_t      step,
1479                            bool         bNS,
1480                            const matrix box,
1481                            const rvec   x[],
1482                            double       t)
1483 {
1484     /* IMD at all? */
1485     if (!sessionPossible)
1486     {
1487         return false;
1488     }
1489
1490     wallcycle_start(wcycle, ewcIMD);
1491
1492     /* read command from client and check if new incoming connection */
1493     if (MASTER(cr))
1494     {
1495         /* If not already connected, check for new connections */
1496         if (!clientsocket)
1497         {
1498             if (bWConnect)
1499             {
1500                 blockConnect();
1501             }
1502             else
1503             {
1504                 tryConnect();
1505             }
1506         }
1507
1508         /* Let's see if we have new IMD messages for us */
1509         if (clientsocket)
1510         {
1511             readCommand();
1512         }
1513     }
1514
1515     /* is this an IMD communication step? */
1516     bool imdstep = do_per_step(step, nstimd);
1517
1518     /* OK so this is an IMD step ... */
1519     if (imdstep)
1520     {
1521         /* First we sync all nodes to let everybody know whether we are connected to VMD */
1522         syncNodes(cr, t);
1523     }
1524
1525     /* If a client is connected, we collect the positions
1526      * and put molecules back into the box before transfer */
1527     if ((imdstep && bConnected)
1528         || bNS)            /* independent of imdstep, we communicate positions at each NS step */
1529     {
1530         /* Transfer the IMD positions to the master node. Every node contributes
1531          * its local positions x and stores them in the assembled xa array. */
1532         communicate_group_positions(cr, xa, xa_shifts, xa_eshifts,
1533                                     true, x, nat, nat_loc,
1534                                     ind_loc, xa_ind, xa_old, box);
1535
1536         /* If connected and master -> remove shifts */
1537         if ((imdstep && bConnected) && MASTER(cr))
1538         {
1539             removeMolecularShifts(box);
1540         }
1541     }
1542
1543     wallcycle_stop(wcycle, ewcIMD);
1544
1545     return imdstep;
1546 }
1547
1548 bool ImdSession::run(int64_t      step,
1549                      bool         bNS,
1550                      const matrix box,
1551                      const rvec   x[],
1552                      double       t)
1553 {
1554     return impl_->run(step, bNS, box, x, t);
1555 }
1556
1557 void ImdSession::fillEnergyRecord(int64_t step,
1558                                   bool    bHaveNewEnergies)
1559 {
1560     if (!impl_->sessionPossible || !impl_->clientsocket)
1561     {
1562         return;
1563     }
1564
1565     IMDEnergyBlock *ene = impl_->energies;
1566
1567     ene->tstep = step;
1568
1569     /* In MPI-parallel simulations the energies are not accessible a at every time step.
1570      * We update them if we have new values, otherwise, the energy values from the
1571      * last global communication step are still on display in the viewer. */
1572     if (bHaveNewEnergies)
1573     {
1574         ene->T_abs   = static_cast<float>(impl_->enerd->term[F_TEMP   ]);
1575         ene->E_pot   = static_cast<float>(impl_->enerd->term[F_EPOT   ]);
1576         ene->E_tot   = static_cast<float>(impl_->enerd->term[F_ETOT   ]);
1577         ene->E_bond  = static_cast<float>(impl_->enerd->term[F_BONDS  ]);
1578         ene->E_angle = static_cast<float>(impl_->enerd->term[F_ANGLES ]);
1579         ene->E_dihe  = static_cast<float>(impl_->enerd->term[F_PDIHS  ]);
1580         ene->E_impr  = static_cast<float>(impl_->enerd->term[F_IDIHS  ]);
1581         ene->E_vdw   = static_cast<float>(impl_->enerd->term[F_LJ     ]);
1582         ene->E_coul  = static_cast<float>(impl_->enerd->term[F_COUL_SR]);
1583     }
1584 }
1585
1586
1587 void
1588 ImdSession::sendPositionsAndEnergies()
1589 {
1590     if (!impl_->sessionPossible || !impl_->clientsocket)
1591     {
1592         return;
1593     }
1594
1595     if (imd_send_energies(impl_->clientsocket, impl_->energies, impl_->energysendbuf))
1596     {
1597         impl_->issueFatalError("Error sending updated energies. Disconnecting client.");
1598     }
1599
1600     if (imd_send_rvecs(impl_->clientsocket, impl_->nat, impl_->xa, impl_->coordsendbuf))
1601     {
1602         impl_->issueFatalError("Error sending updated positions. Disconnecting client.");
1603     }
1604 }
1605
1606
1607 void
1608 ImdSession::updateEnergyRecordAndSendPositionsAndEnergies(bool    bIMDstep,
1609                                                           int64_t step,
1610                                                           bool    bHaveNewEnergies)
1611 {
1612     if (!impl_->sessionPossible)
1613     {
1614         return;
1615     }
1616
1617     wallcycle_start(impl_->wcycle, ewcIMD);
1618
1619     /* Update time step for IMD and prepare IMD energy record if we have new energies. */
1620     fillEnergyRecord(step, bHaveNewEnergies);
1621
1622     if (bIMDstep)
1623     {
1624         /* Send positions and energies to VMD client via IMD */
1625         sendPositionsAndEnergies();
1626     }
1627
1628     wallcycle_stop(impl_->wcycle, ewcIMD);
1629 }
1630
1631 void ImdSession::applyForces(rvec *f)
1632 {
1633     if (!impl_->sessionPossible || !impl_->bForceActivated)
1634     {
1635         return;
1636     }
1637
1638     wallcycle_start(impl_->wcycle, ewcIMD);
1639
1640     for (int i = 0; i < impl_->nforces; i++)
1641     {
1642         /* j are the indices in the "System group".*/
1643         int j = impl_->ind[impl_->f_ind[i]];
1644
1645         /* check if this is a local atom and find out locndx */
1646         const int       *locndx;
1647         const t_commrec *cr = impl_->cr;
1648         if (PAR(cr) && (locndx = cr->dd->ga2la->findHome(j)))
1649         {
1650             j = *locndx;
1651         }
1652
1653         rvec_inc(f[j], impl_->f[i]);
1654     }
1655
1656     wallcycle_stop(impl_->wcycle, ewcIMD);
1657 }
1658
1659 ImdSession::ImdSession(const MDLogger &mdlog)
1660     : impl_(new Impl(mdlog))
1661 {}
1662
1663 ImdSession::~ImdSession() = default;
1664
1665 } // namespace gmx