Broke up copyright.*
[gromacs/AngularHB.git] / src / programs / mdrun / repl_ex.cpp
blob05867f110ea0fc3433746b8648b8b9a4f8ef81d9
1 /*
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,2015, 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.
38 #include "gmxpre.h"
40 #include "repl_ex.h"
42 #include "config.h"
44 #include <math.h>
46 #include "gromacs/domdec/domdec.h"
47 #include "gromacs/gmxlib/main.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/math/units.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/random/random.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/pleasecite.h"
56 #include "gromacs/utility/smalloc.h"
58 #define PROBABILITYCUTOFF 100
59 /* we don't bother evaluating if events are more rare than exp(-100) = 3.7x10^-44 */
61 //! Rank in the multisimulaiton
62 #define MSRANK(ms, nodeid) (nodeid)
64 enum {
65 ereTEMP, ereLAMBDA, ereENDSINGLE, ereTL, ereNR
67 const char *erename[ereNR] = { "temperature", "lambda", "end_single_marker", "temperature and lambda"};
68 /* end_single_marker merely notes the end of single variable replica exchange. All types higher than
69 it are multiple replica exchange methods */
70 /* Eventually, should add 'pressure', 'temperature and pressure', 'lambda_and_pressure', 'temperature_lambda_pressure'?;
71 Let's wait until we feel better about the pressure control methods giving exact ensembles. Right now, we assume constant pressure */
73 typedef struct gmx_repl_ex
75 int repl; /* replica ID */
76 int nrepl; /* total number of replica */
77 real temp; /* temperature */
78 int type; /* replica exchange type from ere enum */
79 real **q; /* quantity, e.g. temperature or lambda; first index is ere, second index is replica ID */
80 gmx_bool bNPT; /* use constant pressure and temperature */
81 real *pres; /* replica pressures */
82 int *ind; /* replica indices */
83 int *allswaps; /* used for keeping track of all the replica swaps */
84 int nst; /* replica exchange interval (number of steps) */
85 int nex; /* number of exchanges per interval */
86 int seed; /* random seed */
87 int nattempt[2]; /* number of even and odd replica change attempts */
88 real *prob_sum; /* sum of probabilities */
89 int **nmoves; /* number of moves between replicas i and j */
90 int *nexchange; /* i-th element of the array is the number of exchanges between replica i-1 and i */
91 gmx_rng_t rng; /* random number generator */
93 /* these are helper arrays for replica exchange; allocated here so they
94 don't have to be allocated each time */
95 int *destinations;
96 int **cyclic;
97 int **order;
98 int *tmpswap;
99 gmx_bool *incycle;
100 gmx_bool *bEx;
102 /* helper arrays to hold the quantities that are exchanged */
103 real *prob;
104 real *Epot;
105 real *beta;
106 real *Vol;
107 real **de;
109 } t_gmx_repl_ex;
111 static gmx_bool repl_quantity(const gmx_multisim_t *ms,
112 struct gmx_repl_ex *re, int ere, real q)
114 real *qall;
115 gmx_bool bDiff;
116 int s;
118 snew(qall, ms->nsim);
119 qall[re->repl] = q;
120 gmx_sum_sim(ms->nsim, qall, ms);
122 bDiff = FALSE;
123 for (s = 1; s < ms->nsim; s++)
125 if (qall[s] != qall[0])
127 bDiff = TRUE;
131 if (bDiff)
133 /* Set the replica exchange type and quantities */
134 re->type = ere;
136 snew(re->q[ere], re->nrepl);
137 for (s = 0; s < ms->nsim; s++)
139 re->q[ere][s] = qall[s];
142 sfree(qall);
143 return bDiff;
146 gmx_repl_ex_t init_replica_exchange(FILE *fplog,
147 const gmx_multisim_t *ms,
148 const t_state *state,
149 const t_inputrec *ir,
150 int nst, int nex, int init_seed)
152 real pres;
153 int i, j, k;
154 struct gmx_repl_ex *re;
155 gmx_bool bTemp;
156 gmx_bool bLambda = FALSE;
158 fprintf(fplog, "\nInitializing Replica Exchange\n");
160 if (ms == NULL || ms->nsim == 1)
162 gmx_fatal(FARGS, "Nothing to exchange with only one replica, maybe you forgot to set the -multi option of mdrun?");
164 if (!EI_DYNAMICS(ir->eI))
166 gmx_fatal(FARGS, "Replica exchange is only supported by dynamical simulations");
167 /* Note that PAR(cr) is defined by cr->nnodes > 1, which is
168 * distinct from MULTISIM(cr). A multi-simulation only runs
169 * with real MPI parallelism, but this does not imply PAR(cr)
170 * is true!
172 * Since we are using a dynamical integrator, the only
173 * decomposition is DD, so PAR(cr) and DOMAINDECOMP(cr) are
174 * synonymous. The only way for cr->nnodes > 1 to be true is
175 * if we are using DD. */
178 snew(re, 1);
180 re->repl = ms->sim;
181 re->nrepl = ms->nsim;
182 snew(re->q, ereENDSINGLE);
184 fprintf(fplog, "Repl There are %d replicas:\n", re->nrepl);
186 check_multi_int(fplog, ms, state->natoms, "the number of atoms", FALSE);
187 check_multi_int(fplog, ms, ir->eI, "the integrator", FALSE);
188 check_multi_int64(fplog, ms, ir->init_step+ir->nsteps, "init_step+nsteps", FALSE);
189 check_multi_int64(fplog, ms, (ir->init_step+nst-1)/nst,
190 "first exchange step: init_step/-replex", FALSE);
191 check_multi_int(fplog, ms, ir->etc, "the temperature coupling", FALSE);
192 check_multi_int(fplog, ms, ir->opts.ngtc,
193 "the number of temperature coupling groups", FALSE);
194 check_multi_int(fplog, ms, ir->epc, "the pressure coupling", FALSE);
195 check_multi_int(fplog, ms, ir->efep, "free energy", FALSE);
196 check_multi_int(fplog, ms, ir->fepvals->n_lambda, "number of lambda states", FALSE);
198 re->temp = ir->opts.ref_t[0];
199 for (i = 1; (i < ir->opts.ngtc); i++)
201 if (ir->opts.ref_t[i] != re->temp)
203 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
204 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
208 re->type = -1;
209 bTemp = repl_quantity(ms, re, ereTEMP, re->temp);
210 if (ir->efep != efepNO)
212 bLambda = repl_quantity(ms, re, ereLAMBDA, (real)ir->fepvals->init_fep_state);
214 if (re->type == -1) /* nothing was assigned */
216 gmx_fatal(FARGS, "The properties of the %d systems are all the same, there is nothing to exchange", re->nrepl);
218 if (bLambda && bTemp)
220 re->type = ereTL;
223 if (bTemp)
225 please_cite(fplog, "Sugita1999a");
226 if (ir->epc != epcNO)
228 re->bNPT = TRUE;
229 fprintf(fplog, "Repl Using Constant Pressure REMD.\n");
230 please_cite(fplog, "Okabe2001a");
232 if (ir->etc == etcBERENDSEN)
234 gmx_fatal(FARGS, "REMD with the %s thermostat does not produce correct potential energy distributions, consider using the %s thermostat instead",
235 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
238 if (bLambda)
240 if (ir->fepvals->delta_lambda != 0) /* check this? */
242 gmx_fatal(FARGS, "delta_lambda is not zero");
245 if (re->bNPT)
247 snew(re->pres, re->nrepl);
248 if (ir->epct == epctSURFACETENSION)
250 pres = ir->ref_p[ZZ][ZZ];
252 else
254 pres = 0;
255 j = 0;
256 for (i = 0; i < DIM; i++)
258 if (ir->compress[i][i] != 0)
260 pres += ir->ref_p[i][i];
261 j++;
264 pres /= j;
266 re->pres[re->repl] = pres;
267 gmx_sum_sim(re->nrepl, re->pres, ms);
270 /* Make an index for increasing replica order */
271 /* only makes sense if one or the other is varying, not both!
272 if both are varying, we trust the order the person gave. */
273 snew(re->ind, re->nrepl);
274 for (i = 0; i < re->nrepl; i++)
276 re->ind[i] = i;
279 if (re->type < ereENDSINGLE)
282 for (i = 0; i < re->nrepl; i++)
284 for (j = i+1; j < re->nrepl; j++)
286 if (re->q[re->type][re->ind[j]] < re->q[re->type][re->ind[i]])
288 /* Unordered replicas are supposed to work, but there
289 * is still an issues somewhere.
290 * Note that at this point still re->ind[i]=i.
292 gmx_fatal(FARGS, "Replicas with indices %d < %d have %ss %g > %g, please order your replicas on increasing %s",
293 i, j,
294 erename[re->type],
295 re->q[re->type][i], re->q[re->type][j],
296 erename[re->type]);
298 k = re->ind[i];
299 re->ind[i] = re->ind[j];
300 re->ind[j] = k;
302 else if (re->q[re->type][re->ind[j]] == re->q[re->type][re->ind[i]])
304 gmx_fatal(FARGS, "Two replicas have identical %ss", erename[re->type]);
310 /* keep track of all the swaps, starting with the initial placement. */
311 snew(re->allswaps, re->nrepl);
312 for (i = 0; i < re->nrepl; i++)
314 re->allswaps[i] = re->ind[i];
317 switch (re->type)
319 case ereTEMP:
320 fprintf(fplog, "\nReplica exchange in temperature\n");
321 for (i = 0; i < re->nrepl; i++)
323 fprintf(fplog, " %5.1f", re->q[re->type][re->ind[i]]);
325 fprintf(fplog, "\n");
326 break;
327 case ereLAMBDA:
328 fprintf(fplog, "\nReplica exchange in lambda\n");
329 for (i = 0; i < re->nrepl; i++)
331 fprintf(fplog, " %3d", (int)re->q[re->type][re->ind[i]]);
333 fprintf(fplog, "\n");
334 break;
335 case ereTL:
336 fprintf(fplog, "\nReplica exchange in temperature and lambda state\n");
337 for (i = 0; i < re->nrepl; i++)
339 fprintf(fplog, " %5.1f", re->q[ereTEMP][re->ind[i]]);
341 fprintf(fplog, "\n");
342 for (i = 0; i < re->nrepl; i++)
344 fprintf(fplog, " %5d", (int)re->q[ereLAMBDA][re->ind[i]]);
346 fprintf(fplog, "\n");
347 break;
348 default:
349 gmx_incons("Unknown replica exchange quantity");
351 if (re->bNPT)
353 fprintf(fplog, "\nRepl p");
354 for (i = 0; i < re->nrepl; i++)
356 fprintf(fplog, " %5.2f", re->pres[re->ind[i]]);
359 for (i = 0; i < re->nrepl; i++)
361 if ((i > 0) && (re->pres[re->ind[i]] < re->pres[re->ind[i-1]]))
363 fprintf(fplog, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
364 fprintf(stderr, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
368 re->nst = nst;
369 if (init_seed == -1)
371 if (MASTERSIM(ms))
373 re->seed = (int)gmx_rng_make_seed();
375 else
377 re->seed = 0;
379 gmx_sumi_sim(1, &(re->seed), ms);
381 else
383 re->seed = init_seed;
385 fprintf(fplog, "\nReplica exchange interval: %d\n", re->nst);
386 fprintf(fplog, "\nReplica random seed: %d\n", re->seed);
387 re->rng = gmx_rng_init(re->seed);
389 re->nattempt[0] = 0;
390 re->nattempt[1] = 0;
392 snew(re->prob_sum, re->nrepl);
393 snew(re->nexchange, re->nrepl);
394 snew(re->nmoves, re->nrepl);
395 for (i = 0; i < re->nrepl; i++)
397 snew(re->nmoves[i], re->nrepl);
399 fprintf(fplog, "Replica exchange information below: ex and x = exchange, pr = probability\n");
401 /* generate space for the helper functions so we don't have to snew each time */
403 snew(re->destinations, re->nrepl);
404 snew(re->incycle, re->nrepl);
405 snew(re->tmpswap, re->nrepl);
406 snew(re->cyclic, re->nrepl);
407 snew(re->order, re->nrepl);
408 for (i = 0; i < re->nrepl; i++)
410 snew(re->cyclic[i], re->nrepl+1);
411 snew(re->order[i], re->nrepl);
413 /* allocate space for the functions storing the data for the replicas */
414 /* not all of these arrays needed in all cases, but they don't take
415 up much space, since the max size is nrepl**2 */
416 snew(re->prob, re->nrepl);
417 snew(re->bEx, re->nrepl);
418 snew(re->beta, re->nrepl);
419 snew(re->Vol, re->nrepl);
420 snew(re->Epot, re->nrepl);
421 snew(re->de, re->nrepl);
422 for (i = 0; i < re->nrepl; i++)
424 snew(re->de[i], re->nrepl);
426 re->nex = nex;
427 return re;
430 static void exchange_reals(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, real *v, int n)
432 real *buf;
433 int i;
435 if (v)
437 snew(buf, n);
438 #ifdef GMX_MPI
440 MPI_Sendrecv(v, n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
441 buf,n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
442 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
445 MPI_Request mpi_req;
447 MPI_Isend(v, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
448 ms->mpi_comm_masters, &mpi_req);
449 MPI_Recv(buf, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
450 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
451 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
453 #endif
454 for (i = 0; i < n; i++)
456 v[i] = buf[i];
458 sfree(buf);
463 static void exchange_doubles(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, double *v, int n)
465 double *buf;
466 int i;
468 if (v)
470 snew(buf, n);
471 #ifdef GMX_MPI
473 MPI_Sendrecv(v, n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
474 buf,n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
475 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
478 MPI_Request mpi_req;
480 MPI_Isend(v, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
481 ms->mpi_comm_masters, &mpi_req);
482 MPI_Recv(buf, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
483 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
484 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
486 #endif
487 for (i = 0; i < n; i++)
489 v[i] = buf[i];
491 sfree(buf);
495 static void exchange_rvecs(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, rvec *v, int n)
497 rvec *buf;
498 int i;
500 if (v)
502 snew(buf, n);
503 #ifdef GMX_MPI
505 MPI_Sendrecv(v[0], n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
506 buf[0],n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
507 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
510 MPI_Request mpi_req;
512 MPI_Isend(v[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
513 ms->mpi_comm_masters, &mpi_req);
514 MPI_Recv(buf[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
515 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
516 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
518 #endif
519 for (i = 0; i < n; i++)
521 copy_rvec(buf[i], v[i]);
523 sfree(buf);
527 static void exchange_state(const gmx_multisim_t *ms, int b, t_state *state)
529 /* When t_state changes, this code should be updated. */
530 int ngtc, nnhpres;
531 ngtc = state->ngtc * state->nhchainlength;
532 nnhpres = state->nnhpres* state->nhchainlength;
533 exchange_rvecs(ms, b, state->box, DIM);
534 exchange_rvecs(ms, b, state->box_rel, DIM);
535 exchange_rvecs(ms, b, state->boxv, DIM);
536 exchange_reals(ms, b, &(state->veta), 1);
537 exchange_reals(ms, b, &(state->vol0), 1);
538 exchange_rvecs(ms, b, state->svir_prev, DIM);
539 exchange_rvecs(ms, b, state->fvir_prev, DIM);
540 exchange_rvecs(ms, b, state->pres_prev, DIM);
541 exchange_doubles(ms, b, state->nosehoover_xi, ngtc);
542 exchange_doubles(ms, b, state->nosehoover_vxi, ngtc);
543 exchange_doubles(ms, b, state->nhpres_xi, nnhpres);
544 exchange_doubles(ms, b, state->nhpres_vxi, nnhpres);
545 exchange_doubles(ms, b, state->therm_integral, state->ngtc);
546 exchange_rvecs(ms, b, state->x, state->natoms);
547 exchange_rvecs(ms, b, state->v, state->natoms);
548 exchange_rvecs(ms, b, state->sd_X, state->natoms);
551 static void copy_rvecs(rvec *s, rvec *d, int n)
553 int i;
555 if (d != NULL)
557 for (i = 0; i < n; i++)
559 copy_rvec(s[i], d[i]);
564 static void copy_doubles(const double *s, double *d, int n)
566 int i;
568 if (d != NULL)
570 for (i = 0; i < n; i++)
572 d[i] = s[i];
577 static void copy_reals(const real *s, real *d, int n)
579 int i;
581 if (d != NULL)
583 for (i = 0; i < n; i++)
585 d[i] = s[i];
590 static void copy_ints(const int *s, int *d, int n)
592 int i;
594 if (d != NULL)
596 for (i = 0; i < n; i++)
598 d[i] = s[i];
603 #define scopy_rvecs(v, n) copy_rvecs(state->v, state_local->v, n);
604 #define scopy_doubles(v, n) copy_doubles(state->v, state_local->v, n);
605 #define scopy_reals(v, n) copy_reals(state->v, state_local->v, n);
606 #define scopy_ints(v, n) copy_ints(state->v, state_local->v, n);
608 static void copy_state_nonatomdata(t_state *state, t_state *state_local)
610 /* When t_state changes, this code should be updated. */
611 int ngtc, nnhpres;
612 ngtc = state->ngtc * state->nhchainlength;
613 nnhpres = state->nnhpres* state->nhchainlength;
614 scopy_rvecs(box, DIM);
615 scopy_rvecs(box_rel, DIM);
616 scopy_rvecs(boxv, DIM);
617 state_local->veta = state->veta;
618 state_local->vol0 = state->vol0;
619 scopy_rvecs(svir_prev, DIM);
620 scopy_rvecs(fvir_prev, DIM);
621 scopy_rvecs(pres_prev, DIM);
622 scopy_doubles(nosehoover_xi, ngtc);
623 scopy_doubles(nosehoover_vxi, ngtc);
624 scopy_doubles(nhpres_xi, nnhpres);
625 scopy_doubles(nhpres_vxi, nnhpres);
626 scopy_doubles(therm_integral, state->ngtc);
627 scopy_rvecs(x, state->natoms);
628 scopy_rvecs(v, state->natoms);
629 scopy_rvecs(sd_X, state->natoms);
630 copy_ints(&(state->fep_state), &(state_local->fep_state), 1);
631 scopy_reals(lambda, efptNR);
634 static void scale_velocities(t_state *state, real fac)
636 int i;
638 if (state->v)
640 for (i = 0; i < state->natoms; i++)
642 svmul(fac, state->v[i], state->v[i]);
647 static void print_transition_matrix(FILE *fplog, int n, int **nmoves, int *nattempt)
649 int i, j, ntot;
650 float Tprint;
652 ntot = nattempt[0] + nattempt[1];
653 fprintf(fplog, "\n");
654 fprintf(fplog, "Repl");
655 for (i = 0; i < n; i++)
657 fprintf(fplog, " "); /* put the title closer to the center */
659 fprintf(fplog, "Empirical Transition Matrix\n");
661 fprintf(fplog, "Repl");
662 for (i = 0; i < n; i++)
664 fprintf(fplog, "%8d", (i+1));
666 fprintf(fplog, "\n");
668 for (i = 0; i < n; i++)
670 fprintf(fplog, "Repl");
671 for (j = 0; j < n; j++)
673 Tprint = 0.0;
674 if (nmoves[i][j] > 0)
676 Tprint = nmoves[i][j]/(2.0*ntot);
678 fprintf(fplog, "%8.4f", Tprint);
680 fprintf(fplog, "%3d\n", i);
684 static void print_ind(FILE *fplog, const char *leg, int n, int *ind, gmx_bool *bEx)
686 int i;
688 fprintf(fplog, "Repl %2s %2d", leg, ind[0]);
689 for (i = 1; i < n; i++)
691 fprintf(fplog, " %c %2d", (bEx != 0 && bEx[i]) ? 'x' : ' ', ind[i]);
693 fprintf(fplog, "\n");
696 static void print_allswitchind(FILE *fplog, int n, int *pind, int *allswaps, int *tmpswap)
698 int i;
700 for (i = 0; i < n; i++)
702 tmpswap[i] = allswaps[i];
704 for (i = 0; i < n; i++)
706 allswaps[i] = tmpswap[pind[i]];
709 fprintf(fplog, "\nAccepted Exchanges: ");
710 for (i = 0; i < n; i++)
712 fprintf(fplog, "%d ", pind[i]);
714 fprintf(fplog, "\n");
716 /* the "Order After Exchange" is the state label corresponding to the configuration that
717 started in state listed in order, i.e.
719 3 0 1 2
721 means that the:
722 configuration starting in simulation 3 is now in simulation 0,
723 configuration starting in simulation 0 is now in simulation 1,
724 configuration starting in simulation 1 is now in simulation 2,
725 configuration starting in simulation 2 is now in simulation 3
727 fprintf(fplog, "Order After Exchange: ");
728 for (i = 0; i < n; i++)
730 fprintf(fplog, "%d ", allswaps[i]);
732 fprintf(fplog, "\n\n");
735 static void print_prob(FILE *fplog, const char *leg, int n, real *prob)
737 int i;
738 char buf[8];
740 fprintf(fplog, "Repl %2s ", leg);
741 for (i = 1; i < n; i++)
743 if (prob[i] >= 0)
745 sprintf(buf, "%4.2f", prob[i]);
746 fprintf(fplog, " %3s", buf[0] == '1' ? "1.0" : buf+1);
748 else
750 fprintf(fplog, " ");
753 fprintf(fplog, "\n");
756 static void print_count(FILE *fplog, const char *leg, int n, int *count)
758 int i;
760 fprintf(fplog, "Repl %2s ", leg);
761 for (i = 1; i < n; i++)
763 fprintf(fplog, " %4d", count[i]);
765 fprintf(fplog, "\n");
768 static real calc_delta(FILE *fplog, gmx_bool bPrint, struct gmx_repl_ex *re, int a, int b, int ap, int bp)
771 real ediff, dpV, delta = 0;
772 real *Epot = re->Epot;
773 real *Vol = re->Vol;
774 real **de = re->de;
775 real *beta = re->beta;
777 /* Two cases; we are permuted and not. In all cases, setting ap = a and bp = b will reduce
778 to the non permuted case */
780 switch (re->type)
782 case ereTEMP:
784 * Okabe et. al. Chem. Phys. Lett. 335 (2001) 435-439
786 ediff = Epot[b] - Epot[a];
787 delta = -(beta[bp] - beta[ap])*ediff;
788 break;
789 case ereLAMBDA:
790 /* two cases: when we are permuted, and not. */
791 /* non-permuted:
792 ediff = E_new - E_old
793 = [H_b(x_a) + H_a(x_b)] - [H_b(x_b) + H_a(x_a)]
794 = [H_b(x_a) - H_a(x_a)] + [H_a(x_b) - H_b(x_b)]
795 = de[b][a] + de[a][b] */
797 /* permuted:
798 ediff = E_new - E_old
799 = [H_bp(x_a) + H_ap(x_b)] - [H_bp(x_b) + H_ap(x_a)]
800 = [H_bp(x_a) - H_ap(x_a)] + [H_ap(x_b) - H_bp(x_b)]
801 = [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)]
802 = [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)]
803 = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]) */
804 /* but, in the current code implementation, we flip configurations, not indices . . .
805 So let's examine that.
806 = [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)]
807 = [H_b(x_ap) - H_a(x_ap)] + [H_a(x_bp) - H_b(x_pb)]
808 = (de[b][ap] - de[a][ap]) + (de[a][bp] - de[b][bp]
809 So, if we exchange b<=> bp and a<=> ap, we return to the same result.
810 So the simple solution is to flip the
811 position of perturbed and original indices in the tests.
814 ediff = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]);
815 delta = ediff*beta[a]; /* assume all same temperature in this case */
816 break;
817 case ereTL:
818 /* not permuted: */
819 /* delta = reduced E_new - reduced E_old
820 = [beta_b H_b(x_a) + beta_a H_a(x_b)] - [beta_b H_b(x_b) + beta_a H_a(x_a)]
821 = [beta_b H_b(x_a) - beta_a H_a(x_a)] + [beta_a H_a(x_b) - beta_b H_b(x_b)]
822 = [beta_b dH_b(x_a) + beta_b H_a(x_a) - beta_a H_a(x_a)] +
823 [beta_a dH_a(x_b) + beta_a H_b(x_b) - beta_b H_b(x_b)]
824 = [beta_b dH_b(x_a) + [beta_a dH_a(x_b) +
825 beta_b (H_a(x_a) - H_b(x_b)]) - beta_a (H_a(x_a) - H_b(x_b))
826 = beta_b dH_b(x_a) + beta_a dH_a(x_b) - (beta_b - beta_a)(H_b(x_b) - H_a(x_a) */
827 /* delta = beta[b]*de[b][a] + beta[a]*de[a][b] - (beta[b] - beta[a])*(Epot[b] - Epot[a]; */
828 /* permuted (big breath!) */
829 /* delta = reduced E_new - reduced E_old
830 = [beta_bp H_bp(x_a) + beta_ap H_ap(x_b)] - [beta_bp H_bp(x_b) + beta_ap H_ap(x_a)]
831 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
832 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
833 - beta_pb H_a(x_a) + beta_ap H_a(x_a) + beta_pb H_a(x_a) - beta_ap H_a(x_a)
834 - beta_ap H_b(x_b) + beta_bp H_b(x_b) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
835 = [(beta_bp H_bp(x_a) - beta_bp H_a(x_a)) - (beta_ap H_ap(x_a) - beta_ap H_a(x_a))] +
836 [(beta_ap H_ap(x_b) - beta_ap H_b(x_b)) - (beta_bp H_bp(x_b) - beta_bp H_b(x_b))]
837 + beta_pb H_a(x_a) - beta_ap H_a(x_a) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
838 = [beta_bp (H_bp(x_a) - H_a(x_a)) - beta_ap (H_ap(x_a) - H_a(x_a))] +
839 [beta_ap (H_ap(x_b) - H_b(x_b)) - beta_bp (H_bp(x_b) - H_b(x_b))]
840 + beta_pb (H_a(x_a) - H_b(x_b)) - beta_ap (H_a(x_a) - H_b(x_b))
841 = ([beta_bp de[bp][a] - beta_ap de[ap][a]) + beta_ap de[ap][b] - beta_bp de[bp][b])
842 + (beta_pb-beta_ap)(H_a(x_a) - H_b(x_b)) */
843 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]);
844 break;
845 default:
846 gmx_incons("Unknown replica exchange quantity");
848 if (bPrint)
850 fprintf(fplog, "Repl %d <-> %d dE_term = %10.3e (kT)\n", a, b, delta);
852 if (re->bNPT)
854 /* revist the calculation for 5.0. Might be some improvements. */
855 dpV = (beta[ap]*re->pres[ap]-beta[bp]*re->pres[bp])*(Vol[b]-Vol[a])/PRESFAC;
856 if (bPrint)
858 fprintf(fplog, " dpV = %10.3e d = %10.3e\n", dpV, delta + dpV);
860 delta += dpV;
862 return delta;
865 static void
866 test_for_replica_exchange(FILE *fplog,
867 const gmx_multisim_t *ms,
868 struct gmx_repl_ex *re,
869 gmx_enerdata_t *enerd,
870 real vol,
871 gmx_int64_t step,
872 real time)
874 int m, i, j, a, b, ap, bp, i0, i1, tmp;
875 real delta = 0;
876 gmx_bool bPrint, bMultiEx;
877 gmx_bool *bEx = re->bEx;
878 real *prob = re->prob;
879 int *pind = re->destinations; /* permuted index */
880 gmx_bool bEpot = FALSE;
881 gmx_bool bDLambda = FALSE;
882 gmx_bool bVol = FALSE;
884 bMultiEx = (re->nex > 1); /* multiple exchanges at each state */
885 fprintf(fplog, "Replica exchange at step %" GMX_PRId64 " time %.5f\n", step, time);
887 if (re->bNPT)
889 for (i = 0; i < re->nrepl; i++)
891 re->Vol[i] = 0;
893 bVol = TRUE;
894 re->Vol[re->repl] = vol;
896 if ((re->type == ereTEMP || re->type == ereTL))
898 for (i = 0; i < re->nrepl; i++)
900 re->Epot[i] = 0;
902 bEpot = TRUE;
903 re->Epot[re->repl] = enerd->term[F_EPOT];
904 /* temperatures of different states*/
905 for (i = 0; i < re->nrepl; i++)
907 re->beta[i] = 1.0/(re->q[ereTEMP][i]*BOLTZ);
910 else
912 for (i = 0; i < re->nrepl; i++)
914 re->beta[i] = 1.0/(re->temp*BOLTZ); /* we have a single temperature */
917 if (re->type == ereLAMBDA || re->type == ereTL)
919 bDLambda = TRUE;
920 /* lambda differences. */
921 /* de[i][j] is the energy of the jth simulation in the ith Hamiltonian
922 minus the energy of the jth simulation in the jth Hamiltonian */
923 for (i = 0; i < re->nrepl; i++)
925 for (j = 0; j < re->nrepl; j++)
927 re->de[i][j] = 0;
930 for (i = 0; i < re->nrepl; i++)
932 re->de[i][re->repl] = (enerd->enerpart_lambda[(int)re->q[ereLAMBDA][i]+1]-enerd->enerpart_lambda[0]);
936 /* now actually do the communication */
937 if (bVol)
939 gmx_sum_sim(re->nrepl, re->Vol, ms);
941 if (bEpot)
943 gmx_sum_sim(re->nrepl, re->Epot, ms);
945 if (bDLambda)
947 for (i = 0; i < re->nrepl; i++)
949 gmx_sum_sim(re->nrepl, re->de[i], ms);
953 /* make a duplicate set of indices for shuffling */
954 for (i = 0; i < re->nrepl; i++)
956 pind[i] = re->ind[i];
959 if (bMultiEx)
961 /* multiple random switch exchange */
962 int nself = 0;
963 for (i = 0; i < re->nex + nself; i++)
965 double rnd[2];
967 gmx_rng_cycle_2uniform(step, i*2, re->seed, RND_SEED_REPLEX, rnd);
968 /* randomly select a pair */
969 /* in theory, could reduce this by identifying only which switches had a nonneglibible
970 probability of occurring (log p > -100) and only operate on those switches */
971 /* find out which state it is from, and what label that state currently has. Likely
972 more work that useful. */
973 i0 = (int)(re->nrepl*rnd[0]);
974 i1 = (int)(re->nrepl*rnd[1]);
975 if (i0 == i1)
977 nself++;
978 continue; /* self-exchange, back up and do it again */
981 a = re->ind[i0]; /* what are the indices of these states? */
982 b = re->ind[i1];
983 ap = pind[i0];
984 bp = pind[i1];
986 bPrint = FALSE; /* too noisy */
987 /* calculate the energy difference */
988 /* if the code changes to flip the STATES, rather than the configurations,
989 use the commented version of the code */
990 /* delta = calc_delta(fplog,bPrint,re,a,b,ap,bp); */
991 delta = calc_delta(fplog, bPrint, re, ap, bp, a, b);
993 /* we actually only use the first space in the prob and bEx array,
994 since there are actually many switches between pairs. */
996 if (delta <= 0)
998 /* accepted */
999 prob[0] = 1;
1000 bEx[0] = TRUE;
1002 else
1004 if (delta > PROBABILITYCUTOFF)
1006 prob[0] = 0;
1008 else
1010 prob[0] = exp(-delta);
1012 /* roll a number to determine if accepted */
1013 gmx_rng_cycle_2uniform(step, i*2+1, re->seed, RND_SEED_REPLEX, rnd);
1014 bEx[0] = rnd[0] < prob[0];
1016 re->prob_sum[0] += prob[0];
1018 if (bEx[0])
1020 /* swap the states */
1021 tmp = pind[i0];
1022 pind[i0] = pind[i1];
1023 pind[i1] = tmp;
1026 re->nattempt[0]++; /* keep track of total permutation trials here */
1027 print_allswitchind(fplog, re->nrepl, pind, re->allswaps, re->tmpswap);
1029 else
1031 /* standard nearest neighbor replica exchange */
1033 m = (step / re->nst) % 2;
1034 for (i = 1; i < re->nrepl; i++)
1036 a = re->ind[i-1];
1037 b = re->ind[i];
1039 bPrint = (re->repl == a || re->repl == b);
1040 if (i % 2 == m)
1042 delta = calc_delta(fplog, bPrint, re, a, b, a, b);
1043 if (delta <= 0)
1045 /* accepted */
1046 prob[i] = 1;
1047 bEx[i] = TRUE;
1049 else
1051 double rnd[2];
1053 if (delta > PROBABILITYCUTOFF)
1055 prob[i] = 0;
1057 else
1059 prob[i] = exp(-delta);
1061 /* roll a number to determine if accepted */
1062 gmx_rng_cycle_2uniform(step, i, re->seed, RND_SEED_REPLEX, rnd);
1063 bEx[i] = rnd[0] < prob[i];
1065 re->prob_sum[i] += prob[i];
1067 if (bEx[i])
1069 /* swap these two */
1070 tmp = pind[i-1];
1071 pind[i-1] = pind[i];
1072 pind[i] = tmp;
1073 re->nexchange[i]++; /* statistics for back compatibility */
1076 else
1078 prob[i] = -1;
1079 bEx[i] = FALSE;
1082 /* print some statistics */
1083 print_ind(fplog, "ex", re->nrepl, re->ind, bEx);
1084 print_prob(fplog, "pr", re->nrepl, prob);
1085 fprintf(fplog, "\n");
1086 re->nattempt[m]++;
1089 /* record which moves were made and accepted */
1090 for (i = 0; i < re->nrepl; i++)
1092 re->nmoves[re->ind[i]][pind[i]] += 1;
1093 re->nmoves[pind[i]][re->ind[i]] += 1;
1095 fflush(fplog); /* make sure we can see what the last exchange was */
1098 static void
1099 cyclic_decomposition(const int *destinations,
1100 int **cyclic,
1101 gmx_bool *incycle,
1102 const int nrepl,
1103 int *nswap)
1106 int i, j, c, p;
1107 int maxlen = 1;
1108 for (i = 0; i < nrepl; i++)
1110 incycle[i] = FALSE;
1112 for (i = 0; i < nrepl; i++) /* one cycle for each replica */
1114 if (incycle[i])
1116 cyclic[i][0] = -1;
1117 continue;
1119 cyclic[i][0] = i;
1120 incycle[i] = TRUE;
1121 c = 1;
1122 p = i;
1123 for (j = 0; j < nrepl; j++) /* potentially all cycles are part, but we will break first */
1125 p = destinations[p]; /* start permuting */
1126 if (p == i)
1128 cyclic[i][c] = -1;
1129 if (c > maxlen)
1131 maxlen = c;
1133 break; /* we've reached the original element, the cycle is complete, and we marked the end. */
1135 else
1137 cyclic[i][c] = p; /* each permutation gives a new member of the cycle */
1138 incycle[p] = TRUE;
1139 c++;
1143 *nswap = maxlen - 1;
1145 if (debug)
1147 for (i = 0; i < nrepl; i++)
1149 fprintf(debug, "Cycle %d:", i);
1150 for (j = 0; j < nrepl; j++)
1152 if (cyclic[i][j] < 0)
1154 break;
1156 fprintf(debug, "%2d", cyclic[i][j]);
1158 fprintf(debug, "\n");
1160 fflush(debug);
1164 static void
1165 compute_exchange_order(int **cyclic,
1166 int **order,
1167 const int nrepl,
1168 const int maxswap)
1170 int i, j;
1172 for (j = 0; j < maxswap; j++)
1174 for (i = 0; i < nrepl; i++)
1176 if (cyclic[i][j+1] >= 0)
1178 order[cyclic[i][j+1]][j] = cyclic[i][j];
1179 order[cyclic[i][j]][j] = cyclic[i][j+1];
1182 for (i = 0; i < nrepl; i++)
1184 if (order[i][j] < 0)
1186 order[i][j] = i; /* if it's not exchanging, it should stay this round*/
1191 if (debug)
1193 fprintf(debug, "Replica Exchange Order\n");
1194 for (i = 0; i < nrepl; i++)
1196 fprintf(debug, "Replica %d:", i);
1197 for (j = 0; j < maxswap; j++)
1199 if (order[i][j] < 0)
1201 break;
1203 fprintf(debug, "%2d", order[i][j]);
1205 fprintf(debug, "\n");
1207 fflush(debug);
1211 static void
1212 prepare_to_do_exchange(struct gmx_repl_ex *re,
1213 const int replica_id,
1214 int *maxswap,
1215 gmx_bool *bThisReplicaExchanged)
1217 int i, j;
1218 /* Hold the cyclic decomposition of the (multiple) replica
1219 * exchange. */
1220 gmx_bool bAnyReplicaExchanged = FALSE;
1221 *bThisReplicaExchanged = FALSE;
1223 for (i = 0; i < re->nrepl; i++)
1225 if (re->destinations[i] != re->ind[i])
1227 /* only mark as exchanged if the index has been shuffled */
1228 bAnyReplicaExchanged = TRUE;
1229 break;
1232 if (bAnyReplicaExchanged)
1234 /* reinitialize the placeholder arrays */
1235 for (i = 0; i < re->nrepl; i++)
1237 for (j = 0; j < re->nrepl; j++)
1239 re->cyclic[i][j] = -1;
1240 re->order[i][j] = -1;
1244 /* Identify the cyclic decomposition of the permutation (very
1245 * fast if neighbor replica exchange). */
1246 cyclic_decomposition(re->destinations, re->cyclic, re->incycle, re->nrepl, maxswap);
1248 /* Now translate the decomposition into a replica exchange
1249 * order at each step. */
1250 compute_exchange_order(re->cyclic, re->order, re->nrepl, *maxswap);
1252 /* Did this replica do any exchange at any point? */
1253 for (j = 0; j < *maxswap; j++)
1255 if (replica_id != re->order[replica_id][j])
1257 *bThisReplicaExchanged = TRUE;
1258 break;
1264 gmx_bool replica_exchange(FILE *fplog, const t_commrec *cr, struct gmx_repl_ex *re,
1265 t_state *state, gmx_enerdata_t *enerd,
1266 t_state *state_local, gmx_int64_t step, real time)
1268 int j;
1269 int replica_id = 0;
1270 int exchange_partner;
1271 int maxswap = 0;
1272 /* Number of rounds of exchanges needed to deal with any multiple
1273 * exchanges. */
1274 /* Where each replica ends up after the exchange attempt(s). */
1275 /* The order in which multiple exchanges will occur. */
1276 gmx_bool bThisReplicaExchanged = FALSE;
1278 if (MASTER(cr))
1280 replica_id = re->repl;
1281 test_for_replica_exchange(fplog, cr->ms, re, enerd, det(state_local->box), step, time);
1282 prepare_to_do_exchange(re, replica_id, &maxswap, &bThisReplicaExchanged);
1284 /* Do intra-simulation broadcast so all processors belonging to
1285 * each simulation know whether they need to participate in
1286 * collecting the state. Otherwise, they might as well get on with
1287 * the next thing to do. */
1288 if (DOMAINDECOMP(cr))
1290 #ifdef GMX_MPI
1291 MPI_Bcast(&bThisReplicaExchanged, sizeof(gmx_bool), MPI_BYTE, MASTERRANK(cr),
1292 cr->mpi_comm_mygroup);
1293 #endif
1296 if (bThisReplicaExchanged)
1298 /* Exchange the states */
1299 /* Collect the global state on the master node */
1300 if (DOMAINDECOMP(cr))
1302 dd_collect_state(cr->dd, state_local, state);
1304 else
1306 copy_state_nonatomdata(state_local, state);
1309 if (MASTER(cr))
1311 /* There will be only one swap cycle with standard replica
1312 * exchange, but there may be multiple swap cycles if we
1313 * allow multiple swaps. */
1315 for (j = 0; j < maxswap; j++)
1317 exchange_partner = re->order[replica_id][j];
1319 if (exchange_partner != replica_id)
1321 /* Exchange the global states between the master nodes */
1322 if (debug)
1324 fprintf(debug, "Exchanging %d with %d\n", replica_id, exchange_partner);
1326 exchange_state(cr->ms, exchange_partner, state);
1329 /* For temperature-type replica exchange, we need to scale
1330 * the velocities. */
1331 if (re->type == ereTEMP || re->type == ereTL)
1333 scale_velocities(state, sqrt(re->q[ereTEMP][replica_id]/re->q[ereTEMP][re->destinations[replica_id]]));
1338 /* With domain decomposition the global state is distributed later */
1339 if (!DOMAINDECOMP(cr))
1341 /* Copy the global state to the local state data structure */
1342 copy_state_nonatomdata(state, state_local);
1346 return bThisReplicaExchanged;
1349 void print_replica_exchange_statistics(FILE *fplog, struct gmx_repl_ex *re)
1351 int i;
1353 fprintf(fplog, "\nReplica exchange statistics\n");
1355 if (re->nex == 0)
1357 fprintf(fplog, "Repl %d attempts, %d odd, %d even\n",
1358 re->nattempt[0]+re->nattempt[1], re->nattempt[1], re->nattempt[0]);
1360 fprintf(fplog, "Repl average probabilities:\n");
1361 for (i = 1; i < re->nrepl; i++)
1363 if (re->nattempt[i%2] == 0)
1365 re->prob[i] = 0;
1367 else
1369 re->prob[i] = re->prob_sum[i]/re->nattempt[i%2];
1372 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1373 print_prob(fplog, "", re->nrepl, re->prob);
1375 fprintf(fplog, "Repl number of exchanges:\n");
1376 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1377 print_count(fplog, "", re->nrepl, re->nexchange);
1379 fprintf(fplog, "Repl average number of exchanges:\n");
1380 for (i = 1; i < re->nrepl; i++)
1382 if (re->nattempt[i%2] == 0)
1384 re->prob[i] = 0;
1386 else
1388 re->prob[i] = ((real)re->nexchange[i])/re->nattempt[i%2];
1391 print_ind(fplog, "", re->nrepl, re->ind, NULL);
1392 print_prob(fplog, "", re->nrepl, re->prob);
1394 fprintf(fplog, "\n");
1396 /* print the transition matrix */
1397 print_transition_matrix(fplog, re->nrepl, re->nmoves, re->nattempt);