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 * 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.
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.
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.
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.
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.
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.
56 #define PROBABILITYCUTOFF 100
57 /* we don't bother evaluating if events are more rare than exp(-100) = 3.7x10^-44 */
60 ereTEMP, ereLAMBDA, ereENDSINGLE, ereTL, ereNR
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 */
68 typedef struct gmx_repl_ex
87 /* these are helper arrays for replica exchange; allocated here so they
88 don't have to be allocated each time */
96 /* helper arrays to hold the quantities that are exchanged */
105 static gmx_bool repl_quantity(FILE *fplog, const gmx_multisim_t *ms,
106 struct gmx_repl_ex *re, int ere, real q)
112 snew(qall, ms->nsim);
114 gmx_sum_sim(ms->nsim, qall, ms);
117 for (s = 1; s < ms->nsim; s++)
119 if (qall[s] != qall[0])
127 /* Set the replica exchange type and quantities */
130 snew(re->q[ere], re->nrepl);
131 for (s = 0; s < ms->nsim; s++)
133 re->q[ere][s] = qall[s];
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)
148 struct gmx_repl_ex *re;
150 gmx_bool bLambda = FALSE;
152 fprintf(fplog, "\nInitializing Replica Exchange\n");
154 if (ms == NULL || ms->nsim == 1)
156 gmx_fatal(FARGS, "Nothing to exchange with only one replica, maybe you forgot to set the -multi option of mdrun?");
162 re->nrepl = ms->nsim;
163 snew(re->q, ereENDSINGLE);
165 fprintf(fplog, "Repl There are %d replicas:\n", re->nrepl);
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);
179 re->temp = ir->opts.ref_t[0];
180 for (i = 1; (i < ir->opts.ngtc); i++)
182 if (ir->opts.ref_t[i] != re->temp)
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");
190 bTemp = repl_quantity(fplog, ms, re, ereTEMP, re->temp);
191 if (ir->efep != efepNO)
193 bLambda = repl_quantity(fplog, ms, re, ereLAMBDA, (real)ir->fepvals->init_fep_state);
195 if (re->type == -1) /* nothing was assigned */
197 gmx_fatal(FARGS, "The properties of the %d systems are all the same, there is nothing to exchange", re->nrepl);
199 if (bLambda && bTemp)
206 please_cite(fplog, "Sugita1999a");
207 if (ir->epc != epcNO)
210 fprintf(fplog, "Repl Using Constant Pressure REMD.\n");
211 please_cite(fplog, "Okabe2001a");
213 if (ir->etc == etcBERENDSEN)
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));
221 if (ir->fepvals->delta_lambda != 0) /* check this? */
223 gmx_fatal(FARGS, "delta_lambda is not zero");
228 snew(re->pres, re->nrepl);
229 if (ir->epct == epctSURFACETENSION)
231 pres = ir->ref_p[ZZ][ZZ];
237 for (i = 0; i < DIM; i++)
239 if (ir->compress[i][i] != 0)
241 pres += ir->ref_p[i][i];
247 re->pres[re->repl] = pres;
248 gmx_sum_sim(re->nrepl, re->pres, ms);
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++)
260 if (re->type < ereENDSINGLE)
263 for (i = 0; i < re->nrepl; i++)
265 for (j = i+1; j < re->nrepl; j++)
267 if (re->q[re->type][re->ind[j]] < re->q[re->type][re->ind[i]])
270 re->ind[i] = re->ind[j];
273 else if (re->q[re->type][re->ind[j]] == re->q[re->type][re->ind[i]])
275 gmx_fatal(FARGS, "Two replicas have identical %ss", erename[re->type]);
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++)
285 re->allswaps[i] = re->ind[i];
291 fprintf(fplog, "\nReplica exchange in temperature\n");
292 for (i = 0; i < re->nrepl; i++)
294 fprintf(fplog, " %5.1f", re->q[re->type][re->ind[i]]);
296 fprintf(fplog, "\n");
299 fprintf(fplog, "\nReplica exchange in lambda\n");
300 for (i = 0; i < re->nrepl; i++)
302 fprintf(fplog, " %3d", (int)re->q[re->type][re->ind[i]]);
304 fprintf(fplog, "\n");
307 fprintf(fplog, "\nReplica exchange in temperature and lambda state\n");
308 for (i = 0; i < re->nrepl; i++)
310 fprintf(fplog, " %5.1f", re->q[ereTEMP][re->ind[i]]);
312 fprintf(fplog, "\n");
313 for (i = 0; i < re->nrepl; i++)
315 fprintf(fplog, " %5d", (int)re->q[ereLAMBDA][re->ind[i]]);
317 fprintf(fplog, "\n");
320 gmx_incons("Unknown replica exchange quantity");
324 fprintf(fplog, "\nRepl p");
325 for (i = 0; i < re->nrepl; i++)
327 fprintf(fplog, " %5.2f", re->pres[re->ind[i]]);
330 for (i = 0; i < re->nrepl; i++)
332 if ((i > 0) && (re->pres[re->ind[i]] < re->pres[re->ind[i-1]]))
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");
344 re->seed = make_seed();
350 gmx_sumi_sim(1, &(re->seed), ms);
354 re->seed = init_seed;
356 fprintf(fplog, "\nReplica exchange interval: %d\n", re->nst);
357 fprintf(fplog, "\nReplica random seed: %d\n", re->seed);
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++)
367 snew(re->nmoves[i], re->nrepl);
369 fprintf(fplog, "Replica exchange information below: x=exchange, pr=probability\n");
371 /* generate space for the helper functions so we don't have to snew each time */
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++)
380 snew(re->cyclic[i], re->nrepl);
381 snew(re->order[i], re->nrepl);
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++)
394 snew(re->de[i], re->nrepl);
400 static void exchange_reals(const gmx_multisim_t *ms, int b, real *v, int n)
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);
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);
424 for (i = 0; i < n; i++)
433 static void exchange_ints(const gmx_multisim_t *ms, int b, int *v, int n)
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);
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);
457 for (i = 0; i < n; i++)
465 static void exchange_doubles(const gmx_multisim_t *ms, int b, double *v, int n)
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);
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);
489 for (i = 0; i < n; i++)
497 static void exchange_rvecs(const gmx_multisim_t *ms, int b, rvec *v, int n)
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);
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);
521 for (i = 0; i < n; i++)
523 copy_rvec(buf[i], v[i]);
529 static void exchange_state(const gmx_multisim_t *ms, int b, t_state *state)
531 /* When t_state changes, this code should be updated. */
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);
553 static void copy_rvecs(rvec *s, rvec *d, int n)
559 for (i = 0; i < n; i++)
561 copy_rvec(s[i], d[i]);
566 static void copy_doubles(const double *s, double *d, int n)
572 for (i = 0; i < n; i++)
579 static void copy_reals(const real *s, real *d, int n)
585 for (i = 0; i < n; i++)
592 static void copy_ints(const int *s, int *d, int n)
598 for (i = 0; i < n; i++)
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);
610 static void copy_state_nonatomdata(t_state *state, t_state *state_local)
612 /* When t_state changes, this code should be updated. */
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);
636 static void scale_velocities(t_state *state, real fac)
642 for (i = 0; i < state->natoms; i++)
644 svmul(fac, state->v[i], state->v[i]);
649 static void pd_collect_state(const t_commrec *cr, t_state *state)
655 fprintf(debug, "Collecting state before exchange\n");
657 shift = cr->nnodes - cr->npmenodes - 1;
658 move_rvecs(cr, FALSE, FALSE, GMX_LEFT, GMX_RIGHT, state->x, NULL, shift, NULL);
661 move_rvecs(cr, FALSE, FALSE, GMX_LEFT, GMX_RIGHT, state->v, NULL, shift, NULL);
665 move_rvecs(cr, FALSE, FALSE, GMX_LEFT, GMX_RIGHT, state->sd_X, NULL, shift, NULL);
669 static void print_transition_matrix(FILE *fplog, const char *leg, int n, int **nmoves, int *nattempt)
674 ntot = nattempt[0] + nattempt[1];
675 fprintf(fplog, "\n");
676 fprintf(fplog, "Repl");
677 for (i = 0; i < n; i++)
679 fprintf(fplog, " "); /* put the title closer to the center */
681 fprintf(fplog, "Empirical Transition Matrix\n");
683 fprintf(fplog, "Repl");
684 for (i = 0; i < n; i++)
686 fprintf(fplog, "%8d", (i+1));
688 fprintf(fplog, "\n");
690 for (i = 0; i < n; i++)
692 fprintf(fplog, "Repl");
693 for (j = 0; j < n; j++)
696 if (nmoves[i][j] > 0)
698 Tprint = nmoves[i][j]/(2.0*ntot);
700 fprintf(fplog, "%8.4f", Tprint);
702 fprintf(fplog, "%3d\n", i);
706 static void print_ind(FILE *fplog, const char *leg, int n, int *ind, gmx_bool *bEx)
710 fprintf(fplog, "Repl %2s %2d", leg, ind[0]);
711 for (i = 1; i < n; i++)
713 fprintf(fplog, " %c %2d", (bEx != 0 && bEx[i]) ? 'x' : ' ', ind[i]);
715 fprintf(fplog, "\n");
718 static void print_allswitchind(FILE *fplog, int n, int *ind, int *pind, int *allswaps, int *tmpswap)
722 for (i = 0; i < n; i++)
724 tmpswap[i] = allswaps[i];
726 for (i = 0; i < n; i++)
728 allswaps[i] = tmpswap[pind[i]];
731 fprintf(fplog, "\nAccepted Exchanges: ");
732 for (i = 0; i < n; i++)
734 fprintf(fplog, "%d ", pind[i]);
736 fprintf(fplog, "\n");
738 /* the "Order After Exchange" is the state label corresponding to the configuration that
739 started in state listed in order, i.e.
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
749 fprintf(fplog, "Order After Exchange: ");
750 for (i = 0; i < n; i++)
752 fprintf(fplog, "%d ", allswaps[i]);
754 fprintf(fplog, "\n\n");
757 static void print_prob(FILE *fplog, const char *leg, int n, real *prob)
762 fprintf(fplog, "Repl %2s ", leg);
763 for (i = 1; i < n; i++)
767 sprintf(buf, "%4.2f", prob[i]);
768 fprintf(fplog, " %3s", buf[0] == '1' ? "1.0" : buf+1);
775 fprintf(fplog, "\n");
778 static void print_count(FILE *fplog, const char *leg, int n, int *count)
782 fprintf(fplog, "Repl %2s ", leg);
783 for (i = 1; i < n; i++)
785 fprintf(fplog, " %4d", count[i]);
787 fprintf(fplog, "\n");
790 static real calc_delta(FILE *fplog, gmx_bool bPrint, struct gmx_repl_ex *re, int a, int b, int ap, int bp)
793 real ediff, dpV, delta = 0;
794 real *Epot = re->Epot;
797 real *beta = re->beta;
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 */
806 * Okabe et. al. Chem. Phys. Lett. 335 (2001) 435-439
808 ediff = Epot[b] - Epot[a];
809 delta = -(beta[bp] - beta[ap])*ediff;
812 /* two cases: when we are permuted, and not. */
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] */
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.
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 */
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]);
868 gmx_incons("Unknown replica exchange quantity");
872 fprintf(fplog, "Repl %d <-> %d dE_term = %10.3e (kT)\n", a, b, delta);
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;
880 fprintf(fplog, " dpV = %10.3e d = %10.3e\nb", dpV, delta + dpV);
888 test_for_replica_exchange(FILE *fplog,
889 const gmx_multisim_t *ms,
890 struct gmx_repl_ex *re,
891 gmx_enerdata_t *enerd,
893 gmx_large_int_t step,
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;
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);
911 for (i = 0; i < re->nrepl; i++)
916 re->Vol[re->repl] = vol;
918 if ((re->type == ereTEMP || re->type == ereTL))
920 for (i = 0; i < re->nrepl; i++)
925 re->Epot[re->repl] = enerd->term[F_EPOT];
926 /* temperatures of different states*/
927 for (i = 0; i < re->nrepl; i++)
929 re->beta[i] = 1.0/(re->q[ereTEMP][i]*BOLTZ);
934 for (i = 0; i < re->nrepl; i++)
936 re->beta[i] = 1.0/(re->temp*BOLTZ); /* we have a single temperature */
939 if (re->type == ereLAMBDA || re->type == ereTL)
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++)
947 for (j = 0; j < re->nrepl; j++)
952 for (i = 0; i < re->nrepl; i++)
954 re->de[i][re->repl] = (enerd->enerpart_lambda[(int)re->q[ereLAMBDA][i]+1]-enerd->enerpart_lambda[0]);
958 /* now actually do the communication */
961 gmx_sum_sim(re->nrepl, re->Vol, ms);
965 gmx_sum_sim(re->nrepl, re->Epot, ms);
969 for (i = 0; i < re->nrepl; i++)
971 gmx_sum_sim(re->nrepl, re->de[i], ms);
975 /* make a duplicate set of indices for shuffling */
976 for (i = 0; i < re->nrepl; i++)
978 pind[i] = re->ind[i];
983 /* multiple random switch exchange */
984 for (i = 0; i < re->nex; i++)
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)));
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 bEx[0] = (rando(&(re->seed)) < prob[0]);
1033 re->prob_sum[0] += prob[0];
1037 /* swap the states */
1039 pind[i0] = pind[i1];
1043 re->nattempt[0]++; /* keep track of total permutation trials here */
1044 print_allswitchind(fplog, re->nrepl, re->ind, pind, re->allswaps, re->tmpswap);
1048 /* standard nearest neighbor replica exchange */
1049 m = (step / re->nst) % 2;
1050 for (i = 1; i < re->nrepl; i++)
1055 bPrint = (re->repl == a || re->repl == b);
1058 delta = calc_delta(fplog, bPrint, re, a, b, a, b);
1067 if (delta > PROBABILITYCUTOFF)
1073 prob[i] = exp(-delta);
1075 /* roll a number to determine if accepted */
1076 bEx[i] = (rando(&(re->seed)) < prob[i]);
1078 re->prob_sum[i] += prob[i];
1082 /* swap these two */
1084 pind[i-1] = pind[i];
1086 re->nexchange[i]++; /* statistics for back compatibility */
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");
1102 /* record which moves were made and accepted */
1103 for (i = 0; i < re->nrepl; i++)
1105 re->nmoves[re->ind[i]][pind[i]] += 1;
1106 re->nmoves[pind[i]][re->ind[i]] += 1;
1108 fflush(fplog); /* make sure we can see what the last exchange was */
1111 static void write_debug_x(t_state *state)
1117 for (i = 0; i < state->natoms; i += 10)
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]);
1125 cyclic_decomposition(FILE *fplog,
1126 const int *destinations,
1135 for (i = 0; i < nrepl; i++)
1139 for (i = 0; i < nrepl; i++) /* one cycle for each replica */
1150 for (j = 0; j < nrepl; j++) /* potentially all cycles are part, but we will break first */
1152 p = destinations[p]; /* start permuting */
1160 break; /* we've reached the original element, the cycle is complete, and we marked the end. */
1164 cyclic[i][c] = p; /* each permutation gives a new member of the cycle */
1170 *nswap = maxlen - 1;
1174 for (i = 0; i < nrepl; i++)
1176 fprintf(debug, "Cycle %d:", i);
1177 for (j = 0; j < nrepl; j++)
1179 if (cyclic[i][j] < 0)
1183 fprintf(debug, "%2d", cyclic[i][j]);
1185 fprintf(debug, "\n");
1192 compute_exchange_order(FILE *fplog,
1200 for (j = 0; j < maxswap; j++)
1202 for (i = 0; i < nrepl; i++)
1204 if (cyclic[i][j+1] >= 0)
1206 order[cyclic[i][j+1]][j] = cyclic[i][j];
1207 order[cyclic[i][j]][j] = cyclic[i][j+1];
1210 for (i = 0; i < nrepl; i++)
1212 if (order[i][j] < 0)
1214 order[i][j] = i; /* if it's not exchanging, it should stay this round*/
1221 fprintf(fplog, "Replica Exchange Order\n");
1222 for (i = 0; i < nrepl; i++)
1224 fprintf(fplog, "Replica %d:", i);
1225 for (j = 0; j < maxswap; j++)
1227 if (order[i][j] < 0)
1231 fprintf(debug, "%2d", order[i][j]);
1233 fprintf(fplog, "\n");
1240 prepare_to_do_exchange(FILE *fplog,
1241 const int *destinations,
1242 const int replica_id,
1248 gmx_bool *bThisReplicaExchanged)
1251 /* Hold the cyclic decomposition of the (multiple) replica
1253 gmx_bool bAnyReplicaExchanged = FALSE;
1254 *bThisReplicaExchanged = FALSE;
1256 for (i = 0; i < nrepl; i++)
1258 if (destinations[i] != i)
1260 /* only mark as exchanged if the index has been shuffled */
1261 bAnyReplicaExchanged = TRUE;
1265 if (bAnyReplicaExchanged)
1267 /* reinitialize the placeholder arrays */
1268 for (i = 0; i < nrepl; i++)
1270 for (j = 0; j < nrepl; j++)
1277 /* Identify the cyclic decomposition of the permutation (very
1278 * fast if neighbor replica exchange). */
1279 cyclic_decomposition(fplog, destinations, cyclic, incycle, nrepl, maxswap);
1281 /* Now translate the decomposition into a replica exchange
1282 * order at each step. */
1283 compute_exchange_order(fplog, cyclic, order, nrepl, *maxswap);
1285 /* Did this replica do any exchange at any point? */
1286 for (j = 0; j < *maxswap; j++)
1288 if (replica_id != order[replica_id][j])
1290 *bThisReplicaExchanged = TRUE;
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)
1303 int exchange_partner;
1305 /* Number of rounds of exchanges needed to deal with any multiple
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;
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);
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. */
1325 MPI_Bcast(&bThisReplicaExchanged, sizeof(gmx_bool), MPI_BYTE, MASTERRANK(cr),
1326 cr->mpi_comm_mygroup);
1330 if (bThisReplicaExchanged)
1332 /* Exchange the states */
1336 /* Collect the global state on the master node */
1337 if (DOMAINDECOMP(cr))
1339 dd_collect_state(cr->dd, state_local, state);
1343 pd_collect_state(cr, 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 bcast_state(cr, state, FALSE);
1389 return bThisReplicaExchanged;
1392 void print_replica_exchange_statistics(FILE *fplog, struct gmx_repl_ex *re)
1396 fprintf(fplog, "\nReplica exchange statistics\n");
1400 fprintf(fplog, "Repl %d attempts, %d odd, %d even\n",
1401 re->nattempt[0]+re->nattempt[1], re->nattempt[1], re->nattempt[0]);
1403 fprintf(fplog, "Repl average probabilities:\n");
1404 for (i = 1; i < re->nrepl; i++)
1406 if (re->nattempt[i%2] == 0)
1412 re->prob[i] = re->prob_sum[i]/re->nattempt[i%2];
1415 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1416 print_prob(fplog, "", re->nrepl, re->prob);
1418 fprintf(fplog, "Repl number of exchanges:\n");
1419 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1420 print_count(fplog, "", re->nrepl, re->nexchange);
1422 fprintf(fplog, "Repl average number of exchanges:\n");
1423 for (i = 1; i < re->nrepl; i++)
1425 if (re->nattempt[i%2] == 0)
1431 re->prob[i] = ((real)re->nexchange[i])/re->nattempt[i%2];
1434 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1435 print_prob(fplog, "", re->nrepl, re->prob);
1437 fprintf(fplog, "\n");
1439 /* print the transition matrix */
1440 print_transition_matrix(fplog, "", re->nrepl, re->nmoves, re->nattempt);