Parse user-supplied GPU task assignment only when needed
[gromacs.git] / src / programs / mdrun / repl_ex.cpp
blobc886b180e4e8f5dfcaf8e2fa8e1ec54b42f59fb8
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,2016,2017, 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 <random>
48 #include "gromacs/domdec/domdec.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/mdlib/main.h"
53 #include "gromacs/mdtypes/commrec.h"
54 #include "gromacs/mdtypes/inputrec.h"
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/mdtypes/state.h"
57 #include "gromacs/random/threefry.h"
58 #include "gromacs/random/uniformintdistribution.h"
59 #include "gromacs/random/uniformrealdistribution.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/pleasecite.h"
62 #include "gromacs/utility/smalloc.h"
65 #define PROBABILITYCUTOFF 100
66 /* we don't bother evaluating if events are more rare than exp(-100) = 3.7x10^-44 */
68 //! Rank in the multisimulaiton
69 #define MSRANK(ms, nodeid) (nodeid)
71 enum {
72 ereTEMP, ereLAMBDA, ereENDSINGLE, ereTL, ereNR
74 const char *erename[ereNR] = { "temperature", "lambda", "end_single_marker", "temperature and lambda"};
75 /* end_single_marker merely notes the end of single variable replica exchange. All types higher than
76 it are multiple replica exchange methods */
77 /* Eventually, should add 'pressure', 'temperature and pressure', 'lambda_and_pressure', 'temperature_lambda_pressure'?;
78 Let's wait until we feel better about the pressure control methods giving exact ensembles. Right now, we assume constant pressure */
80 typedef struct gmx_repl_ex
82 int repl; /* replica ID */
83 int nrepl; /* total number of replica */
84 real temp; /* temperature */
85 int type; /* replica exchange type from ere enum */
86 real **q; /* quantity, e.g. temperature or lambda; first index is ere, second index is replica ID */
87 gmx_bool bNPT; /* use constant pressure and temperature */
88 real *pres; /* replica pressures */
89 int *ind; /* replica indices */
90 int *allswaps; /* used for keeping track of all the replica swaps */
91 int nst; /* replica exchange interval (number of steps) */
92 int nex; /* number of exchanges per interval */
93 int seed; /* random seed */
94 int nattempt[2]; /* number of even and odd replica change attempts */
95 real *prob_sum; /* sum of probabilities */
96 int **nmoves; /* number of moves between replicas i and j */
97 int *nexchange; /* i-th element of the array is the number of exchanges between replica i-1 and i */
99 /* these are helper arrays for replica exchange; allocated here so they
100 don't have to be allocated each time */
101 int *destinations;
102 int **cyclic;
103 int **order;
104 int *tmpswap;
105 gmx_bool *incycle;
106 gmx_bool *bEx;
108 /* helper arrays to hold the quantities that are exchanged */
109 real *prob;
110 real *Epot;
111 real *beta;
112 real *Vol;
113 real **de;
115 } t_gmx_repl_ex;
117 static gmx_bool repl_quantity(const gmx_multisim_t *ms,
118 struct gmx_repl_ex *re, int ere, real q)
120 real *qall;
121 gmx_bool bDiff;
122 int s;
124 snew(qall, ms->nsim);
125 qall[re->repl] = q;
126 gmx_sum_sim(ms->nsim, qall, ms);
128 bDiff = FALSE;
129 for (s = 1; s < ms->nsim; s++)
131 if (qall[s] != qall[0])
133 bDiff = TRUE;
137 if (bDiff)
139 /* Set the replica exchange type and quantities */
140 re->type = ere;
142 snew(re->q[ere], re->nrepl);
143 for (s = 0; s < ms->nsim; s++)
145 re->q[ere][s] = qall[s];
148 sfree(qall);
149 return bDiff;
152 gmx_repl_ex_t
153 init_replica_exchange(FILE *fplog,
154 const gmx_multisim_t *ms,
155 const t_state *state,
156 const t_inputrec *ir,
157 const ReplicaExchangeParameters &replExParams)
159 real pres;
160 int i, j, k;
161 struct gmx_repl_ex *re;
162 gmx_bool bTemp;
163 gmx_bool bLambda = FALSE;
165 fprintf(fplog, "\nInitializing Replica Exchange\n");
167 if (ms == nullptr || ms->nsim == 1)
169 gmx_fatal(FARGS, "Nothing to exchange with only one replica, maybe you forgot to set the -multi option of mdrun?");
171 if (!EI_DYNAMICS(ir->eI))
173 gmx_fatal(FARGS, "Replica exchange is only supported by dynamical simulations");
174 /* Note that PAR(cr) is defined by cr->nnodes > 1, which is
175 * distinct from MULTISIM(cr). A multi-simulation only runs
176 * with real MPI parallelism, but this does not imply PAR(cr)
177 * is true!
179 * Since we are using a dynamical integrator, the only
180 * decomposition is DD, so PAR(cr) and DOMAINDECOMP(cr) are
181 * synonymous. The only way for cr->nnodes > 1 to be true is
182 * if we are using DD. */
185 snew(re, 1);
187 re->repl = ms->sim;
188 re->nrepl = ms->nsim;
189 snew(re->q, ereENDSINGLE);
191 fprintf(fplog, "Repl There are %d replicas:\n", re->nrepl);
193 check_multi_int(fplog, ms, state->natoms, "the number of atoms", FALSE);
194 check_multi_int(fplog, ms, ir->eI, "the integrator", FALSE);
195 check_multi_int64(fplog, ms, ir->init_step+ir->nsteps, "init_step+nsteps", FALSE);
196 const int nst = replExParams.exchangeInterval;
197 check_multi_int64(fplog, ms, (ir->init_step+nst-1)/nst,
198 "first exchange step: init_step/-replex", FALSE);
199 check_multi_int(fplog, ms, ir->etc, "the temperature coupling", FALSE);
200 check_multi_int(fplog, ms, ir->opts.ngtc,
201 "the number of temperature coupling groups", FALSE);
202 check_multi_int(fplog, ms, ir->epc, "the pressure coupling", FALSE);
203 check_multi_int(fplog, ms, ir->efep, "free energy", FALSE);
204 check_multi_int(fplog, ms, ir->fepvals->n_lambda, "number of lambda states", FALSE);
206 re->temp = ir->opts.ref_t[0];
207 for (i = 1; (i < ir->opts.ngtc); i++)
209 if (ir->opts.ref_t[i] != re->temp)
211 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
212 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
216 re->type = -1;
217 bTemp = repl_quantity(ms, re, ereTEMP, re->temp);
218 if (ir->efep != efepNO)
220 bLambda = repl_quantity(ms, re, ereLAMBDA, (real)ir->fepvals->init_fep_state);
222 if (re->type == -1) /* nothing was assigned */
224 gmx_fatal(FARGS, "The properties of the %d systems are all the same, there is nothing to exchange", re->nrepl);
226 if (bLambda && bTemp)
228 re->type = ereTL;
231 if (bTemp)
233 please_cite(fplog, "Sugita1999a");
234 if (ir->epc != epcNO)
236 re->bNPT = TRUE;
237 fprintf(fplog, "Repl Using Constant Pressure REMD.\n");
238 please_cite(fplog, "Okabe2001a");
240 if (ir->etc == etcBERENDSEN)
242 gmx_fatal(FARGS, "REMD with the %s thermostat does not produce correct potential energy distributions, consider using the %s thermostat instead",
243 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
246 if (bLambda)
248 if (ir->fepvals->delta_lambda != 0) /* check this? */
250 gmx_fatal(FARGS, "delta_lambda is not zero");
253 if (re->bNPT)
255 snew(re->pres, re->nrepl);
256 if (ir->epct == epctSURFACETENSION)
258 pres = ir->ref_p[ZZ][ZZ];
260 else
262 pres = 0;
263 j = 0;
264 for (i = 0; i < DIM; i++)
266 if (ir->compress[i][i] != 0)
268 pres += ir->ref_p[i][i];
269 j++;
272 pres /= j;
274 re->pres[re->repl] = pres;
275 gmx_sum_sim(re->nrepl, re->pres, ms);
278 /* Make an index for increasing replica order */
279 /* only makes sense if one or the other is varying, not both!
280 if both are varying, we trust the order the person gave. */
281 snew(re->ind, re->nrepl);
282 for (i = 0; i < re->nrepl; i++)
284 re->ind[i] = i;
287 if (re->type < ereENDSINGLE)
290 for (i = 0; i < re->nrepl; i++)
292 for (j = i+1; j < re->nrepl; j++)
294 if (re->q[re->type][re->ind[j]] < re->q[re->type][re->ind[i]])
296 /* Unordered replicas are supposed to work, but there
297 * is still an issues somewhere.
298 * Note that at this point still re->ind[i]=i.
300 gmx_fatal(FARGS, "Replicas with indices %d < %d have %ss %g > %g, please order your replicas on increasing %s",
301 i, j,
302 erename[re->type],
303 re->q[re->type][i], re->q[re->type][j],
304 erename[re->type]);
306 k = re->ind[i];
307 re->ind[i] = re->ind[j];
308 re->ind[j] = k;
310 else if (re->q[re->type][re->ind[j]] == re->q[re->type][re->ind[i]])
312 gmx_fatal(FARGS, "Two replicas have identical %ss", erename[re->type]);
318 /* keep track of all the swaps, starting with the initial placement. */
319 snew(re->allswaps, re->nrepl);
320 for (i = 0; i < re->nrepl; i++)
322 re->allswaps[i] = re->ind[i];
325 switch (re->type)
327 case ereTEMP:
328 fprintf(fplog, "\nReplica exchange in temperature\n");
329 for (i = 0; i < re->nrepl; i++)
331 fprintf(fplog, " %5.1f", re->q[re->type][re->ind[i]]);
333 fprintf(fplog, "\n");
334 break;
335 case ereLAMBDA:
336 fprintf(fplog, "\nReplica exchange in lambda\n");
337 for (i = 0; i < re->nrepl; i++)
339 fprintf(fplog, " %3d", (int)re->q[re->type][re->ind[i]]);
341 fprintf(fplog, "\n");
342 break;
343 case ereTL:
344 fprintf(fplog, "\nReplica exchange in temperature and lambda state\n");
345 for (i = 0; i < re->nrepl; i++)
347 fprintf(fplog, " %5.1f", re->q[ereTEMP][re->ind[i]]);
349 fprintf(fplog, "\n");
350 for (i = 0; i < re->nrepl; i++)
352 fprintf(fplog, " %5d", (int)re->q[ereLAMBDA][re->ind[i]]);
354 fprintf(fplog, "\n");
355 break;
356 default:
357 gmx_incons("Unknown replica exchange quantity");
359 if (re->bNPT)
361 fprintf(fplog, "\nRepl p");
362 for (i = 0; i < re->nrepl; i++)
364 fprintf(fplog, " %5.2f", re->pres[re->ind[i]]);
367 for (i = 0; i < re->nrepl; i++)
369 if ((i > 0) && (re->pres[re->ind[i]] < re->pres[re->ind[i-1]]))
371 fprintf(fplog, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
372 fprintf(stderr, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
376 re->nst = nst;
377 if (replExParams.randomSeed == -1)
379 if (MASTERSIM(ms))
381 re->seed = static_cast<int>(gmx::makeRandomSeed());
383 else
385 re->seed = 0;
387 gmx_sumi_sim(1, &(re->seed), ms);
389 else
391 re->seed = replExParams.randomSeed;
393 fprintf(fplog, "\nReplica exchange interval: %d\n", re->nst);
394 fprintf(fplog, "\nReplica random seed: %d\n", re->seed);
396 re->nattempt[0] = 0;
397 re->nattempt[1] = 0;
399 snew(re->prob_sum, re->nrepl);
400 snew(re->nexchange, re->nrepl);
401 snew(re->nmoves, re->nrepl);
402 for (i = 0; i < re->nrepl; i++)
404 snew(re->nmoves[i], re->nrepl);
406 fprintf(fplog, "Replica exchange information below: ex and x = exchange, pr = probability\n");
408 /* generate space for the helper functions so we don't have to snew each time */
410 snew(re->destinations, re->nrepl);
411 snew(re->incycle, re->nrepl);
412 snew(re->tmpswap, re->nrepl);
413 snew(re->cyclic, re->nrepl);
414 snew(re->order, re->nrepl);
415 for (i = 0; i < re->nrepl; i++)
417 snew(re->cyclic[i], re->nrepl+1);
418 snew(re->order[i], re->nrepl);
420 /* allocate space for the functions storing the data for the replicas */
421 /* not all of these arrays needed in all cases, but they don't take
422 up much space, since the max size is nrepl**2 */
423 snew(re->prob, re->nrepl);
424 snew(re->bEx, re->nrepl);
425 snew(re->beta, re->nrepl);
426 snew(re->Vol, re->nrepl);
427 snew(re->Epot, re->nrepl);
428 snew(re->de, re->nrepl);
429 for (i = 0; i < re->nrepl; i++)
431 snew(re->de[i], re->nrepl);
433 re->nex = replExParams.numExchanges;
434 return re;
437 static void exchange_reals(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, real *v, int n)
439 real *buf;
440 int i;
442 if (v)
444 snew(buf, n);
445 #if GMX_MPI
447 MPI_Sendrecv(v, n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
448 buf,n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
449 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
452 MPI_Request mpi_req;
454 MPI_Isend(v, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
455 ms->mpi_comm_masters, &mpi_req);
456 MPI_Recv(buf, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
457 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
458 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
460 #endif
461 for (i = 0; i < n; i++)
463 v[i] = buf[i];
465 sfree(buf);
470 static void exchange_doubles(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, double *v, int n)
472 double *buf;
473 int i;
475 if (v)
477 snew(buf, n);
478 #if GMX_MPI
480 MPI_Sendrecv(v, n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
481 buf,n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
482 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
485 MPI_Request mpi_req;
487 MPI_Isend(v, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
488 ms->mpi_comm_masters, &mpi_req);
489 MPI_Recv(buf, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
490 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
491 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
493 #endif
494 for (i = 0; i < n; i++)
496 v[i] = buf[i];
498 sfree(buf);
502 static void exchange_rvecs(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, rvec *v, int n)
504 rvec *buf;
505 int i;
507 if (v)
509 snew(buf, n);
510 #if GMX_MPI
512 MPI_Sendrecv(v[0], n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
513 buf[0],n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
514 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
517 MPI_Request mpi_req;
519 MPI_Isend(v[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
520 ms->mpi_comm_masters, &mpi_req);
521 MPI_Recv(buf[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
522 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
523 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
525 #endif
526 for (i = 0; i < n; i++)
528 copy_rvec(buf[i], v[i]);
530 sfree(buf);
534 static void exchange_state(const gmx_multisim_t *ms, int b, t_state *state)
536 /* When t_state changes, this code should be updated. */
537 int ngtc, nnhpres;
538 ngtc = state->ngtc * state->nhchainlength;
539 nnhpres = state->nnhpres* state->nhchainlength;
540 exchange_rvecs(ms, b, state->box, DIM);
541 exchange_rvecs(ms, b, state->box_rel, DIM);
542 exchange_rvecs(ms, b, state->boxv, DIM);
543 exchange_reals(ms, b, &(state->veta), 1);
544 exchange_reals(ms, b, &(state->vol0), 1);
545 exchange_rvecs(ms, b, state->svir_prev, DIM);
546 exchange_rvecs(ms, b, state->fvir_prev, DIM);
547 exchange_rvecs(ms, b, state->pres_prev, DIM);
548 exchange_doubles(ms, b, state->nosehoover_xi.data(), ngtc);
549 exchange_doubles(ms, b, state->nosehoover_vxi.data(), ngtc);
550 exchange_doubles(ms, b, state->nhpres_xi.data(), nnhpres);
551 exchange_doubles(ms, b, state->nhpres_vxi.data(), nnhpres);
552 exchange_doubles(ms, b, state->therm_integral.data(), state->ngtc);
553 exchange_doubles(ms, b, &state->baros_integral, 1);
554 exchange_rvecs(ms, b, as_rvec_array(state->x.data()), state->natoms);
555 exchange_rvecs(ms, b, as_rvec_array(state->v.data()), state->natoms);
558 static void copy_state_serial(const t_state *src, t_state *dest)
560 if (dest != src)
562 /* Currently the local state is always a pointer to the global
563 * in serial, so we should never end up here.
564 * TODO: Implement a (trivial) t_state copy once converted to C++.
566 GMX_RELEASE_ASSERT(false, "State copying is currently not implemented in replica exchange");
570 static void scale_velocities(t_state *state, real fac)
572 int i;
574 if (as_rvec_array(state->v.data()))
576 for (i = 0; i < state->natoms; i++)
578 svmul(fac, state->v[i], state->v[i]);
583 static void print_transition_matrix(FILE *fplog, int n, int **nmoves, int *nattempt)
585 int i, j, ntot;
586 float Tprint;
588 ntot = nattempt[0] + nattempt[1];
589 fprintf(fplog, "\n");
590 fprintf(fplog, "Repl");
591 for (i = 0; i < n; i++)
593 fprintf(fplog, " "); /* put the title closer to the center */
595 fprintf(fplog, "Empirical Transition Matrix\n");
597 fprintf(fplog, "Repl");
598 for (i = 0; i < n; i++)
600 fprintf(fplog, "%8d", (i+1));
602 fprintf(fplog, "\n");
604 for (i = 0; i < n; i++)
606 fprintf(fplog, "Repl");
607 for (j = 0; j < n; j++)
609 Tprint = 0.0;
610 if (nmoves[i][j] > 0)
612 Tprint = nmoves[i][j]/(2.0*ntot);
614 fprintf(fplog, "%8.4f", Tprint);
616 fprintf(fplog, "%3d\n", i);
620 static void print_ind(FILE *fplog, const char *leg, int n, int *ind, gmx_bool *bEx)
622 int i;
624 fprintf(fplog, "Repl %2s %2d", leg, ind[0]);
625 for (i = 1; i < n; i++)
627 fprintf(fplog, " %c %2d", (bEx != nullptr && bEx[i]) ? 'x' : ' ', ind[i]);
629 fprintf(fplog, "\n");
632 static void print_allswitchind(FILE *fplog, int n, int *pind, int *allswaps, int *tmpswap)
634 int i;
636 for (i = 0; i < n; i++)
638 tmpswap[i] = allswaps[i];
640 for (i = 0; i < n; i++)
642 allswaps[i] = tmpswap[pind[i]];
645 fprintf(fplog, "\nAccepted Exchanges: ");
646 for (i = 0; i < n; i++)
648 fprintf(fplog, "%d ", pind[i]);
650 fprintf(fplog, "\n");
652 /* the "Order After Exchange" is the state label corresponding to the configuration that
653 started in state listed in order, i.e.
655 3 0 1 2
657 means that the:
658 configuration starting in simulation 3 is now in simulation 0,
659 configuration starting in simulation 0 is now in simulation 1,
660 configuration starting in simulation 1 is now in simulation 2,
661 configuration starting in simulation 2 is now in simulation 3
663 fprintf(fplog, "Order After Exchange: ");
664 for (i = 0; i < n; i++)
666 fprintf(fplog, "%d ", allswaps[i]);
668 fprintf(fplog, "\n\n");
671 static void print_prob(FILE *fplog, const char *leg, int n, real *prob)
673 int i;
674 char buf[8];
676 fprintf(fplog, "Repl %2s ", leg);
677 for (i = 1; i < n; i++)
679 if (prob[i] >= 0)
681 sprintf(buf, "%4.2f", prob[i]);
682 fprintf(fplog, " %3s", buf[0] == '1' ? "1.0" : buf+1);
684 else
686 fprintf(fplog, " ");
689 fprintf(fplog, "\n");
692 static void print_count(FILE *fplog, const char *leg, int n, int *count)
694 int i;
696 fprintf(fplog, "Repl %2s ", leg);
697 for (i = 1; i < n; i++)
699 fprintf(fplog, " %4d", count[i]);
701 fprintf(fplog, "\n");
704 static real calc_delta(FILE *fplog, gmx_bool bPrint, struct gmx_repl_ex *re, int a, int b, int ap, int bp)
707 real ediff, dpV, delta = 0;
708 real *Epot = re->Epot;
709 real *Vol = re->Vol;
710 real **de = re->de;
711 real *beta = re->beta;
713 /* Two cases; we are permuted and not. In all cases, setting ap = a and bp = b will reduce
714 to the non permuted case */
716 switch (re->type)
718 case ereTEMP:
720 * Okabe et. al. Chem. Phys. Lett. 335 (2001) 435-439
722 ediff = Epot[b] - Epot[a];
723 delta = -(beta[bp] - beta[ap])*ediff;
724 break;
725 case ereLAMBDA:
726 /* two cases: when we are permuted, and not. */
727 /* non-permuted:
728 ediff = E_new - E_old
729 = [H_b(x_a) + H_a(x_b)] - [H_b(x_b) + H_a(x_a)]
730 = [H_b(x_a) - H_a(x_a)] + [H_a(x_b) - H_b(x_b)]
731 = de[b][a] + de[a][b] */
733 /* permuted:
734 ediff = E_new - E_old
735 = [H_bp(x_a) + H_ap(x_b)] - [H_bp(x_b) + H_ap(x_a)]
736 = [H_bp(x_a) - H_ap(x_a)] + [H_ap(x_b) - H_bp(x_b)]
737 = [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)]
738 = [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)]
739 = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]) */
740 /* but, in the current code implementation, we flip configurations, not indices . . .
741 So let's examine that.
742 = [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)]
743 = [H_b(x_ap) - H_a(x_ap)] + [H_a(x_bp) - H_b(x_pb)]
744 = (de[b][ap] - de[a][ap]) + (de[a][bp] - de[b][bp]
745 So, if we exchange b<=> bp and a<=> ap, we return to the same result.
746 So the simple solution is to flip the
747 position of perturbed and original indices in the tests.
750 ediff = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]);
751 delta = ediff*beta[a]; /* assume all same temperature in this case */
752 break;
753 case ereTL:
754 /* not permuted: */
755 /* delta = reduced E_new - reduced E_old
756 = [beta_b H_b(x_a) + beta_a H_a(x_b)] - [beta_b H_b(x_b) + beta_a H_a(x_a)]
757 = [beta_b H_b(x_a) - beta_a H_a(x_a)] + [beta_a H_a(x_b) - beta_b H_b(x_b)]
758 = [beta_b dH_b(x_a) + beta_b H_a(x_a) - beta_a H_a(x_a)] +
759 [beta_a dH_a(x_b) + beta_a H_b(x_b) - beta_b H_b(x_b)]
760 = [beta_b dH_b(x_a) + [beta_a dH_a(x_b) +
761 beta_b (H_a(x_a) - H_b(x_b)]) - beta_a (H_a(x_a) - H_b(x_b))
762 = beta_b dH_b(x_a) + beta_a dH_a(x_b) - (beta_b - beta_a)(H_b(x_b) - H_a(x_a) */
763 /* delta = beta[b]*de[b][a] + beta[a]*de[a][b] - (beta[b] - beta[a])*(Epot[b] - Epot[a]; */
764 /* permuted (big breath!) */
765 /* delta = reduced E_new - reduced E_old
766 = [beta_bp H_bp(x_a) + beta_ap H_ap(x_b)] - [beta_bp H_bp(x_b) + beta_ap H_ap(x_a)]
767 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
768 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
769 - beta_pb H_a(x_a) + beta_ap H_a(x_a) + beta_pb H_a(x_a) - beta_ap H_a(x_a)
770 - beta_ap H_b(x_b) + beta_bp H_b(x_b) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
771 = [(beta_bp H_bp(x_a) - beta_bp H_a(x_a)) - (beta_ap H_ap(x_a) - beta_ap H_a(x_a))] +
772 [(beta_ap H_ap(x_b) - beta_ap H_b(x_b)) - (beta_bp H_bp(x_b) - beta_bp H_b(x_b))]
773 + beta_pb H_a(x_a) - beta_ap H_a(x_a) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
774 = [beta_bp (H_bp(x_a) - H_a(x_a)) - beta_ap (H_ap(x_a) - H_a(x_a))] +
775 [beta_ap (H_ap(x_b) - H_b(x_b)) - beta_bp (H_bp(x_b) - H_b(x_b))]
776 + beta_pb (H_a(x_a) - H_b(x_b)) - beta_ap (H_a(x_a) - H_b(x_b))
777 = ([beta_bp de[bp][a] - beta_ap de[ap][a]) + beta_ap de[ap][b] - beta_bp de[bp][b])
778 + (beta_pb-beta_ap)(H_a(x_a) - H_b(x_b)) */
779 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]);
780 break;
781 default:
782 gmx_incons("Unknown replica exchange quantity");
784 if (bPrint)
786 fprintf(fplog, "Repl %d <-> %d dE_term = %10.3e (kT)\n", a, b, delta);
788 if (re->bNPT)
790 /* revist the calculation for 5.0. Might be some improvements. */
791 dpV = (beta[ap]*re->pres[ap]-beta[bp]*re->pres[bp])*(Vol[b]-Vol[a])/PRESFAC;
792 if (bPrint)
794 fprintf(fplog, " dpV = %10.3e d = %10.3e\n", dpV, delta + dpV);
796 delta += dpV;
798 return delta;
801 static void
802 test_for_replica_exchange(FILE *fplog,
803 const gmx_multisim_t *ms,
804 struct gmx_repl_ex *re,
805 gmx_enerdata_t *enerd,
806 real vol,
807 gmx_int64_t step,
808 real time)
810 int m, i, j, a, b, ap, bp, i0, i1, tmp;
811 real delta = 0;
812 gmx_bool bPrint, bMultiEx;
813 gmx_bool *bEx = re->bEx;
814 real *prob = re->prob;
815 int *pind = re->destinations; /* permuted index */
816 gmx_bool bEpot = FALSE;
817 gmx_bool bDLambda = FALSE;
818 gmx_bool bVol = FALSE;
819 gmx::ThreeFry2x64<64> rng(re->seed, gmx::RandomDomain::ReplicaExchange);
820 gmx::UniformRealDistribution<real> uniformRealDist;
821 gmx::UniformIntDistribution<int> uniformNreplDist(0, re->nrepl-1);
823 bMultiEx = (re->nex > 1); /* multiple exchanges at each state */
824 fprintf(fplog, "Replica exchange at step %" GMX_PRId64 " time %.5f\n", step, time);
826 if (re->bNPT)
828 for (i = 0; i < re->nrepl; i++)
830 re->Vol[i] = 0;
832 bVol = TRUE;
833 re->Vol[re->repl] = vol;
835 if ((re->type == ereTEMP || re->type == ereTL))
837 for (i = 0; i < re->nrepl; i++)
839 re->Epot[i] = 0;
841 bEpot = TRUE;
842 re->Epot[re->repl] = enerd->term[F_EPOT];
843 /* temperatures of different states*/
844 for (i = 0; i < re->nrepl; i++)
846 re->beta[i] = 1.0/(re->q[ereTEMP][i]*BOLTZ);
849 else
851 for (i = 0; i < re->nrepl; i++)
853 re->beta[i] = 1.0/(re->temp*BOLTZ); /* we have a single temperature */
856 if (re->type == ereLAMBDA || re->type == ereTL)
858 bDLambda = TRUE;
859 /* lambda differences. */
860 /* de[i][j] is the energy of the jth simulation in the ith Hamiltonian
861 minus the energy of the jth simulation in the jth Hamiltonian */
862 for (i = 0; i < re->nrepl; i++)
864 for (j = 0; j < re->nrepl; j++)
866 re->de[i][j] = 0;
869 for (i = 0; i < re->nrepl; i++)
871 re->de[i][re->repl] = (enerd->enerpart_lambda[(int)re->q[ereLAMBDA][i]+1]-enerd->enerpart_lambda[0]);
875 /* now actually do the communication */
876 if (bVol)
878 gmx_sum_sim(re->nrepl, re->Vol, ms);
880 if (bEpot)
882 gmx_sum_sim(re->nrepl, re->Epot, ms);
884 if (bDLambda)
886 for (i = 0; i < re->nrepl; i++)
888 gmx_sum_sim(re->nrepl, re->de[i], ms);
892 /* make a duplicate set of indices for shuffling */
893 for (i = 0; i < re->nrepl; i++)
895 pind[i] = re->ind[i];
898 rng.restart( step, 0 );
900 if (bMultiEx)
902 /* multiple random switch exchange */
903 int nself = 0;
906 for (i = 0; i < re->nex + nself; i++)
908 // For now this is superfluous, but just in case we ever add more
909 // calls in different branches it is safer to always reset the distribution.
910 uniformNreplDist.reset();
912 /* randomly select a pair */
913 /* in theory, could reduce this by identifying only which switches had a nonneglibible
914 probability of occurring (log p > -100) and only operate on those switches */
915 /* find out which state it is from, and what label that state currently has. Likely
916 more work that useful. */
917 i0 = uniformNreplDist(rng);
918 i1 = uniformNreplDist(rng);
919 if (i0 == i1)
921 nself++;
922 continue; /* self-exchange, back up and do it again */
925 a = re->ind[i0]; /* what are the indices of these states? */
926 b = re->ind[i1];
927 ap = pind[i0];
928 bp = pind[i1];
930 bPrint = FALSE; /* too noisy */
931 /* calculate the energy difference */
932 /* if the code changes to flip the STATES, rather than the configurations,
933 use the commented version of the code */
934 /* delta = calc_delta(fplog,bPrint,re,a,b,ap,bp); */
935 delta = calc_delta(fplog, bPrint, re, ap, bp, a, b);
937 /* we actually only use the first space in the prob and bEx array,
938 since there are actually many switches between pairs. */
940 if (delta <= 0)
942 /* accepted */
943 prob[0] = 1;
944 bEx[0] = TRUE;
946 else
948 if (delta > PROBABILITYCUTOFF)
950 prob[0] = 0;
952 else
954 prob[0] = exp(-delta);
956 // roll a number to determine if accepted. For now it is superfluous to
957 // reset, but just in case we ever add more calls in different branches
958 // it is safer to always reset the distribution.
959 uniformRealDist.reset();
960 bEx[0] = uniformRealDist(rng) < prob[0];
962 re->prob_sum[0] += prob[0];
964 if (bEx[0])
966 /* swap the states */
967 tmp = pind[i0];
968 pind[i0] = pind[i1];
969 pind[i1] = tmp;
972 re->nattempt[0]++; /* keep track of total permutation trials here */
973 print_allswitchind(fplog, re->nrepl, pind, re->allswaps, re->tmpswap);
975 else
977 /* standard nearest neighbor replica exchange */
979 m = (step / re->nst) % 2;
980 for (i = 1; i < re->nrepl; i++)
982 a = re->ind[i-1];
983 b = re->ind[i];
985 bPrint = (re->repl == a || re->repl == b);
986 if (i % 2 == m)
988 delta = calc_delta(fplog, bPrint, re, a, b, a, b);
989 if (delta <= 0)
991 /* accepted */
992 prob[i] = 1;
993 bEx[i] = TRUE;
995 else
997 if (delta > PROBABILITYCUTOFF)
999 prob[i] = 0;
1001 else
1003 prob[i] = exp(-delta);
1005 // roll a number to determine if accepted. For now it is superfluous to
1006 // reset, but just in case we ever add more calls in different branches
1007 // it is safer to always reset the distribution.
1008 uniformRealDist.reset();
1009 bEx[i] = uniformRealDist(rng) < prob[i];
1011 re->prob_sum[i] += prob[i];
1013 if (bEx[i])
1015 /* swap these two */
1016 tmp = pind[i-1];
1017 pind[i-1] = pind[i];
1018 pind[i] = tmp;
1019 re->nexchange[i]++; /* statistics for back compatibility */
1022 else
1024 prob[i] = -1;
1025 bEx[i] = FALSE;
1028 /* print some statistics */
1029 print_ind(fplog, "ex", re->nrepl, re->ind, bEx);
1030 print_prob(fplog, "pr", re->nrepl, prob);
1031 fprintf(fplog, "\n");
1032 re->nattempt[m]++;
1035 /* record which moves were made and accepted */
1036 for (i = 0; i < re->nrepl; i++)
1038 re->nmoves[re->ind[i]][pind[i]] += 1;
1039 re->nmoves[pind[i]][re->ind[i]] += 1;
1041 fflush(fplog); /* make sure we can see what the last exchange was */
1044 static void
1045 cyclic_decomposition(const int *destinations,
1046 int **cyclic,
1047 gmx_bool *incycle,
1048 const int nrepl,
1049 int *nswap)
1052 int i, j, c, p;
1053 int maxlen = 1;
1054 for (i = 0; i < nrepl; i++)
1056 incycle[i] = FALSE;
1058 for (i = 0; i < nrepl; i++) /* one cycle for each replica */
1060 if (incycle[i])
1062 cyclic[i][0] = -1;
1063 continue;
1065 cyclic[i][0] = i;
1066 incycle[i] = TRUE;
1067 c = 1;
1068 p = i;
1069 for (j = 0; j < nrepl; j++) /* potentially all cycles are part, but we will break first */
1071 p = destinations[p]; /* start permuting */
1072 if (p == i)
1074 cyclic[i][c] = -1;
1075 if (c > maxlen)
1077 maxlen = c;
1079 break; /* we've reached the original element, the cycle is complete, and we marked the end. */
1081 else
1083 cyclic[i][c] = p; /* each permutation gives a new member of the cycle */
1084 incycle[p] = TRUE;
1085 c++;
1089 *nswap = maxlen - 1;
1091 if (debug)
1093 for (i = 0; i < nrepl; i++)
1095 fprintf(debug, "Cycle %d:", i);
1096 for (j = 0; j < nrepl; j++)
1098 if (cyclic[i][j] < 0)
1100 break;
1102 fprintf(debug, "%2d", cyclic[i][j]);
1104 fprintf(debug, "\n");
1106 fflush(debug);
1110 static void
1111 compute_exchange_order(int **cyclic,
1112 int **order,
1113 const int nrepl,
1114 const int maxswap)
1116 int i, j;
1118 for (j = 0; j < maxswap; j++)
1120 for (i = 0; i < nrepl; i++)
1122 if (cyclic[i][j+1] >= 0)
1124 order[cyclic[i][j+1]][j] = cyclic[i][j];
1125 order[cyclic[i][j]][j] = cyclic[i][j+1];
1128 for (i = 0; i < nrepl; i++)
1130 if (order[i][j] < 0)
1132 order[i][j] = i; /* if it's not exchanging, it should stay this round*/
1137 if (debug)
1139 fprintf(debug, "Replica Exchange Order\n");
1140 for (i = 0; i < nrepl; i++)
1142 fprintf(debug, "Replica %d:", i);
1143 for (j = 0; j < maxswap; j++)
1145 if (order[i][j] < 0)
1147 break;
1149 fprintf(debug, "%2d", order[i][j]);
1151 fprintf(debug, "\n");
1153 fflush(debug);
1157 static void
1158 prepare_to_do_exchange(struct gmx_repl_ex *re,
1159 const int replica_id,
1160 int *maxswap,
1161 gmx_bool *bThisReplicaExchanged)
1163 int i, j;
1164 /* Hold the cyclic decomposition of the (multiple) replica
1165 * exchange. */
1166 gmx_bool bAnyReplicaExchanged = FALSE;
1167 *bThisReplicaExchanged = FALSE;
1169 for (i = 0; i < re->nrepl; i++)
1171 if (re->destinations[i] != re->ind[i])
1173 /* only mark as exchanged if the index has been shuffled */
1174 bAnyReplicaExchanged = TRUE;
1175 break;
1178 if (bAnyReplicaExchanged)
1180 /* reinitialize the placeholder arrays */
1181 for (i = 0; i < re->nrepl; i++)
1183 for (j = 0; j < re->nrepl; j++)
1185 re->cyclic[i][j] = -1;
1186 re->order[i][j] = -1;
1190 /* Identify the cyclic decomposition of the permutation (very
1191 * fast if neighbor replica exchange). */
1192 cyclic_decomposition(re->destinations, re->cyclic, re->incycle, re->nrepl, maxswap);
1194 /* Now translate the decomposition into a replica exchange
1195 * order at each step. */
1196 compute_exchange_order(re->cyclic, re->order, re->nrepl, *maxswap);
1198 /* Did this replica do any exchange at any point? */
1199 for (j = 0; j < *maxswap; j++)
1201 if (replica_id != re->order[replica_id][j])
1203 *bThisReplicaExchanged = TRUE;
1204 break;
1210 gmx_bool replica_exchange(FILE *fplog, const t_commrec *cr, struct gmx_repl_ex *re,
1211 t_state *state, gmx_enerdata_t *enerd,
1212 t_state *state_local, gmx_int64_t step, real time)
1214 int j;
1215 int replica_id = 0;
1216 int exchange_partner;
1217 int maxswap = 0;
1218 /* Number of rounds of exchanges needed to deal with any multiple
1219 * exchanges. */
1220 /* Where each replica ends up after the exchange attempt(s). */
1221 /* The order in which multiple exchanges will occur. */
1222 gmx_bool bThisReplicaExchanged = FALSE;
1224 if (MASTER(cr))
1226 replica_id = re->repl;
1227 test_for_replica_exchange(fplog, cr->ms, re, enerd, det(state_local->box), step, time);
1228 prepare_to_do_exchange(re, replica_id, &maxswap, &bThisReplicaExchanged);
1230 /* Do intra-simulation broadcast so all processors belonging to
1231 * each simulation know whether they need to participate in
1232 * collecting the state. Otherwise, they might as well get on with
1233 * the next thing to do. */
1234 if (DOMAINDECOMP(cr))
1236 #if GMX_MPI
1237 MPI_Bcast(&bThisReplicaExchanged, sizeof(gmx_bool), MPI_BYTE, MASTERRANK(cr),
1238 cr->mpi_comm_mygroup);
1239 #endif
1242 if (bThisReplicaExchanged)
1244 /* Exchange the states */
1245 /* Collect the global state on the master node */
1246 if (DOMAINDECOMP(cr))
1248 dd_collect_state(cr->dd, state_local, state);
1250 else
1252 copy_state_serial(state_local, state);
1255 if (MASTER(cr))
1257 /* There will be only one swap cycle with standard replica
1258 * exchange, but there may be multiple swap cycles if we
1259 * allow multiple swaps. */
1261 for (j = 0; j < maxswap; j++)
1263 exchange_partner = re->order[replica_id][j];
1265 if (exchange_partner != replica_id)
1267 /* Exchange the global states between the master nodes */
1268 if (debug)
1270 fprintf(debug, "Exchanging %d with %d\n", replica_id, exchange_partner);
1272 exchange_state(cr->ms, exchange_partner, state);
1275 /* For temperature-type replica exchange, we need to scale
1276 * the velocities. */
1277 if (re->type == ereTEMP || re->type == ereTL)
1279 scale_velocities(state, sqrt(re->q[ereTEMP][replica_id]/re->q[ereTEMP][re->destinations[replica_id]]));
1284 /* With domain decomposition the global state is distributed later */
1285 if (!DOMAINDECOMP(cr))
1287 /* Copy the global state to the local state data structure */
1288 copy_state_serial(state, state_local);
1292 return bThisReplicaExchanged;
1295 void print_replica_exchange_statistics(FILE *fplog, struct gmx_repl_ex *re)
1297 int i;
1299 fprintf(fplog, "\nReplica exchange statistics\n");
1301 if (re->nex == 0)
1303 fprintf(fplog, "Repl %d attempts, %d odd, %d even\n",
1304 re->nattempt[0]+re->nattempt[1], re->nattempt[1], re->nattempt[0]);
1306 fprintf(fplog, "Repl average probabilities:\n");
1307 for (i = 1; i < re->nrepl; i++)
1309 if (re->nattempt[i%2] == 0)
1311 re->prob[i] = 0;
1313 else
1315 re->prob[i] = re->prob_sum[i]/re->nattempt[i%2];
1318 print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1319 print_prob(fplog, "", re->nrepl, re->prob);
1321 fprintf(fplog, "Repl number of exchanges:\n");
1322 print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1323 print_count(fplog, "", re->nrepl, re->nexchange);
1325 fprintf(fplog, "Repl average number of exchanges:\n");
1326 for (i = 1; i < re->nrepl; i++)
1328 if (re->nattempt[i%2] == 0)
1330 re->prob[i] = 0;
1332 else
1334 re->prob[i] = ((real)re->nexchange[i])/re->nattempt[i%2];
1337 print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1338 print_prob(fplog, "", re->nrepl, re->prob);
1340 fprintf(fplog, "\n");
1342 /* print the transition matrix */
1343 print_transition_matrix(fplog, re->nrepl, re->nmoves, re->nattempt);