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