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