2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
47 #include "gromacs/random/random.h"
48 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/math/units.h"
51 #include "gromacs/math/vec.h"
55 #include "gromacs/random/random.h"
57 #define PROBABILITYCUTOFF 100
58 /* we don't bother evaluating if events are more rare than exp(-100) = 3.7x10^-44 */
61 ereTEMP, ereLAMBDA, ereENDSINGLE, ereTL, ereNR
63 const char *erename[ereNR] = { "temperature", "lambda", "end_single_marker", "temperature and lambda"};
64 /* end_single_marker merely notes the end of single variable replica exchange. All types higher than
65 it are multiple replica exchange methods */
66 /* Eventually, should add 'pressure', 'temperature and pressure', 'lambda_and_pressure', 'temperature_lambda_pressure'?;
67 Let's wait until we feel better about the pressure control methods giving exact ensembles. Right now, we assume constant pressure */
69 typedef struct gmx_repl_ex
89 /* these are helper arrays for replica exchange; allocated here so they
90 don't have to be allocated each time */
98 /* helper arrays to hold the quantities that are exchanged */
107 static gmx_bool repl_quantity(const gmx_multisim_t *ms,
108 struct gmx_repl_ex *re, int ere, real q)
114 snew(qall, ms->nsim);
116 gmx_sum_sim(ms->nsim, qall, ms);
119 for (s = 1; s < ms->nsim; s++)
121 if (qall[s] != qall[0])
129 /* Set the replica exchange type and quantities */
132 snew(re->q[ere], re->nrepl);
133 for (s = 0; s < ms->nsim; s++)
135 re->q[ere][s] = qall[s];
142 gmx_repl_ex_t init_replica_exchange(FILE *fplog,
143 const gmx_multisim_t *ms,
144 const t_state *state,
145 const t_inputrec *ir,
146 int nst, int nex, int init_seed)
150 struct gmx_repl_ex *re;
152 gmx_bool bLambda = FALSE;
154 fprintf(fplog, "\nInitializing Replica Exchange\n");
156 if (ms == NULL || ms->nsim == 1)
158 gmx_fatal(FARGS, "Nothing to exchange with only one replica, maybe you forgot to set the -multi option of mdrun?");
160 if (!EI_DYNAMICS(ir->eI))
162 gmx_fatal(FARGS, "Replica exchange is only supported by dynamical simulations");
163 /* Note that PAR(cr) is defined by cr->nnodes > 1, which is
164 * distinct from MULTISIM(cr). A multi-simulation only runs
165 * with real MPI parallelism, but this does not imply PAR(cr)
168 * Since we are using a dynamical integrator, the only
169 * decomposition is DD, so PAR(cr) and DOMAINDECOMP(cr) are
170 * synonymous. The only way for cr->nnodes > 1 to be true is
171 * if we are using DD. */
177 re->nrepl = ms->nsim;
178 snew(re->q, ereENDSINGLE);
180 fprintf(fplog, "Repl There are %d replicas:\n", re->nrepl);
182 check_multi_int(fplog, ms, state->natoms, "the number of atoms", FALSE);
183 check_multi_int(fplog, ms, ir->eI, "the integrator", FALSE);
184 check_multi_int64(fplog, ms, ir->init_step+ir->nsteps, "init_step+nsteps", FALSE);
185 check_multi_int64(fplog, ms, (ir->init_step+nst-1)/nst,
186 "first exchange step: init_step/-replex", FALSE);
187 check_multi_int(fplog, ms, ir->etc, "the temperature coupling", FALSE);
188 check_multi_int(fplog, ms, ir->opts.ngtc,
189 "the number of temperature coupling groups", FALSE);
190 check_multi_int(fplog, ms, ir->epc, "the pressure coupling", FALSE);
191 check_multi_int(fplog, ms, ir->efep, "free energy", FALSE);
192 check_multi_int(fplog, ms, ir->fepvals->n_lambda, "number of lambda states", FALSE);
194 re->temp = ir->opts.ref_t[0];
195 for (i = 1; (i < ir->opts.ngtc); i++)
197 if (ir->opts.ref_t[i] != re->temp)
199 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
200 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
205 bTemp = repl_quantity(ms, re, ereTEMP, re->temp);
206 if (ir->efep != efepNO)
208 bLambda = repl_quantity(ms, re, ereLAMBDA, (real)ir->fepvals->init_fep_state);
210 if (re->type == -1) /* nothing was assigned */
212 gmx_fatal(FARGS, "The properties of the %d systems are all the same, there is nothing to exchange", re->nrepl);
214 if (bLambda && bTemp)
221 please_cite(fplog, "Sugita1999a");
222 if (ir->epc != epcNO)
225 fprintf(fplog, "Repl Using Constant Pressure REMD.\n");
226 please_cite(fplog, "Okabe2001a");
228 if (ir->etc == etcBERENDSEN)
230 gmx_fatal(FARGS, "REMD with the %s thermostat does not produce correct potential energy distributions, consider using the %s thermostat instead",
231 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
236 if (ir->fepvals->delta_lambda != 0) /* check this? */
238 gmx_fatal(FARGS, "delta_lambda is not zero");
243 snew(re->pres, re->nrepl);
244 if (ir->epct == epctSURFACETENSION)
246 pres = ir->ref_p[ZZ][ZZ];
252 for (i = 0; i < DIM; i++)
254 if (ir->compress[i][i] != 0)
256 pres += ir->ref_p[i][i];
262 re->pres[re->repl] = pres;
263 gmx_sum_sim(re->nrepl, re->pres, ms);
266 /* Make an index for increasing replica order */
267 /* only makes sense if one or the other is varying, not both!
268 if both are varying, we trust the order the person gave. */
269 snew(re->ind, re->nrepl);
270 for (i = 0; i < re->nrepl; i++)
275 if (re->type < ereENDSINGLE)
278 for (i = 0; i < re->nrepl; i++)
280 for (j = i+1; j < re->nrepl; j++)
282 if (re->q[re->type][re->ind[j]] < re->q[re->type][re->ind[i]])
285 re->ind[i] = re->ind[j];
288 else if (re->q[re->type][re->ind[j]] == re->q[re->type][re->ind[i]])
290 gmx_fatal(FARGS, "Two replicas have identical %ss", erename[re->type]);
296 /* keep track of all the swaps, starting with the initial placement. */
297 snew(re->allswaps, re->nrepl);
298 for (i = 0; i < re->nrepl; i++)
300 re->allswaps[i] = re->ind[i];
306 fprintf(fplog, "\nReplica exchange in temperature\n");
307 for (i = 0; i < re->nrepl; i++)
309 fprintf(fplog, " %5.1f", re->q[re->type][re->ind[i]]);
311 fprintf(fplog, "\n");
314 fprintf(fplog, "\nReplica exchange in lambda\n");
315 for (i = 0; i < re->nrepl; i++)
317 fprintf(fplog, " %3d", (int)re->q[re->type][re->ind[i]]);
319 fprintf(fplog, "\n");
322 fprintf(fplog, "\nReplica exchange in temperature and lambda state\n");
323 for (i = 0; i < re->nrepl; i++)
325 fprintf(fplog, " %5.1f", re->q[ereTEMP][re->ind[i]]);
327 fprintf(fplog, "\n");
328 for (i = 0; i < re->nrepl; i++)
330 fprintf(fplog, " %5d", (int)re->q[ereLAMBDA][re->ind[i]]);
332 fprintf(fplog, "\n");
335 gmx_incons("Unknown replica exchange quantity");
339 fprintf(fplog, "\nRepl p");
340 for (i = 0; i < re->nrepl; i++)
342 fprintf(fplog, " %5.2f", re->pres[re->ind[i]]);
345 for (i = 0; i < re->nrepl; i++)
347 if ((i > 0) && (re->pres[re->ind[i]] < re->pres[re->ind[i-1]]))
349 fprintf(fplog, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
350 fprintf(stderr, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
359 re->seed = (int)gmx_rng_make_seed();
365 gmx_sumi_sim(1, &(re->seed), ms);
369 re->seed = init_seed;
371 fprintf(fplog, "\nReplica exchange interval: %d\n", re->nst);
372 fprintf(fplog, "\nReplica random seed: %d\n", re->seed);
373 re->rng = gmx_rng_init(re->seed);
378 snew(re->prob_sum, re->nrepl);
379 snew(re->nexchange, re->nrepl);
380 snew(re->nmoves, re->nrepl);
381 for (i = 0; i < re->nrepl; i++)
383 snew(re->nmoves[i], re->nrepl);
385 fprintf(fplog, "Replica exchange information below: x=exchange, pr=probability\n");
387 /* generate space for the helper functions so we don't have to snew each time */
389 snew(re->destinations, re->nrepl);
390 snew(re->incycle, re->nrepl);
391 snew(re->tmpswap, re->nrepl);
392 snew(re->cyclic, re->nrepl);
393 snew(re->order, re->nrepl);
394 for (i = 0; i < re->nrepl; i++)
396 snew(re->cyclic[i], re->nrepl);
397 snew(re->order[i], re->nrepl);
399 /* allocate space for the functions storing the data for the replicas */
400 /* not all of these arrays needed in all cases, but they don't take
401 up much space, since the max size is nrepl**2 */
402 snew(re->prob, re->nrepl);
403 snew(re->bEx, re->nrepl);
404 snew(re->beta, re->nrepl);
405 snew(re->Vol, re->nrepl);
406 snew(re->Epot, re->nrepl);
407 snew(re->de, re->nrepl);
408 for (i = 0; i < re->nrepl; i++)
410 snew(re->de[i], re->nrepl);
416 static void exchange_reals(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, real *v, int n)
426 MPI_Sendrecv(v, n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
427 buf,n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
428 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
433 MPI_Isend(v, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
434 ms->mpi_comm_masters, &mpi_req);
435 MPI_Recv(buf, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
436 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
437 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
440 for (i = 0; i < n; i++)
449 static void exchange_ints(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, int *v, int n)
459 MPI_Sendrecv(v, n*sizeof(int),MPI_BYTE,MSRANK(ms,b),0,
460 buf,n*sizeof(int),MPI_BYTE,MSRANK(ms,b),0,
461 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
466 MPI_Isend(v, n*sizeof(int), MPI_BYTE, MSRANK(ms, b), 0,
467 ms->mpi_comm_masters, &mpi_req);
468 MPI_Recv(buf, n*sizeof(int), MPI_BYTE, MSRANK(ms, b), 0,
469 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
470 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
473 for (i = 0; i < n; i++)
481 static void exchange_doubles(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, double *v, int n)
491 MPI_Sendrecv(v, n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
492 buf,n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
493 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
498 MPI_Isend(v, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
499 ms->mpi_comm_masters, &mpi_req);
500 MPI_Recv(buf, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
501 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
502 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
505 for (i = 0; i < n; i++)
513 static void exchange_rvecs(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, rvec *v, int n)
523 MPI_Sendrecv(v[0], n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
524 buf[0],n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
525 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
530 MPI_Isend(v[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
531 ms->mpi_comm_masters, &mpi_req);
532 MPI_Recv(buf[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
533 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
534 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
537 for (i = 0; i < n; i++)
539 copy_rvec(buf[i], v[i]);
545 static void exchange_state(const gmx_multisim_t *ms, int b, t_state *state)
547 /* When t_state changes, this code should be updated. */
549 ngtc = state->ngtc * state->nhchainlength;
550 nnhpres = state->nnhpres* state->nhchainlength;
551 exchange_rvecs(ms, b, state->box, DIM);
552 exchange_rvecs(ms, b, state->box_rel, DIM);
553 exchange_rvecs(ms, b, state->boxv, DIM);
554 exchange_reals(ms, b, &(state->veta), 1);
555 exchange_reals(ms, b, &(state->vol0), 1);
556 exchange_rvecs(ms, b, state->svir_prev, DIM);
557 exchange_rvecs(ms, b, state->fvir_prev, DIM);
558 exchange_rvecs(ms, b, state->pres_prev, DIM);
559 exchange_doubles(ms, b, state->nosehoover_xi, ngtc);
560 exchange_doubles(ms, b, state->nosehoover_vxi, ngtc);
561 exchange_doubles(ms, b, state->nhpres_xi, nnhpres);
562 exchange_doubles(ms, b, state->nhpres_vxi, nnhpres);
563 exchange_doubles(ms, b, state->therm_integral, state->ngtc);
564 exchange_rvecs(ms, b, state->x, state->natoms);
565 exchange_rvecs(ms, b, state->v, state->natoms);
566 exchange_rvecs(ms, b, state->sd_X, state->natoms);
569 static void copy_rvecs(rvec *s, rvec *d, int n)
575 for (i = 0; i < n; i++)
577 copy_rvec(s[i], d[i]);
582 static void copy_doubles(const double *s, double *d, int n)
588 for (i = 0; i < n; i++)
595 static void copy_reals(const real *s, real *d, int n)
601 for (i = 0; i < n; i++)
608 static void copy_ints(const int *s, int *d, int n)
614 for (i = 0; i < n; i++)
621 #define scopy_rvecs(v, n) copy_rvecs(state->v, state_local->v, n);
622 #define scopy_doubles(v, n) copy_doubles(state->v, state_local->v, n);
623 #define scopy_reals(v, n) copy_reals(state->v, state_local->v, n);
624 #define scopy_ints(v, n) copy_ints(state->v, state_local->v, n);
626 static void copy_state_nonatomdata(t_state *state, t_state *state_local)
628 /* When t_state changes, this code should be updated. */
630 ngtc = state->ngtc * state->nhchainlength;
631 nnhpres = state->nnhpres* state->nhchainlength;
632 scopy_rvecs(box, DIM);
633 scopy_rvecs(box_rel, DIM);
634 scopy_rvecs(boxv, DIM);
635 state_local->veta = state->veta;
636 state_local->vol0 = state->vol0;
637 scopy_rvecs(svir_prev, DIM);
638 scopy_rvecs(fvir_prev, DIM);
639 scopy_rvecs(pres_prev, DIM);
640 scopy_doubles(nosehoover_xi, ngtc);
641 scopy_doubles(nosehoover_vxi, ngtc);
642 scopy_doubles(nhpres_xi, nnhpres);
643 scopy_doubles(nhpres_vxi, nnhpres);
644 scopy_doubles(therm_integral, state->ngtc);
645 scopy_rvecs(x, state->natoms);
646 scopy_rvecs(v, state->natoms);
647 scopy_rvecs(sd_X, state->natoms);
648 copy_ints(&(state->fep_state), &(state_local->fep_state), 1);
649 scopy_reals(lambda, efptNR);
652 static void scale_velocities(t_state *state, real fac)
658 for (i = 0; i < state->natoms; i++)
660 svmul(fac, state->v[i], state->v[i]);
665 static void print_transition_matrix(FILE *fplog, int n, int **nmoves, int *nattempt)
670 ntot = nattempt[0] + nattempt[1];
671 fprintf(fplog, "\n");
672 fprintf(fplog, "Repl");
673 for (i = 0; i < n; i++)
675 fprintf(fplog, " "); /* put the title closer to the center */
677 fprintf(fplog, "Empirical Transition Matrix\n");
679 fprintf(fplog, "Repl");
680 for (i = 0; i < n; i++)
682 fprintf(fplog, "%8d", (i+1));
684 fprintf(fplog, "\n");
686 for (i = 0; i < n; i++)
688 fprintf(fplog, "Repl");
689 for (j = 0; j < n; j++)
692 if (nmoves[i][j] > 0)
694 Tprint = nmoves[i][j]/(2.0*ntot);
696 fprintf(fplog, "%8.4f", Tprint);
698 fprintf(fplog, "%3d\n", i);
702 static void print_ind(FILE *fplog, const char *leg, int n, int *ind, gmx_bool *bEx)
706 fprintf(fplog, "Repl %2s %2d", leg, ind[0]);
707 for (i = 1; i < n; i++)
709 fprintf(fplog, " %c %2d", (bEx != 0 && bEx[i]) ? 'x' : ' ', ind[i]);
711 fprintf(fplog, "\n");
714 static void print_allswitchind(FILE *fplog, int n, int *pind, int *allswaps, int *tmpswap)
718 for (i = 0; i < n; i++)
720 tmpswap[i] = allswaps[i];
722 for (i = 0; i < n; i++)
724 allswaps[i] = tmpswap[pind[i]];
727 fprintf(fplog, "\nAccepted Exchanges: ");
728 for (i = 0; i < n; i++)
730 fprintf(fplog, "%d ", pind[i]);
732 fprintf(fplog, "\n");
734 /* the "Order After Exchange" is the state label corresponding to the configuration that
735 started in state listed in order, i.e.
740 configuration starting in simulation 3 is now in simulation 0,
741 configuration starting in simulation 0 is now in simulation 1,
742 configuration starting in simulation 1 is now in simulation 2,
743 configuration starting in simulation 2 is now in simulation 3
745 fprintf(fplog, "Order After Exchange: ");
746 for (i = 0; i < n; i++)
748 fprintf(fplog, "%d ", allswaps[i]);
750 fprintf(fplog, "\n\n");
753 static void print_prob(FILE *fplog, const char *leg, int n, real *prob)
758 fprintf(fplog, "Repl %2s ", leg);
759 for (i = 1; i < n; i++)
763 sprintf(buf, "%4.2f", prob[i]);
764 fprintf(fplog, " %3s", buf[0] == '1' ? "1.0" : buf+1);
771 fprintf(fplog, "\n");
774 static void print_count(FILE *fplog, const char *leg, int n, int *count)
778 fprintf(fplog, "Repl %2s ", leg);
779 for (i = 1; i < n; i++)
781 fprintf(fplog, " %4d", count[i]);
783 fprintf(fplog, "\n");
786 static real calc_delta(FILE *fplog, gmx_bool bPrint, struct gmx_repl_ex *re, int a, int b, int ap, int bp)
789 real ediff, dpV, delta = 0;
790 real *Epot = re->Epot;
793 real *beta = re->beta;
795 /* Two cases; we are permuted and not. In all cases, setting ap = a and bp = b will reduce
796 to the non permuted case */
802 * Okabe et. al. Chem. Phys. Lett. 335 (2001) 435-439
804 ediff = Epot[b] - Epot[a];
805 delta = -(beta[bp] - beta[ap])*ediff;
808 /* two cases: when we are permuted, and not. */
810 ediff = E_new - E_old
811 = [H_b(x_a) + H_a(x_b)] - [H_b(x_b) + H_a(x_a)]
812 = [H_b(x_a) - H_a(x_a)] + [H_a(x_b) - H_b(x_b)]
813 = de[b][a] + de[a][b] */
816 ediff = E_new - E_old
817 = [H_bp(x_a) + H_ap(x_b)] - [H_bp(x_b) + H_ap(x_a)]
818 = [H_bp(x_a) - H_ap(x_a)] + [H_ap(x_b) - H_bp(x_b)]
819 = [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)]
820 = [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)]
821 = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]) */
822 /* but, in the current code implementation, we flip configurations, not indices . . .
823 So let's examine that.
824 = [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)]
825 = [H_b(x_ap) - H_a(x_ap)] + [H_a(x_bp) - H_b(x_pb)]
826 = (de[b][ap] - de[a][ap]) + (de[a][bp] - de[b][bp]
827 So, if we exchange b<=> bp and a<=> ap, we return to the same result.
828 So the simple solution is to flip the
829 position of perturbed and original indices in the tests.
832 ediff = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]);
833 delta = ediff*beta[a]; /* assume all same temperature in this case */
837 /* delta = reduced E_new - reduced E_old
838 = [beta_b H_b(x_a) + beta_a H_a(x_b)] - [beta_b H_b(x_b) + beta_a H_a(x_a)]
839 = [beta_b H_b(x_a) - beta_a H_a(x_a)] + [beta_a H_a(x_b) - beta_b H_b(x_b)]
840 = [beta_b dH_b(x_a) + beta_b H_a(x_a) - beta_a H_a(x_a)] +
841 [beta_a dH_a(x_b) + beta_a H_b(x_b) - beta_b H_b(x_b)]
842 = [beta_b dH_b(x_a) + [beta_a dH_a(x_b) +
843 beta_b (H_a(x_a) - H_b(x_b)]) - beta_a (H_a(x_a) - H_b(x_b))
844 = beta_b dH_b(x_a) + beta_a dH_a(x_b) - (beta_b - beta_a)(H_b(x_b) - H_a(x_a) */
845 /* delta = beta[b]*de[b][a] + beta[a]*de[a][b] - (beta[b] - beta[a])*(Epot[b] - Epot[a]; */
846 /* permuted (big breath!) */
847 /* delta = reduced E_new - reduced E_old
848 = [beta_bp H_bp(x_a) + beta_ap H_ap(x_b)] - [beta_bp H_bp(x_b) + beta_ap H_ap(x_a)]
849 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
850 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
851 - beta_pb H_a(x_a) + beta_ap H_a(x_a) + beta_pb H_a(x_a) - beta_ap H_a(x_a)
852 - beta_ap H_b(x_b) + beta_bp H_b(x_b) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
853 = [(beta_bp H_bp(x_a) - beta_bp H_a(x_a)) - (beta_ap H_ap(x_a) - beta_ap H_a(x_a))] +
854 [(beta_ap H_ap(x_b) - beta_ap H_b(x_b)) - (beta_bp H_bp(x_b) - beta_bp H_b(x_b))]
855 + beta_pb H_a(x_a) - beta_ap H_a(x_a) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
856 = [beta_bp (H_bp(x_a) - H_a(x_a)) - beta_ap (H_ap(x_a) - H_a(x_a))] +
857 [beta_ap (H_ap(x_b) - H_b(x_b)) - beta_bp (H_bp(x_b) - H_b(x_b))]
858 + beta_pb (H_a(x_a) - H_b(x_b)) - beta_ap (H_a(x_a) - H_b(x_b))
859 = ([beta_bp de[bp][a] - beta_ap de[ap][a]) + beta_ap de[ap][b] - beta_bp de[bp][b])
860 + (beta_pb-beta_ap)(H_a(x_a) - H_b(x_b)) */
861 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]);
864 gmx_incons("Unknown replica exchange quantity");
868 fprintf(fplog, "Repl %d <-> %d dE_term = %10.3e (kT)\n", a, b, delta);
872 /* revist the calculation for 5.0. Might be some improvements. */
873 dpV = (beta[ap]*re->pres[ap]-beta[bp]*re->pres[bp])*(Vol[b]-Vol[a])/PRESFAC;
876 fprintf(fplog, " dpV = %10.3e d = %10.3e\nb", dpV, delta + dpV);
884 test_for_replica_exchange(FILE *fplog,
885 const gmx_multisim_t *ms,
886 struct gmx_repl_ex *re,
887 gmx_enerdata_t *enerd,
892 int m, i, j, a, b, ap, bp, i0, i1, tmp;
894 gmx_bool bPrint, bMultiEx;
895 gmx_bool *bEx = re->bEx;
896 real *prob = re->prob;
897 int *pind = re->destinations; /* permuted index */
898 gmx_bool bEpot = FALSE;
899 gmx_bool bDLambda = FALSE;
900 gmx_bool bVol = FALSE;
902 bMultiEx = (re->nex > 1); /* multiple exchanges at each state */
903 fprintf(fplog, "Replica exchange at step " "%" GMX_PRId64 " time %g\n", step, time);
907 for (i = 0; i < re->nrepl; i++)
912 re->Vol[re->repl] = vol;
914 if ((re->type == ereTEMP || re->type == ereTL))
916 for (i = 0; i < re->nrepl; i++)
921 re->Epot[re->repl] = enerd->term[F_EPOT];
922 /* temperatures of different states*/
923 for (i = 0; i < re->nrepl; i++)
925 re->beta[i] = 1.0/(re->q[ereTEMP][i]*BOLTZ);
930 for (i = 0; i < re->nrepl; i++)
932 re->beta[i] = 1.0/(re->temp*BOLTZ); /* we have a single temperature */
935 if (re->type == ereLAMBDA || re->type == ereTL)
938 /* lambda differences. */
939 /* de[i][j] is the energy of the jth simulation in the ith Hamiltonian
940 minus the energy of the jth simulation in the jth Hamiltonian */
941 for (i = 0; i < re->nrepl; i++)
943 for (j = 0; j < re->nrepl; j++)
948 for (i = 0; i < re->nrepl; i++)
950 re->de[i][re->repl] = (enerd->enerpart_lambda[(int)re->q[ereLAMBDA][i]+1]-enerd->enerpart_lambda[0]);
954 /* now actually do the communication */
957 gmx_sum_sim(re->nrepl, re->Vol, ms);
961 gmx_sum_sim(re->nrepl, re->Epot, ms);
965 for (i = 0; i < re->nrepl; i++)
967 gmx_sum_sim(re->nrepl, re->de[i], ms);
971 /* make a duplicate set of indices for shuffling */
972 for (i = 0; i < re->nrepl; i++)
974 pind[i] = re->ind[i];
979 /* multiple random switch exchange */
981 for (i = 0; i < re->nex + nself; i++)
985 gmx_rng_cycle_2uniform(step, i*2, re->seed, RND_SEED_REPLEX, rnd);
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*rnd[0]);
992 i1 = (int)(re->nrepl*rnd[1]);
996 continue; /* self-exchange, back up and do it again */
999 a = re->ind[i0]; /* what are the indices of these states? */
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);
1011 /* we actually only use the first space in the prob and bEx array,
1012 since there are actually many switches between pairs. */
1022 if (delta > PROBABILITYCUTOFF)
1028 prob[0] = exp(-delta);
1030 /* roll a number to determine if accepted */
1031 gmx_rng_cycle_2uniform(step, i*2+1, re->seed, RND_SEED_REPLEX, rnd);
1032 bEx[0] = rnd[0] < prob[0];
1034 re->prob_sum[0] += prob[0];
1038 /* swap the states */
1040 pind[i0] = pind[i1];
1044 re->nattempt[0]++; /* keep track of total permutation trials here */
1045 print_allswitchind(fplog, re->nrepl, pind, re->allswaps, re->tmpswap);
1049 /* standard nearest neighbor replica exchange */
1051 m = (step / re->nst) % 2;
1052 for (i = 1; i < re->nrepl; i++)
1057 bPrint = (re->repl == a || re->repl == b);
1060 delta = calc_delta(fplog, bPrint, re, a, b, a, b);
1071 if (delta > PROBABILITYCUTOFF)
1077 prob[i] = exp(-delta);
1079 /* roll a number to determine if accepted */
1080 gmx_rng_cycle_2uniform(step, i, re->seed, RND_SEED_REPLEX, rnd);
1081 bEx[i] = rnd[0] < prob[i];
1083 re->prob_sum[i] += prob[i];
1087 /* swap these two */
1089 pind[i-1] = pind[i];
1091 re->nexchange[i]++; /* statistics for back compatibility */
1100 /* print some statistics */
1101 print_ind(fplog, "ex", re->nrepl, re->ind, bEx);
1102 print_prob(fplog, "pr", re->nrepl, prob);
1103 fprintf(fplog, "\n");
1107 /* record which moves were made and accepted */
1108 for (i = 0; i < re->nrepl; i++)
1110 re->nmoves[re->ind[i]][pind[i]] += 1;
1111 re->nmoves[pind[i]][re->ind[i]] += 1;
1113 fflush(fplog); /* make sure we can see what the last exchange was */
1116 static void write_debug_x(t_state *state)
1122 for (i = 0; i < state->natoms; i += 10)
1124 fprintf(debug, "dx %5d %10.5f %10.5f %10.5f\n", i, state->x[i][XX], state->x[i][YY], state->x[i][ZZ]);
1130 cyclic_decomposition(const int *destinations,
1139 for (i = 0; i < nrepl; i++)
1143 for (i = 0; i < nrepl; i++) /* one cycle for each replica */
1154 for (j = 0; j < nrepl; j++) /* potentially all cycles are part, but we will break first */
1156 p = destinations[p]; /* start permuting */
1164 break; /* we've reached the original element, the cycle is complete, and we marked the end. */
1168 cyclic[i][c] = p; /* each permutation gives a new member of the cycle */
1174 *nswap = maxlen - 1;
1178 for (i = 0; i < nrepl; i++)
1180 fprintf(debug, "Cycle %d:", i);
1181 for (j = 0; j < nrepl; j++)
1183 if (cyclic[i][j] < 0)
1187 fprintf(debug, "%2d", cyclic[i][j]);
1189 fprintf(debug, "\n");
1196 compute_exchange_order(FILE *fplog,
1204 for (j = 0; j < maxswap; j++)
1206 for (i = 0; i < nrepl; i++)
1208 if (cyclic[i][j+1] >= 0)
1210 order[cyclic[i][j+1]][j] = cyclic[i][j];
1211 order[cyclic[i][j]][j] = cyclic[i][j+1];
1214 for (i = 0; i < nrepl; i++)
1216 if (order[i][j] < 0)
1218 order[i][j] = i; /* if it's not exchanging, it should stay this round*/
1225 fprintf(fplog, "Replica Exchange Order\n");
1226 for (i = 0; i < nrepl; i++)
1228 fprintf(fplog, "Replica %d:", i);
1229 for (j = 0; j < maxswap; j++)
1231 if (order[i][j] < 0)
1235 fprintf(debug, "%2d", order[i][j]);
1237 fprintf(fplog, "\n");
1244 prepare_to_do_exchange(FILE *fplog,
1245 const int *destinations,
1246 const int replica_id,
1252 gmx_bool *bThisReplicaExchanged)
1255 /* Hold the cyclic decomposition of the (multiple) replica
1257 gmx_bool bAnyReplicaExchanged = FALSE;
1258 *bThisReplicaExchanged = FALSE;
1260 for (i = 0; i < nrepl; i++)
1262 if (destinations[i] != i)
1264 /* only mark as exchanged if the index has been shuffled */
1265 bAnyReplicaExchanged = TRUE;
1269 if (bAnyReplicaExchanged)
1271 /* reinitialize the placeholder arrays */
1272 for (i = 0; i < nrepl; i++)
1274 for (j = 0; j < nrepl; j++)
1281 /* Identify the cyclic decomposition of the permutation (very
1282 * fast if neighbor replica exchange). */
1283 cyclic_decomposition(destinations, cyclic, incycle, nrepl, maxswap);
1285 /* Now translate the decomposition into a replica exchange
1286 * order at each step. */
1287 compute_exchange_order(fplog, cyclic, order, nrepl, *maxswap);
1289 /* Did this replica do any exchange at any point? */
1290 for (j = 0; j < *maxswap; j++)
1292 if (replica_id != order[replica_id][j])
1294 *bThisReplicaExchanged = TRUE;
1301 gmx_bool replica_exchange(FILE *fplog, const t_commrec *cr, struct gmx_repl_ex *re,
1302 t_state *state, gmx_enerdata_t *enerd,
1303 t_state *state_local, gmx_int64_t step, real time)
1307 int exchange_partner;
1309 /* Number of rounds of exchanges needed to deal with any multiple
1311 /* Where each replica ends up after the exchange attempt(s). */
1312 /* The order in which multiple exchanges will occur. */
1313 gmx_bool bThisReplicaExchanged = FALSE;
1317 replica_id = re->repl;
1318 test_for_replica_exchange(fplog, cr->ms, re, enerd, det(state_local->box), step, time);
1319 prepare_to_do_exchange(fplog, re->destinations, replica_id, re->nrepl, &maxswap,
1320 re->order, re->cyclic, re->incycle, &bThisReplicaExchanged);
1322 /* Do intra-simulation broadcast so all processors belonging to
1323 * each simulation know whether they need to participate in
1324 * collecting the state. Otherwise, they might as well get on with
1325 * the next thing to do. */
1326 if (DOMAINDECOMP(cr))
1329 MPI_Bcast(&bThisReplicaExchanged, sizeof(gmx_bool), MPI_BYTE, MASTERRANK(cr),
1330 cr->mpi_comm_mygroup);
1334 if (bThisReplicaExchanged)
1336 /* Exchange the states */
1337 /* Collect the global state on the master node */
1338 if (DOMAINDECOMP(cr))
1340 dd_collect_state(cr->dd, state_local, state);
1344 copy_state_nonatomdata(state_local, state);
1349 /* There will be only one swap cycle with standard replica
1350 * exchange, but there may be multiple swap cycles if we
1351 * allow multiple swaps. */
1353 for (j = 0; j < maxswap; j++)
1355 exchange_partner = re->order[replica_id][j];
1357 if (exchange_partner != replica_id)
1359 /* Exchange the global states between the master nodes */
1362 fprintf(debug, "Exchanging %d with %d\n", replica_id, exchange_partner);
1364 exchange_state(cr->ms, exchange_partner, state);
1367 /* For temperature-type replica exchange, we need to scale
1368 * the velocities. */
1369 if (re->type == ereTEMP || re->type == ereTL)
1371 scale_velocities(state, sqrt(re->q[ereTEMP][replica_id]/re->q[ereTEMP][re->destinations[replica_id]]));
1376 /* With domain decomposition the global state is distributed later */
1377 if (!DOMAINDECOMP(cr))
1379 /* Copy the global state to the local state data structure */
1380 copy_state_nonatomdata(state, state_local);
1384 return bThisReplicaExchanged;
1387 void print_replica_exchange_statistics(FILE *fplog, struct gmx_repl_ex *re)
1391 fprintf(fplog, "\nReplica exchange statistics\n");
1395 fprintf(fplog, "Repl %d attempts, %d odd, %d even\n",
1396 re->nattempt[0]+re->nattempt[1], re->nattempt[1], re->nattempt[0]);
1398 fprintf(fplog, "Repl average probabilities:\n");
1399 for (i = 1; i < re->nrepl; i++)
1401 if (re->nattempt[i%2] == 0)
1407 re->prob[i] = re->prob_sum[i]/re->nattempt[i%2];
1410 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1411 print_prob(fplog, "", re->nrepl, re->prob);
1413 fprintf(fplog, "Repl number of exchanges:\n");
1414 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1415 print_count(fplog, "", re->nrepl, re->nexchange);
1417 fprintf(fplog, "Repl average number of exchanges:\n");
1418 for (i = 1; i < re->nrepl; i++)
1420 if (re->nattempt[i%2] == 0)
1426 re->prob[i] = ((real)re->nexchange[i])/re->nattempt[i%2];
1429 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1430 print_prob(fplog, "", re->nrepl, re->prob);
1432 fprintf(fplog, "\n");
1434 /* print the transition matrix */
1435 print_transition_matrix(fplog, re->nrepl, re->nmoves, re->nattempt);