c40161d9ef31b4a88521ca77c7e75e32dd207389
[alexxy/gromacs.git] / src / gromacs / mdrun / replicaexchange.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2011-2019,2020, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 /*! \internal \file
39  *
40  * \brief Implements the replica exchange routines.
41  *
42  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43  * \author Mark Abraham <mark.j.abraham@gmail.com>
44  * \ingroup module_mdrun
45  */
46 #include "gmxpre.h"
47
48 #include "replicaexchange.h"
49
50 #include "config.h"
51
52 #include <cmath>
53
54 #include <random>
55
56 #include "gromacs/domdec/collect.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdrunutility/multisim.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/mdtypes/enerdata.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/mdtypes/state.h"
66 #include "gromacs/random/threefry.h"
67 #include "gromacs/random/uniformintdistribution.h"
68 #include "gromacs/random/uniformrealdistribution.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/pleasecite.h"
71 #include "gromacs/utility/smalloc.h"
72
73 //! Helps cut off probability values.
74 constexpr int c_probabilityCutoff = 100;
75
76 /* we don't bother evaluating if events are more rare than exp(-100) = 3.7x10^-44 */
77
78 //! Rank in the multisimulation
79 #define MSRANK(ms, nodeid) (nodeid)
80
81 //! Enum for replica exchange flavours
82 enum
83 {
84     ereTEMP,
85     ereLAMBDA,
86     ereENDSINGLE,
87     ereTL,
88     ereNR
89 };
90 /*! \brief Strings describing replica exchange flavours.
91  *
92  *  end_single_marker merely notes the end of single variable replica
93  *  exchange. All types higher than it are multiple replica exchange
94  *  methods.
95  *
96  * Eventually, should add 'pressure', 'temperature and pressure',
97  *  'lambda_and_pressure', 'temperature_lambda_pressure'?; Let's wait
98  *  until we feel better about the pressure control methods giving
99  *  exact ensembles.  Right now, we assume constant pressure */
100 static const char* erename[ereNR] = { "temperature", "lambda", "end_single_marker",
101                                       "temperature and lambda" };
102
103 //! Working data for replica exchange.
104 struct gmx_repl_ex
105 {
106     //! Replica ID
107     int repl;
108     //! Total number of replica
109     int nrepl;
110     //! Temperature
111     real temp;
112     //! Replica exchange type from ere enum
113     int type;
114     //! Quantity, e.g. temperature or lambda; first index is ere, second index is replica ID
115     real** q;
116     //! Use constant pressure and temperature
117     gmx_bool bNPT;
118     //! Replica pressures
119     real* pres;
120     //! Replica indices
121     int* ind;
122     //! Used for keeping track of all the replica swaps
123     int* allswaps;
124     //! Replica exchange interval (number of steps)
125     int nst;
126     //! Number of exchanges per interval
127     int nex;
128     //! Random seed
129     int seed;
130     //! Number of even and odd replica change attempts
131     int nattempt[2];
132     //! Sum of probabilities
133     real* prob_sum;
134     //! Number of moves between replicas i and j
135     int** nmoves;
136     //! i-th element of the array is the number of exchanges between replica i-1 and i
137     int* nexchange;
138
139     /*! \brief Helper arrays for replica exchange; allocated here
140      * so they don't have to be allocated each time */
141     //! \{
142     int*      destinations;
143     int**     cyclic;
144     int**     order;
145     int*      tmpswap;
146     gmx_bool* incycle;
147     gmx_bool* bEx;
148     //! \}
149
150     //! Helper arrays to hold the quantities that are exchanged.
151     //! \{
152     real*  prob;
153     real*  Epot;
154     real*  beta;
155     real*  Vol;
156     real** de;
157     //! \}
158 };
159
160 // TODO We should add Doxygen here some time.
161 //! \cond
162
163 static gmx_bool repl_quantity(const gmx_multisim_t* ms, struct gmx_repl_ex* re, int ere, real q)
164 {
165     real*    qall;
166     gmx_bool bDiff;
167     int      s;
168
169     snew(qall, ms->numSimulations_);
170     qall[re->repl] = q;
171     gmx_sum_sim(ms->numSimulations_, qall, ms);
172
173     bDiff = FALSE;
174     for (s = 1; s < ms->numSimulations_; s++)
175     {
176         if (qall[s] != qall[0])
177         {
178             bDiff = TRUE;
179         }
180     }
181
182     if (bDiff)
183     {
184         /* Set the replica exchange type and quantities */
185         re->type = ere;
186
187         snew(re->q[ere], re->nrepl);
188         for (s = 0; s < ms->numSimulations_; s++)
189         {
190             re->q[ere][s] = qall[s];
191         }
192     }
193     sfree(qall);
194     return bDiff;
195 }
196
197 gmx_repl_ex_t init_replica_exchange(FILE*                            fplog,
198                                     const gmx_multisim_t*            ms,
199                                     int                              numAtomsInSystem,
200                                     const t_inputrec*                ir,
201                                     const ReplicaExchangeParameters& replExParams)
202 {
203     real                pres;
204     int                 i, j;
205     struct gmx_repl_ex* re;
206     gmx_bool            bTemp;
207     gmx_bool            bLambda = FALSE;
208
209     fprintf(fplog, "\nInitializing Replica Exchange\n");
210
211     if (!isMultiSim(ms) || ms->numSimulations_ == 1)
212     {
213         gmx_fatal(FARGS,
214                   "Nothing to exchange with only one replica, maybe you forgot to set the "
215                   "-multidir option of mdrun?");
216     }
217     if (replExParams.numExchanges < 0)
218     {
219         gmx_fatal(FARGS, "Replica exchange number of exchanges needs to be positive");
220     }
221
222     if (!EI_DYNAMICS(ir->eI))
223     {
224         gmx_fatal(FARGS, "Replica exchange is only supported by dynamical simulations");
225         /* Note that PAR(cr) is defined by cr->nnodes > 1, which is
226          * distinct from isMultiSim(ms). A multi-simulation only runs
227          * with real MPI parallelism, but this does not imply PAR(cr)
228          * is true!
229          *
230          * Since we are using a dynamical integrator, the only
231          * decomposition is DD, so PAR(cr) and DOMAINDECOMP(cr) are
232          * synonymous. The only way for cr->nnodes > 1 to be true is
233          * if we are using DD. */
234     }
235
236     snew(re, 1);
237
238     re->repl  = ms->simulationIndex_;
239     re->nrepl = ms->numSimulations_;
240     snew(re->q, ereENDSINGLE);
241
242     fprintf(fplog, "Repl  There are %d replicas:\n", re->nrepl);
243
244     /* We only check that the number of atoms in the systms match.
245      * This, of course, do not guarantee that the systems are the same,
246      * but it does guarantee that we can perform replica exchange.
247      */
248     check_multi_int(fplog, ms, numAtomsInSystem, "the number of atoms", FALSE);
249     check_multi_int(fplog, ms, ir->eI, "the integrator", FALSE);
250     check_multi_int64(fplog, ms, ir->init_step + ir->nsteps, "init_step+nsteps", FALSE);
251     const int nst = replExParams.exchangeInterval;
252     check_multi_int64(fplog, ms, (ir->init_step + nst - 1) / nst,
253                       "first exchange step: init_step/-replex", FALSE);
254     check_multi_int(fplog, ms, ir->etc, "the temperature coupling", FALSE);
255     check_multi_int(fplog, ms, ir->opts.ngtc, "the number of temperature coupling groups", FALSE);
256     check_multi_int(fplog, ms, ir->epc, "the pressure coupling", FALSE);
257     check_multi_int(fplog, ms, ir->efep, "free energy", FALSE);
258     check_multi_int(fplog, ms, ir->fepvals->n_lambda, "number of lambda states", FALSE);
259
260     re->temp = ir->opts.ref_t[0];
261     for (i = 1; (i < ir->opts.ngtc); i++)
262     {
263         if (ir->opts.ref_t[i] != re->temp)
264         {
265             fprintf(fplog,
266                     "\nWARNING: The temperatures of the different temperature coupling groups are "
267                     "not identical\n\n");
268             fprintf(stderr,
269                     "\nWARNING: The temperatures of the different temperature coupling groups are "
270                     "not identical\n\n");
271         }
272     }
273
274     re->type = -1;
275     bTemp    = repl_quantity(ms, re, ereTEMP, re->temp);
276     if (ir->efep != efepNO)
277     {
278         bLambda = repl_quantity(ms, re, ereLAMBDA, static_cast<real>(ir->fepvals->init_fep_state));
279     }
280     if (re->type == -1) /* nothing was assigned */
281     {
282         gmx_fatal(FARGS,
283                   "The properties of the %d systems are all the same, there is nothing to exchange",
284                   re->nrepl);
285     }
286     if (bLambda && bTemp)
287     {
288         re->type = ereTL;
289     }
290
291     if (bTemp)
292     {
293         please_cite(fplog, "Sugita1999a");
294         if (ir->epc != epcNO)
295         {
296             re->bNPT = TRUE;
297             fprintf(fplog, "Repl  Using Constant Pressure REMD.\n");
298             please_cite(fplog, "Okabe2001a");
299         }
300         if (ir->etc == etcBERENDSEN)
301         {
302             gmx_fatal(FARGS,
303                       "REMD with the %s thermostat does not produce correct potential energy "
304                       "distributions, consider using the %s thermostat instead",
305                       ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
306         }
307     }
308     if (bLambda)
309     {
310         if (ir->fepvals->delta_lambda != 0) /* check this? */
311         {
312             gmx_fatal(FARGS, "delta_lambda is not zero");
313         }
314     }
315     if (re->bNPT)
316     {
317         snew(re->pres, re->nrepl);
318         if (ir->epct == epctSURFACETENSION)
319         {
320             pres = ir->ref_p[ZZ][ZZ];
321         }
322         else
323         {
324             pres = 0;
325             j    = 0;
326             for (i = 0; i < DIM; i++)
327             {
328                 if (ir->compress[i][i] != 0)
329                 {
330                     pres += ir->ref_p[i][i];
331                     j++;
332                 }
333             }
334             pres /= j;
335         }
336         re->pres[re->repl] = pres;
337         gmx_sum_sim(re->nrepl, re->pres, ms);
338     }
339
340     /* Make an index for increasing replica order */
341     /* only makes sense if one or the other is varying, not both!
342        if both are varying, we trust the order the person gave. */
343     snew(re->ind, re->nrepl);
344     for (i = 0; i < re->nrepl; i++)
345     {
346         re->ind[i] = i;
347     }
348
349     if (re->type < ereENDSINGLE)
350     {
351
352         for (i = 0; i < re->nrepl; i++)
353         {
354             for (j = i + 1; j < re->nrepl; j++)
355             {
356                 if (re->q[re->type][re->ind[j]] < re->q[re->type][re->ind[i]])
357                 {
358                     /* Unordered replicas are supposed to work, but there
359                      * is still an issues somewhere.
360                      * Note that at this point still re->ind[i]=i.
361                      */
362                     gmx_fatal(FARGS,
363                               "Replicas with indices %d < %d have %ss %g > %g, please order your "
364                               "replicas on increasing %s",
365                               i, j, erename[re->type], re->q[re->type][i], re->q[re->type][j],
366                               erename[re->type]);
367                 }
368                 else if (re->q[re->type][re->ind[j]] == re->q[re->type][re->ind[i]])
369                 {
370                     gmx_fatal(FARGS, "Two replicas have identical %ss", erename[re->type]);
371                 }
372             }
373         }
374     }
375
376     /* keep track of all the swaps, starting with the initial placement. */
377     snew(re->allswaps, re->nrepl);
378     for (i = 0; i < re->nrepl; i++)
379     {
380         re->allswaps[i] = re->ind[i];
381     }
382
383     switch (re->type)
384     {
385         case ereTEMP:
386             fprintf(fplog, "\nReplica exchange in temperature\n");
387             for (i = 0; i < re->nrepl; i++)
388             {
389                 fprintf(fplog, " %5.1f", re->q[re->type][re->ind[i]]);
390             }
391             fprintf(fplog, "\n");
392             break;
393         case ereLAMBDA:
394             fprintf(fplog, "\nReplica exchange in lambda\n");
395             for (i = 0; i < re->nrepl; i++)
396             {
397                 fprintf(fplog, " %3d", static_cast<int>(re->q[re->type][re->ind[i]]));
398             }
399             fprintf(fplog, "\n");
400             break;
401         case ereTL:
402             fprintf(fplog, "\nReplica exchange in temperature and lambda state\n");
403             for (i = 0; i < re->nrepl; i++)
404             {
405                 fprintf(fplog, " %5.1f", re->q[ereTEMP][re->ind[i]]);
406             }
407             fprintf(fplog, "\n");
408             for (i = 0; i < re->nrepl; i++)
409             {
410                 fprintf(fplog, " %5d", static_cast<int>(re->q[ereLAMBDA][re->ind[i]]));
411             }
412             fprintf(fplog, "\n");
413             break;
414         default: gmx_incons("Unknown replica exchange quantity");
415     }
416     if (re->bNPT)
417     {
418         fprintf(fplog, "\nRepl  p");
419         for (i = 0; i < re->nrepl; i++)
420         {
421             fprintf(fplog, " %5.2f", re->pres[re->ind[i]]);
422         }
423
424         for (i = 0; i < re->nrepl; i++)
425         {
426             if ((i > 0) && (re->pres[re->ind[i]] < re->pres[re->ind[i - 1]]))
427             {
428                 fprintf(fplog,
429                         "\nWARNING: The reference pressures decrease with increasing "
430                         "temperatures\n\n");
431                 fprintf(stderr,
432                         "\nWARNING: The reference pressures decrease with increasing "
433                         "temperatures\n\n");
434             }
435         }
436     }
437     re->nst = nst;
438     if (replExParams.randomSeed == -1)
439     {
440         if (isMasterSim(ms))
441         {
442             re->seed = static_cast<int>(gmx::makeRandomSeed());
443         }
444         else
445         {
446             re->seed = 0;
447         }
448         gmx_sumi_sim(1, &(re->seed), ms);
449     }
450     else
451     {
452         re->seed = replExParams.randomSeed;
453     }
454     fprintf(fplog, "\nReplica exchange interval: %d\n", re->nst);
455     fprintf(fplog, "\nReplica random seed: %d\n", re->seed);
456
457     re->nattempt[0] = 0;
458     re->nattempt[1] = 0;
459
460     snew(re->prob_sum, re->nrepl);
461     snew(re->nexchange, re->nrepl);
462     snew(re->nmoves, re->nrepl);
463     for (i = 0; i < re->nrepl; i++)
464     {
465         snew(re->nmoves[i], re->nrepl);
466     }
467     fprintf(fplog, "Replica exchange information below: ex and x = exchange, pr = probability\n");
468
469     /* generate space for the helper functions so we don't have to snew each time */
470
471     snew(re->destinations, re->nrepl);
472     snew(re->incycle, re->nrepl);
473     snew(re->tmpswap, re->nrepl);
474     snew(re->cyclic, re->nrepl);
475     snew(re->order, re->nrepl);
476     for (i = 0; i < re->nrepl; i++)
477     {
478         snew(re->cyclic[i], re->nrepl + 1);
479         snew(re->order[i], re->nrepl);
480     }
481     /* allocate space for the functions storing the data for the replicas */
482     /* not all of these arrays needed in all cases, but they don't take
483        up much space, since the max size is nrepl**2 */
484     snew(re->prob, re->nrepl);
485     snew(re->bEx, re->nrepl);
486     snew(re->beta, re->nrepl);
487     snew(re->Vol, re->nrepl);
488     snew(re->Epot, re->nrepl);
489     snew(re->de, re->nrepl);
490     for (i = 0; i < re->nrepl; i++)
491     {
492         snew(re->de[i], re->nrepl);
493     }
494     re->nex = replExParams.numExchanges;
495     return re;
496 }
497
498 static void exchange_reals(const gmx_multisim_t gmx_unused* ms, int gmx_unused b, real* v, int n)
499 {
500     real* buf;
501     int   i;
502
503     if (v)
504     {
505         snew(buf, n);
506 #if GMX_MPI
507         /*
508            MPI_Sendrecv(v,  n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
509            buf,n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
510            ms->mastersComm_,MPI_STATUS_IGNORE);
511          */
512         {
513             MPI_Request mpi_req;
514
515             MPI_Isend(v, n * sizeof(real), MPI_BYTE, MSRANK(ms, b), 0, ms->mastersComm_, &mpi_req);
516             MPI_Recv(buf, n * sizeof(real), MPI_BYTE, MSRANK(ms, b), 0, ms->mastersComm_, MPI_STATUS_IGNORE);
517             MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
518         }
519 #endif
520         for (i = 0; i < n; i++)
521         {
522             v[i] = buf[i];
523         }
524         sfree(buf);
525     }
526 }
527
528
529 static void exchange_doubles(const gmx_multisim_t gmx_unused* ms, int gmx_unused b, double* v, int n)
530 {
531     double* buf;
532     int     i;
533
534     if (v)
535     {
536         snew(buf, n);
537 #if GMX_MPI
538         /*
539            MPI_Sendrecv(v,  n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
540            buf,n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
541            ms->mastersComm_,MPI_STATUS_IGNORE);
542          */
543         {
544             MPI_Request mpi_req;
545
546             MPI_Isend(v, n * sizeof(double), MPI_BYTE, MSRANK(ms, b), 0, ms->mastersComm_, &mpi_req);
547             MPI_Recv(buf, n * sizeof(double), MPI_BYTE, MSRANK(ms, b), 0, ms->mastersComm_,
548                      MPI_STATUS_IGNORE);
549             MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
550         }
551 #endif
552         for (i = 0; i < n; i++)
553         {
554             v[i] = buf[i];
555         }
556         sfree(buf);
557     }
558 }
559
560 static void exchange_rvecs(const gmx_multisim_t gmx_unused* ms, int gmx_unused b, rvec* v, int n)
561 {
562     rvec* buf;
563     int   i;
564
565     if (v)
566     {
567         snew(buf, n);
568 #if GMX_MPI
569         /*
570            MPI_Sendrecv(v[0],  n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
571            buf[0],n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
572            ms->mastersComm_,MPI_STATUS_IGNORE);
573          */
574         {
575             MPI_Request mpi_req;
576
577             MPI_Isend(v[0], n * sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0, ms->mastersComm_, &mpi_req);
578             MPI_Recv(buf[0], n * sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0, ms->mastersComm_,
579                      MPI_STATUS_IGNORE);
580             MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
581         }
582 #endif
583         for (i = 0; i < n; i++)
584         {
585             copy_rvec(buf[i], v[i]);
586         }
587         sfree(buf);
588     }
589 }
590
591 static void exchange_state(const gmx_multisim_t* ms, int b, t_state* state)
592 {
593     /* When t_state changes, this code should be updated. */
594     int ngtc, nnhpres;
595     ngtc    = state->ngtc * state->nhchainlength;
596     nnhpres = state->nnhpres * state->nhchainlength;
597     exchange_rvecs(ms, b, state->box, DIM);
598     exchange_rvecs(ms, b, state->box_rel, DIM);
599     exchange_rvecs(ms, b, state->boxv, DIM);
600     exchange_reals(ms, b, &(state->veta), 1);
601     exchange_reals(ms, b, &(state->vol0), 1);
602     exchange_rvecs(ms, b, state->svir_prev, DIM);
603     exchange_rvecs(ms, b, state->fvir_prev, DIM);
604     exchange_rvecs(ms, b, state->pres_prev, DIM);
605     exchange_doubles(ms, b, state->nosehoover_xi.data(), ngtc);
606     exchange_doubles(ms, b, state->nosehoover_vxi.data(), ngtc);
607     exchange_doubles(ms, b, state->nhpres_xi.data(), nnhpres);
608     exchange_doubles(ms, b, state->nhpres_vxi.data(), nnhpres);
609     exchange_doubles(ms, b, state->therm_integral.data(), state->ngtc);
610     exchange_doubles(ms, b, &state->baros_integral, 1);
611     exchange_rvecs(ms, b, state->x.rvec_array(), state->natoms);
612     exchange_rvecs(ms, b, state->v.rvec_array(), state->natoms);
613 }
614
615 static void copy_state_serial(const t_state* src, t_state* dest)
616 {
617     if (dest != src)
618     {
619         /* Currently the local state is always a pointer to the global
620          * in serial, so we should never end up here.
621          * TODO: Implement a (trivial) t_state copy once converted to C++.
622          */
623         GMX_RELEASE_ASSERT(false, "State copying is currently not implemented in replica exchange");
624     }
625 }
626
627 static void scale_velocities(gmx::ArrayRef<gmx::RVec> velocities, real fac)
628 {
629     for (auto& v : velocities)
630     {
631         v *= fac;
632     }
633 }
634
635 static void print_transition_matrix(FILE* fplog, int n, int** nmoves, const int* nattempt)
636 {
637     int   i, j, ntot;
638     float Tprint;
639
640     ntot = nattempt[0] + nattempt[1];
641     fprintf(fplog, "\n");
642     fprintf(fplog, "Repl");
643     for (i = 0; i < n; i++)
644     {
645         fprintf(fplog, "    "); /* put the title closer to the center */
646     }
647     fprintf(fplog, "Empirical Transition Matrix\n");
648
649     fprintf(fplog, "Repl");
650     for (i = 0; i < n; i++)
651     {
652         fprintf(fplog, "%8d", (i + 1));
653     }
654     fprintf(fplog, "\n");
655
656     for (i = 0; i < n; i++)
657     {
658         fprintf(fplog, "Repl");
659         for (j = 0; j < n; j++)
660         {
661             Tprint = 0.0;
662             if (nmoves[i][j] > 0)
663             {
664                 Tprint = nmoves[i][j] / (2.0 * ntot);
665             }
666             fprintf(fplog, "%8.4f", Tprint);
667         }
668         fprintf(fplog, "%3d\n", i);
669     }
670 }
671
672 static void print_ind(FILE* fplog, const char* leg, int n, int* ind, const gmx_bool* bEx)
673 {
674     int i;
675
676     fprintf(fplog, "Repl %2s %2d", leg, ind[0]);
677     for (i = 1; i < n; i++)
678     {
679         fprintf(fplog, " %c %2d", (bEx != nullptr && bEx[i]) ? 'x' : ' ', ind[i]);
680     }
681     fprintf(fplog, "\n");
682 }
683
684 static void print_allswitchind(FILE* fplog, int n, int* pind, int* allswaps, int* tmpswap)
685 {
686     int i;
687
688     for (i = 0; i < n; i++)
689     {
690         tmpswap[i] = allswaps[i];
691     }
692     for (i = 0; i < n; i++)
693     {
694         allswaps[i] = tmpswap[pind[i]];
695     }
696
697     fprintf(fplog, "\nAccepted Exchanges:   ");
698     for (i = 0; i < n; i++)
699     {
700         fprintf(fplog, "%d ", pind[i]);
701     }
702     fprintf(fplog, "\n");
703
704     /* the "Order After Exchange" is the state label corresponding to the configuration that
705        started in state listed in order, i.e.
706
707        3 0 1 2
708
709        means that the:
710        configuration starting in simulation 3 is now in simulation 0,
711        configuration starting in simulation 0 is now in simulation 1,
712        configuration starting in simulation 1 is now in simulation 2,
713        configuration starting in simulation 2 is now in simulation 3
714      */
715     fprintf(fplog, "Order After Exchange: ");
716     for (i = 0; i < n; i++)
717     {
718         fprintf(fplog, "%d ", allswaps[i]);
719     }
720     fprintf(fplog, "\n\n");
721 }
722
723 static void print_prob(FILE* fplog, const char* leg, int n, real* prob)
724 {
725     int  i;
726     char buf[8];
727
728     fprintf(fplog, "Repl %2s ", leg);
729     for (i = 1; i < n; i++)
730     {
731         if (prob[i] >= 0)
732         {
733             sprintf(buf, "%4.2f", prob[i]);
734             fprintf(fplog, "  %3s", buf[0] == '1' ? "1.0" : buf + 1);
735         }
736         else
737         {
738             fprintf(fplog, "     ");
739         }
740     }
741     fprintf(fplog, "\n");
742 }
743
744 static void print_count(FILE* fplog, const char* leg, int n, int* count)
745 {
746     int i;
747
748     fprintf(fplog, "Repl %2s ", leg);
749     for (i = 1; i < n; i++)
750     {
751         fprintf(fplog, " %4d", count[i]);
752     }
753     fprintf(fplog, "\n");
754 }
755
756 static real calc_delta(FILE* fplog, gmx_bool bPrint, struct gmx_repl_ex* re, int a, int b, int ap, int bp)
757 {
758
759     real   ediff, dpV, delta = 0;
760     real*  Epot = re->Epot;
761     real*  Vol  = re->Vol;
762     real** de   = re->de;
763     real*  beta = re->beta;
764
765     /* Two cases; we are permuted and not.  In all cases, setting ap = a and bp = b will reduce
766        to the non permuted case */
767
768     switch (re->type)
769     {
770         case ereTEMP:
771             /*
772              * Okabe et. al. Chem. Phys. Lett. 335 (2001) 435-439
773              */
774             ediff = Epot[b] - Epot[a];
775             delta = -(beta[bp] - beta[ap]) * ediff;
776             break;
777         case ereLAMBDA:
778             /* two cases:  when we are permuted, and not.  */
779             /* non-permuted:
780                ediff =  E_new - E_old
781                      =  [H_b(x_a) + H_a(x_b)] - [H_b(x_b) + H_a(x_a)]
782                      =  [H_b(x_a) - H_a(x_a)] + [H_a(x_b) - H_b(x_b)]
783                      =  de[b][a] + de[a][b] */
784
785             /* permuted:
786                ediff =  E_new - E_old
787                      =  [H_bp(x_a) + H_ap(x_b)] - [H_bp(x_b) + H_ap(x_a)]
788                      =  [H_bp(x_a) - H_ap(x_a)] + [H_ap(x_b) - H_bp(x_b)]
789                      =  [H_bp(x_a) - H_a(x_a) + H_a(x_a) - H_ap(x_a)] + [H_ap(x_b) - H_b(x_b) + H_b(x_b) - H_bp(x_b)]
790                      =  [H_bp(x_a) - H_a(x_a)] - [H_ap(x_a) - H_a(x_a)] + [H_ap(x_b) - H_b(x_b)] - H_bp(x_b) - H_b(x_b)]
791                      =  (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b])    */
792             /* but, in the current code implementation, we flip configurations, not indices . . .
793                So let's examine that.
794                      =  [H_b(x_ap) - H_a(x_a)] - [H_a(x_ap) - H_a(x_a)] + [H_a(x_bp) - H_b(x_b)] - H_b(x_bp) - H_b(x_b)]
795                      =  [H_b(x_ap) - H_a(x_ap)]  + [H_a(x_bp) - H_b(x_pb)]
796                      = (de[b][ap] - de[a][ap]) + (de[a][bp] - de[b][bp]
797                      So, if we exchange b<=> bp and a<=> ap, we return to the same result.
798                      So the simple solution is to flip the
799                      position of perturbed and original indices in the tests.
800              */
801
802             ediff = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]);
803             delta = ediff * beta[a]; /* assume all same temperature in this case */
804             break;
805         case ereTL:
806             /* not permuted:  */
807             /* delta =  reduced E_new - reduced E_old
808                      =  [beta_b H_b(x_a) + beta_a H_a(x_b)] - [beta_b H_b(x_b) + beta_a H_a(x_a)]
809                      =  [beta_b H_b(x_a) - beta_a H_a(x_a)] + [beta_a H_a(x_b) - beta_b H_b(x_b)]
810                      =  [beta_b dH_b(x_a) + beta_b H_a(x_a) - beta_a H_a(x_a)] +
811                         [beta_a dH_a(x_b) + beta_a H_b(x_b) - beta_b H_b(x_b)]
812                      =  [beta_b dH_b(x_a) + [beta_a dH_a(x_b) +
813                         beta_b (H_a(x_a) - H_b(x_b)]) - beta_a (H_a(x_a) - H_b(x_b))
814                      =  beta_b dH_b(x_a) + beta_a dH_a(x_b) - (beta_b - beta_a)(H_b(x_b) - H_a(x_a) */
815             /* delta = beta[b]*de[b][a] + beta[a]*de[a][b] - (beta[b] - beta[a])*(Epot[b] - Epot[a]; */
816             /* permuted (big breath!) */
817             /*   delta =  reduced E_new - reduced E_old
818                      =  [beta_bp H_bp(x_a) + beta_ap H_ap(x_b)] - [beta_bp H_bp(x_b) + beta_ap H_ap(x_a)]
819                      =  [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
820                      =  [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
821                         - beta_pb H_a(x_a) + beta_ap H_a(x_a) + beta_pb H_a(x_a) - beta_ap H_a(x_a)
822                         - beta_ap H_b(x_b) + beta_bp H_b(x_b) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
823                      =  [(beta_bp H_bp(x_a) - beta_bp H_a(x_a)) - (beta_ap H_ap(x_a) - beta_ap H_a(x_a))] +
824                         [(beta_ap H_ap(x_b)  - beta_ap H_b(x_b)) - (beta_bp H_bp(x_b) - beta_bp H_b(x_b))]
825              + beta_pb H_a(x_a) - beta_ap H_a(x_a) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
826                      =  [beta_bp (H_bp(x_a) - H_a(x_a)) - beta_ap (H_ap(x_a) - H_a(x_a))] +
827                         [beta_ap (H_ap(x_b) - H_b(x_b)) - beta_bp (H_bp(x_b) - H_b(x_b))]
828              + beta_pb (H_a(x_a) - H_b(x_b))  - beta_ap (H_a(x_a) - H_b(x_b))
829                      =  ([beta_bp de[bp][a] - beta_ap de[ap][a]) + beta_ap de[ap][b]  - beta_bp de[bp][b])
830              + (beta_pb-beta_ap)(H_a(x_a) - H_b(x_b))  */
831             delta = beta[bp] * (de[bp][a] - de[bp][b]) + beta[ap] * (de[ap][b] - de[ap][a])
832                     - (beta[bp] - beta[ap]) * (Epot[b] - Epot[a]);
833             break;
834         default: gmx_incons("Unknown replica exchange quantity");
835     }
836     if (bPrint)
837     {
838         fprintf(fplog, "Repl %d <-> %d  dE_term = %10.3e (kT)\n", a, b, delta);
839     }
840     if (re->bNPT)
841     {
842         /* revist the calculation for 5.0.  Might be some improvements. */
843         dpV = (beta[ap] * re->pres[ap] - beta[bp] * re->pres[bp]) * (Vol[b] - Vol[a]) / PRESFAC;
844         if (bPrint)
845         {
846             fprintf(fplog, "  dpV = %10.3e  d = %10.3e\n", dpV, delta + dpV);
847         }
848         delta += dpV;
849     }
850     return delta;
851 }
852
853 static void test_for_replica_exchange(FILE*                 fplog,
854                                       const gmx_multisim_t* ms,
855                                       struct gmx_repl_ex*   re,
856                                       const gmx_enerdata_t* enerd,
857                                       real                  vol,
858                                       int64_t               step,
859                                       real                  time)
860 {
861     int                                m, i, j, a, b, ap, bp, i0, i1, tmp;
862     real                               delta = 0;
863     gmx_bool                           bPrint, bMultiEx;
864     gmx_bool*                          bEx      = re->bEx;
865     real*                              prob     = re->prob;
866     int*                               pind     = re->destinations; /* permuted index */
867     gmx_bool                           bEpot    = FALSE;
868     gmx_bool                           bDLambda = FALSE;
869     gmx_bool                           bVol     = FALSE;
870     gmx::ThreeFry2x64<64>              rng(re->seed, gmx::RandomDomain::ReplicaExchange);
871     gmx::UniformRealDistribution<real> uniformRealDist;
872     gmx::UniformIntDistribution<int>   uniformNreplDist(0, re->nrepl - 1);
873
874     bMultiEx = (re->nex > 1); /* multiple exchanges at each state */
875     fprintf(fplog, "Replica exchange at step %" PRId64 " time %.5f\n", step, time);
876
877     if (re->bNPT)
878     {
879         for (i = 0; i < re->nrepl; i++)
880         {
881             re->Vol[i] = 0;
882         }
883         bVol              = TRUE;
884         re->Vol[re->repl] = vol;
885     }
886     if ((re->type == ereTEMP || re->type == ereTL))
887     {
888         for (i = 0; i < re->nrepl; i++)
889         {
890             re->Epot[i] = 0;
891         }
892         bEpot              = TRUE;
893         re->Epot[re->repl] = enerd->term[F_EPOT];
894         /* temperatures of different states*/
895         for (i = 0; i < re->nrepl; i++)
896         {
897             re->beta[i] = 1.0 / (re->q[ereTEMP][i] * BOLTZ);
898         }
899     }
900     else
901     {
902         for (i = 0; i < re->nrepl; i++)
903         {
904             re->beta[i] = 1.0 / (re->temp * BOLTZ); /* we have a single temperature */
905         }
906     }
907     if (re->type == ereLAMBDA || re->type == ereTL)
908     {
909         bDLambda = TRUE;
910         /* lambda differences. */
911         /* de[i][j] is the energy of the jth simulation in the ith Hamiltonian
912            minus the energy of the jth simulation in the jth Hamiltonian */
913         for (i = 0; i < re->nrepl; i++)
914         {
915             for (j = 0; j < re->nrepl; j++)
916             {
917                 re->de[i][j] = 0;
918             }
919         }
920         for (i = 0; i < re->nrepl; i++)
921         {
922             re->de[i][re->repl] = enerd->foreignLambdaTerms.deltaH(re->q[ereLAMBDA][i]);
923         }
924     }
925
926     /* now actually do the communication */
927     if (bVol)
928     {
929         gmx_sum_sim(re->nrepl, re->Vol, ms);
930     }
931     if (bEpot)
932     {
933         gmx_sum_sim(re->nrepl, re->Epot, ms);
934     }
935     if (bDLambda)
936     {
937         for (i = 0; i < re->nrepl; i++)
938         {
939             gmx_sum_sim(re->nrepl, re->de[i], ms);
940         }
941     }
942
943     /* make a duplicate set of indices for shuffling */
944     for (i = 0; i < re->nrepl; i++)
945     {
946         pind[i] = re->ind[i];
947     }
948
949     rng.restart(step, 0);
950
951     if (bMultiEx)
952     {
953         /* multiple random switch exchange */
954         int nself = 0;
955
956
957         for (i = 0; i < re->nex + nself; i++)
958         {
959             // For now this is superfluous, but just in case we ever add more
960             // calls in different branches it is safer to always reset the distribution.
961             uniformNreplDist.reset();
962
963             /* randomly select a pair  */
964             /* in theory, could reduce this by identifying only which switches had a nonneglibible
965                probability of occurring (log p > -100) and only operate on those switches */
966             /* find out which state it is from, and what label that state currently has. Likely
967                more work that useful. */
968             i0 = uniformNreplDist(rng);
969             i1 = uniformNreplDist(rng);
970             if (i0 == i1)
971             {
972                 nself++;
973                 continue; /* self-exchange, back up and do it again */
974             }
975
976             a  = re->ind[i0]; /* what are the indices of these states? */
977             b  = re->ind[i1];
978             ap = pind[i0];
979             bp = pind[i1];
980
981             bPrint = FALSE; /* too noisy */
982             /* calculate the energy difference */
983             /* if the code changes to flip the STATES, rather than the configurations,
984                use the commented version of the code */
985             /* delta = calc_delta(fplog,bPrint,re,a,b,ap,bp); */
986             delta = calc_delta(fplog, bPrint, re, ap, bp, a, b);
987
988             /* we actually only use the first space in the prob and bEx array,
989                since there are actually many switches between pairs. */
990
991             if (delta <= 0)
992             {
993                 /* accepted */
994                 prob[0] = 1;
995                 bEx[0]  = TRUE;
996             }
997             else
998             {
999                 if (delta > c_probabilityCutoff)
1000                 {
1001                     prob[0] = 0;
1002                 }
1003                 else
1004                 {
1005                     prob[0] = exp(-delta);
1006                 }
1007                 // roll a number to determine if accepted. For now it is superfluous to
1008                 // reset, but just in case we ever add more calls in different branches
1009                 // it is safer to always reset the distribution.
1010                 uniformRealDist.reset();
1011                 bEx[0] = uniformRealDist(rng) < prob[0];
1012             }
1013             re->prob_sum[0] += prob[0];
1014
1015             if (bEx[0])
1016             {
1017                 /* swap the states */
1018                 tmp      = pind[i0];
1019                 pind[i0] = pind[i1];
1020                 pind[i1] = tmp;
1021             }
1022         }
1023         re->nattempt[0]++; /* keep track of total permutation trials here */
1024         print_allswitchind(fplog, re->nrepl, pind, re->allswaps, re->tmpswap);
1025     }
1026     else
1027     {
1028         /* standard nearest neighbor replica exchange */
1029
1030         m = (step / re->nst) % 2;
1031         for (i = 1; i < re->nrepl; i++)
1032         {
1033             a = re->ind[i - 1];
1034             b = re->ind[i];
1035
1036             bPrint = (re->repl == a || re->repl == b);
1037             if (i % 2 == m)
1038             {
1039                 delta = calc_delta(fplog, bPrint, re, a, b, a, b);
1040                 if (delta <= 0)
1041                 {
1042                     /* accepted */
1043                     prob[i] = 1;
1044                     bEx[i]  = TRUE;
1045                 }
1046                 else
1047                 {
1048                     if (delta > c_probabilityCutoff)
1049                     {
1050                         prob[i] = 0;
1051                     }
1052                     else
1053                     {
1054                         prob[i] = exp(-delta);
1055                     }
1056                     // roll a number to determine if accepted. For now it is superfluous to
1057                     // reset, but just in case we ever add more calls in different branches
1058                     // it is safer to always reset the distribution.
1059                     uniformRealDist.reset();
1060                     bEx[i] = uniformRealDist(rng) < prob[i];
1061                 }
1062                 re->prob_sum[i] += prob[i];
1063
1064                 if (bEx[i])
1065                 {
1066                     /* swap these two */
1067                     tmp         = pind[i - 1];
1068                     pind[i - 1] = pind[i];
1069                     pind[i]     = tmp;
1070                     re->nexchange[i]++; /* statistics for back compatibility */
1071                 }
1072             }
1073             else
1074             {
1075                 prob[i] = -1;
1076                 bEx[i]  = FALSE;
1077             }
1078         }
1079         /* print some statistics */
1080         print_ind(fplog, "ex", re->nrepl, re->ind, bEx);
1081         print_prob(fplog, "pr", re->nrepl, prob);
1082         fprintf(fplog, "\n");
1083         re->nattempt[m]++;
1084     }
1085
1086     /* record which moves were made and accepted */
1087     for (i = 0; i < re->nrepl; i++)
1088     {
1089         re->nmoves[re->ind[i]][pind[i]] += 1;
1090         re->nmoves[pind[i]][re->ind[i]] += 1;
1091     }
1092     fflush(fplog); /* make sure we can see what the last exchange was */
1093 }
1094
1095 static void cyclic_decomposition(const int* destinations, int** cyclic, gmx_bool* incycle, const int nrepl, int* nswap)
1096 {
1097
1098     int i, j, c, p;
1099     int maxlen = 1;
1100     for (i = 0; i < nrepl; i++)
1101     {
1102         incycle[i] = FALSE;
1103     }
1104     for (i = 0; i < nrepl; i++) /* one cycle for each replica */
1105     {
1106         if (incycle[i])
1107         {
1108             cyclic[i][0] = -1;
1109             continue;
1110         }
1111         cyclic[i][0] = i;
1112         incycle[i]   = TRUE;
1113         c            = 1;
1114         p            = i;
1115         for (j = 0; j < nrepl; j++) /* potentially all cycles are part, but we will break first */
1116         {
1117             p = destinations[p]; /* start permuting */
1118             if (p == i)
1119             {
1120                 cyclic[i][c] = -1;
1121                 if (c > maxlen)
1122                 {
1123                     maxlen = c;
1124                 }
1125                 break; /* we've reached the original element, the cycle is complete, and we marked the end. */
1126             }
1127             else
1128             {
1129                 cyclic[i][c] = p; /* each permutation gives a new member of the cycle */
1130                 incycle[p]   = TRUE;
1131                 c++;
1132             }
1133         }
1134     }
1135     *nswap = maxlen - 1;
1136
1137     if (debug)
1138     {
1139         for (i = 0; i < nrepl; i++)
1140         {
1141             fprintf(debug, "Cycle %d:", i);
1142             for (j = 0; j < nrepl; j++)
1143             {
1144                 if (cyclic[i][j] < 0)
1145                 {
1146                     break;
1147                 }
1148                 fprintf(debug, "%2d", cyclic[i][j]);
1149             }
1150             fprintf(debug, "\n");
1151         }
1152         fflush(debug);
1153     }
1154 }
1155
1156 static void compute_exchange_order(int** cyclic, int** order, const int nrepl, const int maxswap)
1157 {
1158     int i, j;
1159
1160     for (j = 0; j < maxswap; j++)
1161     {
1162         for (i = 0; i < nrepl; i++)
1163         {
1164             if (cyclic[i][j + 1] >= 0)
1165             {
1166                 order[cyclic[i][j + 1]][j] = cyclic[i][j];
1167                 order[cyclic[i][j]][j]     = cyclic[i][j + 1];
1168             }
1169         }
1170         for (i = 0; i < nrepl; i++)
1171         {
1172             if (order[i][j] < 0)
1173             {
1174                 order[i][j] = i; /* if it's not exchanging, it should stay this round*/
1175             }
1176         }
1177     }
1178
1179     if (debug)
1180     {
1181         fprintf(debug, "Replica Exchange Order\n");
1182         for (i = 0; i < nrepl; i++)
1183         {
1184             fprintf(debug, "Replica %d:", i);
1185             for (j = 0; j < maxswap; j++)
1186             {
1187                 if (order[i][j] < 0)
1188                 {
1189                     break;
1190                 }
1191                 fprintf(debug, "%2d", order[i][j]);
1192             }
1193             fprintf(debug, "\n");
1194         }
1195         fflush(debug);
1196     }
1197 }
1198
1199 static void prepare_to_do_exchange(struct gmx_repl_ex* re, const int replica_id, int* maxswap, gmx_bool* bThisReplicaExchanged)
1200 {
1201     int i, j;
1202     /* Hold the cyclic decomposition of the (multiple) replica
1203      * exchange. */
1204     gmx_bool bAnyReplicaExchanged = FALSE;
1205     *bThisReplicaExchanged        = FALSE;
1206
1207     for (i = 0; i < re->nrepl; i++)
1208     {
1209         if (re->destinations[i] != re->ind[i])
1210         {
1211             /* only mark as exchanged if the index has been shuffled */
1212             bAnyReplicaExchanged = TRUE;
1213             break;
1214         }
1215     }
1216     if (bAnyReplicaExchanged)
1217     {
1218         /* reinitialize the placeholder arrays */
1219         for (i = 0; i < re->nrepl; i++)
1220         {
1221             for (j = 0; j < re->nrepl; j++)
1222             {
1223                 re->cyclic[i][j] = -1;
1224                 re->order[i][j]  = -1;
1225             }
1226         }
1227
1228         /* Identify the cyclic decomposition of the permutation (very
1229          * fast if neighbor replica exchange). */
1230         cyclic_decomposition(re->destinations, re->cyclic, re->incycle, re->nrepl, maxswap);
1231
1232         /* Now translate the decomposition into a replica exchange
1233          * order at each step. */
1234         compute_exchange_order(re->cyclic, re->order, re->nrepl, *maxswap);
1235
1236         /* Did this replica do any exchange at any point? */
1237         for (j = 0; j < *maxswap; j++)
1238         {
1239             if (replica_id != re->order[replica_id][j])
1240             {
1241                 *bThisReplicaExchanged = TRUE;
1242                 break;
1243             }
1244         }
1245     }
1246 }
1247
1248 gmx_bool replica_exchange(FILE*                 fplog,
1249                           const t_commrec*      cr,
1250                           const gmx_multisim_t* ms,
1251                           struct gmx_repl_ex*   re,
1252                           t_state*              state,
1253                           const gmx_enerdata_t* enerd,
1254                           t_state*              state_local,
1255                           int64_t               step,
1256                           real                  time)
1257 {
1258     int j;
1259     int replica_id = 0;
1260     int exchange_partner;
1261     int maxswap = 0;
1262     /* Number of rounds of exchanges needed to deal with any multiple
1263      * exchanges. */
1264     /* Where each replica ends up after the exchange attempt(s). */
1265     /* The order in which multiple exchanges will occur. */
1266     gmx_bool bThisReplicaExchanged = FALSE;
1267
1268     if (MASTER(cr))
1269     {
1270         replica_id = re->repl;
1271         test_for_replica_exchange(fplog, ms, re, enerd, det(state_local->box), step, time);
1272         prepare_to_do_exchange(re, replica_id, &maxswap, &bThisReplicaExchanged);
1273     }
1274     /* Do intra-simulation broadcast so all processors belonging to
1275      * each simulation know whether they need to participate in
1276      * collecting the state. Otherwise, they might as well get on with
1277      * the next thing to do. */
1278     if (DOMAINDECOMP(cr))
1279     {
1280 #if GMX_MPI
1281         MPI_Bcast(&bThisReplicaExchanged, sizeof(gmx_bool), MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
1282 #endif
1283     }
1284
1285     if (bThisReplicaExchanged)
1286     {
1287         /* Exchange the states */
1288         /* Collect the global state on the master node */
1289         if (DOMAINDECOMP(cr))
1290         {
1291             dd_collect_state(cr->dd, state_local, state);
1292         }
1293         else
1294         {
1295             copy_state_serial(state_local, state);
1296         }
1297
1298         if (MASTER(cr))
1299         {
1300             /* There will be only one swap cycle with standard replica
1301              * exchange, but there may be multiple swap cycles if we
1302              * allow multiple swaps. */
1303
1304             for (j = 0; j < maxswap; j++)
1305             {
1306                 exchange_partner = re->order[replica_id][j];
1307
1308                 if (exchange_partner != replica_id)
1309                 {
1310                     /* Exchange the global states between the master nodes */
1311                     if (debug)
1312                     {
1313                         fprintf(debug, "Exchanging %d with %d\n", replica_id, exchange_partner);
1314                     }
1315                     exchange_state(ms, exchange_partner, state);
1316                 }
1317             }
1318             /* For temperature-type replica exchange, we need to scale
1319              * the velocities. */
1320             if (re->type == ereTEMP || re->type == ereTL)
1321             {
1322                 scale_velocities(state->v, std::sqrt(re->q[ereTEMP][replica_id]
1323                                                      / re->q[ereTEMP][re->destinations[replica_id]]));
1324             }
1325         }
1326
1327         /* With domain decomposition the global state is distributed later */
1328         if (!DOMAINDECOMP(cr))
1329         {
1330             /* Copy the global state to the local state data structure */
1331             copy_state_serial(state, state_local);
1332         }
1333     }
1334
1335     return bThisReplicaExchanged;
1336 }
1337
1338 void print_replica_exchange_statistics(FILE* fplog, struct gmx_repl_ex* re)
1339 {
1340     int i;
1341
1342     fprintf(fplog, "\nReplica exchange statistics\n");
1343
1344     if (re->nex == 0)
1345     {
1346         fprintf(fplog, "Repl  %d attempts, %d odd, %d even\n", re->nattempt[0] + re->nattempt[1],
1347                 re->nattempt[1], re->nattempt[0]);
1348
1349         fprintf(fplog, "Repl  average probabilities:\n");
1350         for (i = 1; i < re->nrepl; i++)
1351         {
1352             if (re->nattempt[i % 2] == 0)
1353             {
1354                 re->prob[i] = 0;
1355             }
1356             else
1357             {
1358                 re->prob[i] = re->prob_sum[i] / re->nattempt[i % 2];
1359             }
1360         }
1361         print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1362         print_prob(fplog, "", re->nrepl, re->prob);
1363
1364         fprintf(fplog, "Repl  number of exchanges:\n");
1365         print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1366         print_count(fplog, "", re->nrepl, re->nexchange);
1367
1368         fprintf(fplog, "Repl  average number of exchanges:\n");
1369         for (i = 1; i < re->nrepl; i++)
1370         {
1371             if (re->nattempt[i % 2] == 0)
1372             {
1373                 re->prob[i] = 0;
1374             }
1375             else
1376             {
1377                 re->prob[i] = (static_cast<real>(re->nexchange[i])) / re->nattempt[i % 2];
1378             }
1379         }
1380         print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1381         print_prob(fplog, "", re->nrepl, re->prob);
1382
1383         fprintf(fplog, "\n");
1384     }
1385     /* print the transition matrix */
1386     print_transition_matrix(fplog, re->nrepl, re->nmoves, re->nattempt);
1387 }
1388
1389 //! \endcond