Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / mdrun / replicaexchange.cpp
blob9ff7fa478b3ec93add58a617a34824f1ad68e818
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,2018,2019, 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 /*! \internal \file
40 * \brief Implements the replica exchange routines.
42 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
43 * \author Mark Abraham <mark.j.abraham@gmail.com>
44 * \ingroup module_mdrun
46 #include "gmxpre.h"
48 #include "replicaexchange.h"
50 #include "config.h"
52 #include <cmath>
54 #include <random>
56 #include "gromacs/domdec/collect.h"
57 #include "gromacs/gmxlib/network.h"
58 #include "gromacs/math/units.h"
59 #include "gromacs/math/vec.h"
60 #include "gromacs/mdrunutility/multisim.h"
61 #include "gromacs/mdtypes/commrec.h"
62 #include "gromacs/mdtypes/enerdata.h"
63 #include "gromacs/mdtypes/inputrec.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/mdtypes/state.h"
66 #include "gromacs/random/threefry.h"
67 #include "gromacs/random/uniformintdistribution.h"
68 #include "gromacs/random/uniformrealdistribution.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/pleasecite.h"
71 #include "gromacs/utility/smalloc.h"
73 //! Helps cut off probability values.
74 constexpr int c_probabilityCutoff = 100;
76 /* we don't bother evaluating if events are more rare than exp(-100) = 3.7x10^-44 */
78 //! Rank in the multisimulation
79 #define MSRANK(ms, nodeid) (nodeid)
81 //! Enum for replica exchange flavours
82 enum {
83 ereTEMP, ereLAMBDA, ereENDSINGLE, ereTL, ereNR
85 /*! \brief Strings describing replica exchange flavours.
87 * end_single_marker merely notes the end of single variable replica
88 * exchange. All types higher than it are multiple replica exchange
89 * methods.
91 * Eventually, should add 'pressure', 'temperature and pressure',
92 * 'lambda_and_pressure', 'temperature_lambda_pressure'?; Let's wait
93 * until we feel better about the pressure control methods giving
94 * exact ensembles. Right now, we assume constant pressure */
95 static const char *erename[ereNR] = { "temperature", "lambda", "end_single_marker", "temperature and lambda"};
97 //! Working data for replica exchange.
98 struct gmx_repl_ex
100 //! Replica ID
101 int repl;
102 //! Total number of replica
103 int nrepl;
104 //! Temperature
105 real temp;
106 //! Replica exchange type from ere enum
107 int type;
108 //! Quantity, e.g. temperature or lambda; first index is ere, second index is replica ID
109 real **q;
110 //! Use constant pressure and temperature
111 gmx_bool bNPT;
112 //! Replica pressures
113 real *pres;
114 //! Replica indices
115 int *ind;
116 //! Used for keeping track of all the replica swaps
117 int *allswaps;
118 //! Replica exchange interval (number of steps)
119 int nst;
120 //! Number of exchanges per interval
121 int nex;
122 //! Random seed
123 int seed;
124 //! Number of even and odd replica change attempts
125 int nattempt[2];
126 //! Sum of probabilities
127 real *prob_sum;
128 //! Number of moves between replicas i and j
129 int **nmoves;
130 //! i-th element of the array is the number of exchanges between replica i-1 and i
131 int *nexchange;
133 /*! \brief Helper arrays for replica exchange; allocated here
134 * so they don't have to be allocated each time */
135 //! \{
136 int *destinations;
137 int **cyclic;
138 int **order;
139 int *tmpswap;
140 gmx_bool *incycle;
141 gmx_bool *bEx;
142 //! \}
144 //! Helper arrays to hold the quantities that are exchanged.
145 //! \{
146 real *prob;
147 real *Epot;
148 real *beta;
149 real *Vol;
150 real **de;
151 //! \}
154 // TODO We should add Doxygen here some time.
155 //! \cond
157 static gmx_bool repl_quantity(const gmx_multisim_t *ms,
158 struct gmx_repl_ex *re, int ere, real q)
160 real *qall;
161 gmx_bool bDiff;
162 int s;
164 snew(qall, ms->nsim);
165 qall[re->repl] = q;
166 gmx_sum_sim(ms->nsim, qall, ms);
168 bDiff = FALSE;
169 for (s = 1; s < ms->nsim; s++)
171 if (qall[s] != qall[0])
173 bDiff = TRUE;
177 if (bDiff)
179 /* Set the replica exchange type and quantities */
180 re->type = ere;
182 snew(re->q[ere], re->nrepl);
183 for (s = 0; s < ms->nsim; s++)
185 re->q[ere][s] = qall[s];
188 sfree(qall);
189 return bDiff;
192 gmx_repl_ex_t
193 init_replica_exchange(FILE *fplog,
194 const gmx_multisim_t *ms,
195 int numAtomsInSystem,
196 const t_inputrec *ir,
197 const ReplicaExchangeParameters &replExParams)
199 real pres;
200 int i, j;
201 struct gmx_repl_ex *re;
202 gmx_bool bTemp;
203 gmx_bool bLambda = FALSE;
205 fprintf(fplog, "\nInitializing Replica Exchange\n");
207 if (!isMultiSim(ms) || ms->nsim == 1)
209 gmx_fatal(FARGS, "Nothing to exchange with only one replica, maybe you forgot to set the -multidir option of mdrun?");
211 if (replExParams.numExchanges < 0)
213 gmx_fatal(FARGS, "Replica exchange number of exchanges needs to be positive");
216 if (!EI_DYNAMICS(ir->eI))
218 gmx_fatal(FARGS, "Replica exchange is only supported by dynamical simulations");
219 /* Note that PAR(cr) is defined by cr->nnodes > 1, which is
220 * distinct from isMultiSim(ms). A multi-simulation only runs
221 * with real MPI parallelism, but this does not imply PAR(cr)
222 * is true!
224 * Since we are using a dynamical integrator, the only
225 * decomposition is DD, so PAR(cr) and DOMAINDECOMP(cr) are
226 * synonymous. The only way for cr->nnodes > 1 to be true is
227 * if we are using DD. */
230 snew(re, 1);
232 re->repl = ms->sim;
233 re->nrepl = ms->nsim;
234 snew(re->q, ereENDSINGLE);
236 fprintf(fplog, "Repl There are %d replicas:\n", re->nrepl);
238 /* We only check that the number of atoms in the systms match.
239 * This, of course, do not guarantee that the systems are the same,
240 * but it does guarantee that we can perform replica exchange.
242 check_multi_int(fplog, ms, numAtomsInSystem, "the number of atoms", FALSE);
243 check_multi_int(fplog, ms, ir->eI, "the integrator", FALSE);
244 check_multi_int64(fplog, ms, ir->init_step+ir->nsteps, "init_step+nsteps", FALSE);
245 const int nst = replExParams.exchangeInterval;
246 check_multi_int64(fplog, ms, (ir->init_step+nst-1)/nst,
247 "first exchange step: init_step/-replex", FALSE);
248 check_multi_int(fplog, ms, ir->etc, "the temperature coupling", FALSE);
249 check_multi_int(fplog, ms, ir->opts.ngtc,
250 "the number of temperature coupling groups", FALSE);
251 check_multi_int(fplog, ms, ir->epc, "the pressure coupling", FALSE);
252 check_multi_int(fplog, ms, ir->efep, "free energy", FALSE);
253 check_multi_int(fplog, ms, ir->fepvals->n_lambda, "number of lambda states", FALSE);
255 re->temp = ir->opts.ref_t[0];
256 for (i = 1; (i < ir->opts.ngtc); i++)
258 if (ir->opts.ref_t[i] != re->temp)
260 fprintf(fplog, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
261 fprintf(stderr, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
265 re->type = -1;
266 bTemp = repl_quantity(ms, re, ereTEMP, re->temp);
267 if (ir->efep != efepNO)
269 bLambda = repl_quantity(ms, re, ereLAMBDA, static_cast<real>(ir->fepvals->init_fep_state));
271 if (re->type == -1) /* nothing was assigned */
273 gmx_fatal(FARGS, "The properties of the %d systems are all the same, there is nothing to exchange", re->nrepl);
275 if (bLambda && bTemp)
277 re->type = ereTL;
280 if (bTemp)
282 please_cite(fplog, "Sugita1999a");
283 if (ir->epc != epcNO)
285 re->bNPT = TRUE;
286 fprintf(fplog, "Repl Using Constant Pressure REMD.\n");
287 please_cite(fplog, "Okabe2001a");
289 if (ir->etc == etcBERENDSEN)
291 gmx_fatal(FARGS, "REMD with the %s thermostat does not produce correct potential energy distributions, consider using the %s thermostat instead",
292 ETCOUPLTYPE(ir->etc), ETCOUPLTYPE(etcVRESCALE));
295 if (bLambda)
297 if (ir->fepvals->delta_lambda != 0) /* check this? */
299 gmx_fatal(FARGS, "delta_lambda is not zero");
302 if (re->bNPT)
304 snew(re->pres, re->nrepl);
305 if (ir->epct == epctSURFACETENSION)
307 pres = ir->ref_p[ZZ][ZZ];
309 else
311 pres = 0;
312 j = 0;
313 for (i = 0; i < DIM; i++)
315 if (ir->compress[i][i] != 0)
317 pres += ir->ref_p[i][i];
318 j++;
321 pres /= j;
323 re->pres[re->repl] = pres;
324 gmx_sum_sim(re->nrepl, re->pres, ms);
327 /* Make an index for increasing replica order */
328 /* only makes sense if one or the other is varying, not both!
329 if both are varying, we trust the order the person gave. */
330 snew(re->ind, re->nrepl);
331 for (i = 0; i < re->nrepl; i++)
333 re->ind[i] = i;
336 if (re->type < ereENDSINGLE)
339 for (i = 0; i < re->nrepl; i++)
341 for (j = i+1; j < re->nrepl; j++)
343 if (re->q[re->type][re->ind[j]] < re->q[re->type][re->ind[i]])
345 /* Unordered replicas are supposed to work, but there
346 * is still an issues somewhere.
347 * Note that at this point still re->ind[i]=i.
349 gmx_fatal(FARGS, "Replicas with indices %d < %d have %ss %g > %g, please order your replicas on increasing %s",
350 i, j,
351 erename[re->type],
352 re->q[re->type][i], re->q[re->type][j],
353 erename[re->type]);
355 else if (re->q[re->type][re->ind[j]] == re->q[re->type][re->ind[i]])
357 gmx_fatal(FARGS, "Two replicas have identical %ss", erename[re->type]);
363 /* keep track of all the swaps, starting with the initial placement. */
364 snew(re->allswaps, re->nrepl);
365 for (i = 0; i < re->nrepl; i++)
367 re->allswaps[i] = re->ind[i];
370 switch (re->type)
372 case ereTEMP:
373 fprintf(fplog, "\nReplica exchange in temperature\n");
374 for (i = 0; i < re->nrepl; i++)
376 fprintf(fplog, " %5.1f", re->q[re->type][re->ind[i]]);
378 fprintf(fplog, "\n");
379 break;
380 case ereLAMBDA:
381 fprintf(fplog, "\nReplica exchange in lambda\n");
382 for (i = 0; i < re->nrepl; i++)
384 fprintf(fplog, " %3d", static_cast<int>(re->q[re->type][re->ind[i]]));
386 fprintf(fplog, "\n");
387 break;
388 case ereTL:
389 fprintf(fplog, "\nReplica exchange in temperature and lambda state\n");
390 for (i = 0; i < re->nrepl; i++)
392 fprintf(fplog, " %5.1f", re->q[ereTEMP][re->ind[i]]);
394 fprintf(fplog, "\n");
395 for (i = 0; i < re->nrepl; i++)
397 fprintf(fplog, " %5d", static_cast<int>(re->q[ereLAMBDA][re->ind[i]]));
399 fprintf(fplog, "\n");
400 break;
401 default:
402 gmx_incons("Unknown replica exchange quantity");
404 if (re->bNPT)
406 fprintf(fplog, "\nRepl p");
407 for (i = 0; i < re->nrepl; i++)
409 fprintf(fplog, " %5.2f", re->pres[re->ind[i]]);
412 for (i = 0; i < re->nrepl; i++)
414 if ((i > 0) && (re->pres[re->ind[i]] < re->pres[re->ind[i-1]]))
416 fprintf(fplog, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
417 fprintf(stderr, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
421 re->nst = nst;
422 if (replExParams.randomSeed == -1)
424 if (isMasterSim(ms))
426 re->seed = static_cast<int>(gmx::makeRandomSeed());
428 else
430 re->seed = 0;
432 gmx_sumi_sim(1, &(re->seed), ms);
434 else
436 re->seed = replExParams.randomSeed;
438 fprintf(fplog, "\nReplica exchange interval: %d\n", re->nst);
439 fprintf(fplog, "\nReplica random seed: %d\n", re->seed);
441 re->nattempt[0] = 0;
442 re->nattempt[1] = 0;
444 snew(re->prob_sum, re->nrepl);
445 snew(re->nexchange, re->nrepl);
446 snew(re->nmoves, re->nrepl);
447 for (i = 0; i < re->nrepl; i++)
449 snew(re->nmoves[i], re->nrepl);
451 fprintf(fplog, "Replica exchange information below: ex and x = exchange, pr = probability\n");
453 /* generate space for the helper functions so we don't have to snew each time */
455 snew(re->destinations, re->nrepl);
456 snew(re->incycle, re->nrepl);
457 snew(re->tmpswap, re->nrepl);
458 snew(re->cyclic, re->nrepl);
459 snew(re->order, re->nrepl);
460 for (i = 0; i < re->nrepl; i++)
462 snew(re->cyclic[i], re->nrepl+1);
463 snew(re->order[i], re->nrepl);
465 /* allocate space for the functions storing the data for the replicas */
466 /* not all of these arrays needed in all cases, but they don't take
467 up much space, since the max size is nrepl**2 */
468 snew(re->prob, re->nrepl);
469 snew(re->bEx, re->nrepl);
470 snew(re->beta, re->nrepl);
471 snew(re->Vol, re->nrepl);
472 snew(re->Epot, re->nrepl);
473 snew(re->de, re->nrepl);
474 for (i = 0; i < re->nrepl; i++)
476 snew(re->de[i], re->nrepl);
478 re->nex = replExParams.numExchanges;
479 return re;
482 static void exchange_reals(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, real *v, int n)
484 real *buf;
485 int i;
487 if (v)
489 snew(buf, n);
490 #if GMX_MPI
492 MPI_Sendrecv(v, n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
493 buf,n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
494 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
497 MPI_Request mpi_req;
499 MPI_Isend(v, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
500 ms->mpi_comm_masters, &mpi_req);
501 MPI_Recv(buf, n*sizeof(real), MPI_BYTE, MSRANK(ms, b), 0,
502 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
503 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
505 #endif
506 for (i = 0; i < n; i++)
508 v[i] = buf[i];
510 sfree(buf);
515 static void exchange_doubles(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, double *v, int n)
517 double *buf;
518 int i;
520 if (v)
522 snew(buf, n);
523 #if GMX_MPI
525 MPI_Sendrecv(v, n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
526 buf,n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
527 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
530 MPI_Request mpi_req;
532 MPI_Isend(v, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
533 ms->mpi_comm_masters, &mpi_req);
534 MPI_Recv(buf, n*sizeof(double), MPI_BYTE, MSRANK(ms, b), 0,
535 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
536 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
538 #endif
539 for (i = 0; i < n; i++)
541 v[i] = buf[i];
543 sfree(buf);
547 static void exchange_rvecs(const gmx_multisim_t gmx_unused *ms, int gmx_unused b, rvec *v, int n)
549 rvec *buf;
550 int i;
552 if (v)
554 snew(buf, n);
555 #if GMX_MPI
557 MPI_Sendrecv(v[0], n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
558 buf[0],n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
559 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
562 MPI_Request mpi_req;
564 MPI_Isend(v[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
565 ms->mpi_comm_masters, &mpi_req);
566 MPI_Recv(buf[0], n*sizeof(rvec), MPI_BYTE, MSRANK(ms, b), 0,
567 ms->mpi_comm_masters, MPI_STATUS_IGNORE);
568 MPI_Wait(&mpi_req, MPI_STATUS_IGNORE);
570 #endif
571 for (i = 0; i < n; i++)
573 copy_rvec(buf[i], v[i]);
575 sfree(buf);
579 static void exchange_state(const gmx_multisim_t *ms, int b, t_state *state)
581 /* When t_state changes, this code should be updated. */
582 int ngtc, nnhpres;
583 ngtc = state->ngtc * state->nhchainlength;
584 nnhpres = state->nnhpres* state->nhchainlength;
585 exchange_rvecs(ms, b, state->box, DIM);
586 exchange_rvecs(ms, b, state->box_rel, DIM);
587 exchange_rvecs(ms, b, state->boxv, DIM);
588 exchange_reals(ms, b, &(state->veta), 1);
589 exchange_reals(ms, b, &(state->vol0), 1);
590 exchange_rvecs(ms, b, state->svir_prev, DIM);
591 exchange_rvecs(ms, b, state->fvir_prev, DIM);
592 exchange_rvecs(ms, b, state->pres_prev, DIM);
593 exchange_doubles(ms, b, state->nosehoover_xi.data(), ngtc);
594 exchange_doubles(ms, b, state->nosehoover_vxi.data(), ngtc);
595 exchange_doubles(ms, b, state->nhpres_xi.data(), nnhpres);
596 exchange_doubles(ms, b, state->nhpres_vxi.data(), nnhpres);
597 exchange_doubles(ms, b, state->therm_integral.data(), state->ngtc);
598 exchange_doubles(ms, b, &state->baros_integral, 1);
599 exchange_rvecs(ms, b, state->x.rvec_array(), state->natoms);
600 exchange_rvecs(ms, b, state->v.rvec_array(), state->natoms);
603 static void copy_state_serial(const t_state *src, t_state *dest)
605 if (dest != src)
607 /* Currently the local state is always a pointer to the global
608 * in serial, so we should never end up here.
609 * TODO: Implement a (trivial) t_state copy once converted to C++.
611 GMX_RELEASE_ASSERT(false, "State copying is currently not implemented in replica exchange");
615 static void scale_velocities(gmx::ArrayRef<gmx::RVec> velocities, real fac)
617 for (auto &v : velocities)
619 v *= fac;
623 static void print_transition_matrix(FILE *fplog, int n, int **nmoves, const int *nattempt)
625 int i, j, ntot;
626 float Tprint;
628 ntot = nattempt[0] + nattempt[1];
629 fprintf(fplog, "\n");
630 fprintf(fplog, "Repl");
631 for (i = 0; i < n; i++)
633 fprintf(fplog, " "); /* put the title closer to the center */
635 fprintf(fplog, "Empirical Transition Matrix\n");
637 fprintf(fplog, "Repl");
638 for (i = 0; i < n; i++)
640 fprintf(fplog, "%8d", (i+1));
642 fprintf(fplog, "\n");
644 for (i = 0; i < n; i++)
646 fprintf(fplog, "Repl");
647 for (j = 0; j < n; j++)
649 Tprint = 0.0;
650 if (nmoves[i][j] > 0)
652 Tprint = nmoves[i][j]/(2.0*ntot);
654 fprintf(fplog, "%8.4f", Tprint);
656 fprintf(fplog, "%3d\n", i);
660 static void print_ind(FILE *fplog, const char *leg, int n, int *ind, const gmx_bool *bEx)
662 int i;
664 fprintf(fplog, "Repl %2s %2d", leg, ind[0]);
665 for (i = 1; i < n; i++)
667 fprintf(fplog, " %c %2d", (bEx != nullptr && bEx[i]) ? 'x' : ' ', ind[i]);
669 fprintf(fplog, "\n");
672 static void print_allswitchind(FILE *fplog, int n, int *pind, int *allswaps, int *tmpswap)
674 int i;
676 for (i = 0; i < n; i++)
678 tmpswap[i] = allswaps[i];
680 for (i = 0; i < n; i++)
682 allswaps[i] = tmpswap[pind[i]];
685 fprintf(fplog, "\nAccepted Exchanges: ");
686 for (i = 0; i < n; i++)
688 fprintf(fplog, "%d ", pind[i]);
690 fprintf(fplog, "\n");
692 /* the "Order After Exchange" is the state label corresponding to the configuration that
693 started in state listed in order, i.e.
695 3 0 1 2
697 means that the:
698 configuration starting in simulation 3 is now in simulation 0,
699 configuration starting in simulation 0 is now in simulation 1,
700 configuration starting in simulation 1 is now in simulation 2,
701 configuration starting in simulation 2 is now in simulation 3
703 fprintf(fplog, "Order After Exchange: ");
704 for (i = 0; i < n; i++)
706 fprintf(fplog, "%d ", allswaps[i]);
708 fprintf(fplog, "\n\n");
711 static void print_prob(FILE *fplog, const char *leg, int n, real *prob)
713 int i;
714 char buf[8];
716 fprintf(fplog, "Repl %2s ", leg);
717 for (i = 1; i < n; i++)
719 if (prob[i] >= 0)
721 sprintf(buf, "%4.2f", prob[i]);
722 fprintf(fplog, " %3s", buf[0] == '1' ? "1.0" : buf+1);
724 else
726 fprintf(fplog, " ");
729 fprintf(fplog, "\n");
732 static void print_count(FILE *fplog, const char *leg, int n, int *count)
734 int i;
736 fprintf(fplog, "Repl %2s ", leg);
737 for (i = 1; i < n; i++)
739 fprintf(fplog, " %4d", count[i]);
741 fprintf(fplog, "\n");
744 static real calc_delta(FILE *fplog, gmx_bool bPrint, struct gmx_repl_ex *re, int a, int b, int ap, int bp)
747 real ediff, dpV, delta = 0;
748 real *Epot = re->Epot;
749 real *Vol = re->Vol;
750 real **de = re->de;
751 real *beta = re->beta;
753 /* Two cases; we are permuted and not. In all cases, setting ap = a and bp = b will reduce
754 to the non permuted case */
756 switch (re->type)
758 case ereTEMP:
760 * Okabe et. al. Chem. Phys. Lett. 335 (2001) 435-439
762 ediff = Epot[b] - Epot[a];
763 delta = -(beta[bp] - beta[ap])*ediff;
764 break;
765 case ereLAMBDA:
766 /* two cases: when we are permuted, and not. */
767 /* non-permuted:
768 ediff = E_new - E_old
769 = [H_b(x_a) + H_a(x_b)] - [H_b(x_b) + H_a(x_a)]
770 = [H_b(x_a) - H_a(x_a)] + [H_a(x_b) - H_b(x_b)]
771 = de[b][a] + de[a][b] */
773 /* permuted:
774 ediff = E_new - E_old
775 = [H_bp(x_a) + H_ap(x_b)] - [H_bp(x_b) + H_ap(x_a)]
776 = [H_bp(x_a) - H_ap(x_a)] + [H_ap(x_b) - H_bp(x_b)]
777 = [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)]
778 = [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)]
779 = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]) */
780 /* but, in the current code implementation, we flip configurations, not indices . . .
781 So let's examine that.
782 = [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)]
783 = [H_b(x_ap) - H_a(x_ap)] + [H_a(x_bp) - H_b(x_pb)]
784 = (de[b][ap] - de[a][ap]) + (de[a][bp] - de[b][bp]
785 So, if we exchange b<=> bp and a<=> ap, we return to the same result.
786 So the simple solution is to flip the
787 position of perturbed and original indices in the tests.
790 ediff = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]);
791 delta = ediff*beta[a]; /* assume all same temperature in this case */
792 break;
793 case ereTL:
794 /* not permuted: */
795 /* delta = reduced E_new - reduced E_old
796 = [beta_b H_b(x_a) + beta_a H_a(x_b)] - [beta_b H_b(x_b) + beta_a H_a(x_a)]
797 = [beta_b H_b(x_a) - beta_a H_a(x_a)] + [beta_a H_a(x_b) - beta_b H_b(x_b)]
798 = [beta_b dH_b(x_a) + beta_b H_a(x_a) - beta_a H_a(x_a)] +
799 [beta_a dH_a(x_b) + beta_a H_b(x_b) - beta_b H_b(x_b)]
800 = [beta_b dH_b(x_a) + [beta_a dH_a(x_b) +
801 beta_b (H_a(x_a) - H_b(x_b)]) - beta_a (H_a(x_a) - H_b(x_b))
802 = beta_b dH_b(x_a) + beta_a dH_a(x_b) - (beta_b - beta_a)(H_b(x_b) - H_a(x_a) */
803 /* delta = beta[b]*de[b][a] + beta[a]*de[a][b] - (beta[b] - beta[a])*(Epot[b] - Epot[a]; */
804 /* permuted (big breath!) */
805 /* delta = reduced E_new - reduced E_old
806 = [beta_bp H_bp(x_a) + beta_ap H_ap(x_b)] - [beta_bp H_bp(x_b) + beta_ap H_ap(x_a)]
807 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
808 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
809 - beta_pb H_a(x_a) + beta_ap H_a(x_a) + beta_pb H_a(x_a) - beta_ap H_a(x_a)
810 - beta_ap H_b(x_b) + beta_bp H_b(x_b) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
811 = [(beta_bp H_bp(x_a) - beta_bp H_a(x_a)) - (beta_ap H_ap(x_a) - beta_ap H_a(x_a))] +
812 [(beta_ap H_ap(x_b) - beta_ap H_b(x_b)) - (beta_bp H_bp(x_b) - beta_bp H_b(x_b))]
813 + beta_pb H_a(x_a) - beta_ap H_a(x_a) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
814 = [beta_bp (H_bp(x_a) - H_a(x_a)) - beta_ap (H_ap(x_a) - H_a(x_a))] +
815 [beta_ap (H_ap(x_b) - H_b(x_b)) - beta_bp (H_bp(x_b) - H_b(x_b))]
816 + beta_pb (H_a(x_a) - H_b(x_b)) - beta_ap (H_a(x_a) - H_b(x_b))
817 = ([beta_bp de[bp][a] - beta_ap de[ap][a]) + beta_ap de[ap][b] - beta_bp de[bp][b])
818 + (beta_pb-beta_ap)(H_a(x_a) - H_b(x_b)) */
819 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]);
820 break;
821 default:
822 gmx_incons("Unknown replica exchange quantity");
824 if (bPrint)
826 fprintf(fplog, "Repl %d <-> %d dE_term = %10.3e (kT)\n", a, b, delta);
828 if (re->bNPT)
830 /* revist the calculation for 5.0. Might be some improvements. */
831 dpV = (beta[ap]*re->pres[ap]-beta[bp]*re->pres[bp])*(Vol[b]-Vol[a])/PRESFAC;
832 if (bPrint)
834 fprintf(fplog, " dpV = %10.3e d = %10.3e\n", dpV, delta + dpV);
836 delta += dpV;
838 return delta;
841 static void
842 test_for_replica_exchange(FILE *fplog,
843 const gmx_multisim_t *ms,
844 struct gmx_repl_ex *re,
845 const gmx_enerdata_t *enerd,
846 real vol,
847 int64_t step,
848 real time)
850 int m, i, j, a, b, ap, bp, i0, i1, tmp;
851 real delta = 0;
852 gmx_bool bPrint, bMultiEx;
853 gmx_bool *bEx = re->bEx;
854 real *prob = re->prob;
855 int *pind = re->destinations; /* permuted index */
856 gmx_bool bEpot = FALSE;
857 gmx_bool bDLambda = FALSE;
858 gmx_bool bVol = FALSE;
859 gmx::ThreeFry2x64<64> rng(re->seed, gmx::RandomDomain::ReplicaExchange);
860 gmx::UniformRealDistribution<real> uniformRealDist;
861 gmx::UniformIntDistribution<int> uniformNreplDist(0, re->nrepl-1);
863 bMultiEx = (re->nex > 1); /* multiple exchanges at each state */
864 fprintf(fplog, "Replica exchange at step %" PRId64 " time %.5f\n", step, time);
866 if (re->bNPT)
868 for (i = 0; i < re->nrepl; i++)
870 re->Vol[i] = 0;
872 bVol = TRUE;
873 re->Vol[re->repl] = vol;
875 if ((re->type == ereTEMP || re->type == ereTL))
877 for (i = 0; i < re->nrepl; i++)
879 re->Epot[i] = 0;
881 bEpot = TRUE;
882 re->Epot[re->repl] = enerd->term[F_EPOT];
883 /* temperatures of different states*/
884 for (i = 0; i < re->nrepl; i++)
886 re->beta[i] = 1.0/(re->q[ereTEMP][i]*BOLTZ);
889 else
891 for (i = 0; i < re->nrepl; i++)
893 re->beta[i] = 1.0/(re->temp*BOLTZ); /* we have a single temperature */
896 if (re->type == ereLAMBDA || re->type == ereTL)
898 bDLambda = TRUE;
899 /* lambda differences. */
900 /* de[i][j] is the energy of the jth simulation in the ith Hamiltonian
901 minus the energy of the jth simulation in the jth Hamiltonian */
902 for (i = 0; i < re->nrepl; i++)
904 for (j = 0; j < re->nrepl; j++)
906 re->de[i][j] = 0;
909 for (i = 0; i < re->nrepl; i++)
911 re->de[i][re->repl] = (enerd->enerpart_lambda[static_cast<int>(re->q[ereLAMBDA][i])+1]-enerd->enerpart_lambda[0]);
915 /* now actually do the communication */
916 if (bVol)
918 gmx_sum_sim(re->nrepl, re->Vol, ms);
920 if (bEpot)
922 gmx_sum_sim(re->nrepl, re->Epot, ms);
924 if (bDLambda)
926 for (i = 0; i < re->nrepl; i++)
928 gmx_sum_sim(re->nrepl, re->de[i], ms);
932 /* make a duplicate set of indices for shuffling */
933 for (i = 0; i < re->nrepl; i++)
935 pind[i] = re->ind[i];
938 rng.restart( step, 0 );
940 if (bMultiEx)
942 /* multiple random switch exchange */
943 int nself = 0;
946 for (i = 0; i < re->nex + nself; i++)
948 // For now this is superfluous, but just in case we ever add more
949 // calls in different branches it is safer to always reset the distribution.
950 uniformNreplDist.reset();
952 /* randomly select a pair */
953 /* in theory, could reduce this by identifying only which switches had a nonneglibible
954 probability of occurring (log p > -100) and only operate on those switches */
955 /* find out which state it is from, and what label that state currently has. Likely
956 more work that useful. */
957 i0 = uniformNreplDist(rng);
958 i1 = uniformNreplDist(rng);
959 if (i0 == i1)
961 nself++;
962 continue; /* self-exchange, back up and do it again */
965 a = re->ind[i0]; /* what are the indices of these states? */
966 b = re->ind[i1];
967 ap = pind[i0];
968 bp = pind[i1];
970 bPrint = FALSE; /* too noisy */
971 /* calculate the energy difference */
972 /* if the code changes to flip the STATES, rather than the configurations,
973 use the commented version of the code */
974 /* delta = calc_delta(fplog,bPrint,re,a,b,ap,bp); */
975 delta = calc_delta(fplog, bPrint, re, ap, bp, a, b);
977 /* we actually only use the first space in the prob and bEx array,
978 since there are actually many switches between pairs. */
980 if (delta <= 0)
982 /* accepted */
983 prob[0] = 1;
984 bEx[0] = TRUE;
986 else
988 if (delta > c_probabilityCutoff)
990 prob[0] = 0;
992 else
994 prob[0] = exp(-delta);
996 // roll a number to determine if accepted. For now it is superfluous to
997 // reset, but just in case we ever add more calls in different branches
998 // it is safer to always reset the distribution.
999 uniformRealDist.reset();
1000 bEx[0] = uniformRealDist(rng) < prob[0];
1002 re->prob_sum[0] += prob[0];
1004 if (bEx[0])
1006 /* swap the states */
1007 tmp = pind[i0];
1008 pind[i0] = pind[i1];
1009 pind[i1] = tmp;
1012 re->nattempt[0]++; /* keep track of total permutation trials here */
1013 print_allswitchind(fplog, re->nrepl, pind, re->allswaps, re->tmpswap);
1015 else
1017 /* standard nearest neighbor replica exchange */
1019 m = (step / re->nst) % 2;
1020 for (i = 1; i < re->nrepl; i++)
1022 a = re->ind[i-1];
1023 b = re->ind[i];
1025 bPrint = (re->repl == a || re->repl == b);
1026 if (i % 2 == m)
1028 delta = calc_delta(fplog, bPrint, re, a, b, a, b);
1029 if (delta <= 0)
1031 /* accepted */
1032 prob[i] = 1;
1033 bEx[i] = TRUE;
1035 else
1037 if (delta > c_probabilityCutoff)
1039 prob[i] = 0;
1041 else
1043 prob[i] = exp(-delta);
1045 // roll a number to determine if accepted. For now it is superfluous to
1046 // reset, but just in case we ever add more calls in different branches
1047 // it is safer to always reset the distribution.
1048 uniformRealDist.reset();
1049 bEx[i] = uniformRealDist(rng) < prob[i];
1051 re->prob_sum[i] += prob[i];
1053 if (bEx[i])
1055 /* swap these two */
1056 tmp = pind[i-1];
1057 pind[i-1] = pind[i];
1058 pind[i] = tmp;
1059 re->nexchange[i]++; /* statistics for back compatibility */
1062 else
1064 prob[i] = -1;
1065 bEx[i] = FALSE;
1068 /* print some statistics */
1069 print_ind(fplog, "ex", re->nrepl, re->ind, bEx);
1070 print_prob(fplog, "pr", re->nrepl, prob);
1071 fprintf(fplog, "\n");
1072 re->nattempt[m]++;
1075 /* record which moves were made and accepted */
1076 for (i = 0; i < re->nrepl; i++)
1078 re->nmoves[re->ind[i]][pind[i]] += 1;
1079 re->nmoves[pind[i]][re->ind[i]] += 1;
1081 fflush(fplog); /* make sure we can see what the last exchange was */
1084 static void
1085 cyclic_decomposition(const int *destinations,
1086 int **cyclic,
1087 gmx_bool *incycle,
1088 const int nrepl,
1089 int *nswap)
1092 int i, j, c, p;
1093 int maxlen = 1;
1094 for (i = 0; i < nrepl; i++)
1096 incycle[i] = FALSE;
1098 for (i = 0; i < nrepl; i++) /* one cycle for each replica */
1100 if (incycle[i])
1102 cyclic[i][0] = -1;
1103 continue;
1105 cyclic[i][0] = i;
1106 incycle[i] = TRUE;
1107 c = 1;
1108 p = i;
1109 for (j = 0; j < nrepl; j++) /* potentially all cycles are part, but we will break first */
1111 p = destinations[p]; /* start permuting */
1112 if (p == i)
1114 cyclic[i][c] = -1;
1115 if (c > maxlen)
1117 maxlen = c;
1119 break; /* we've reached the original element, the cycle is complete, and we marked the end. */
1121 else
1123 cyclic[i][c] = p; /* each permutation gives a new member of the cycle */
1124 incycle[p] = TRUE;
1125 c++;
1129 *nswap = maxlen - 1;
1131 if (debug)
1133 for (i = 0; i < nrepl; i++)
1135 fprintf(debug, "Cycle %d:", i);
1136 for (j = 0; j < nrepl; j++)
1138 if (cyclic[i][j] < 0)
1140 break;
1142 fprintf(debug, "%2d", cyclic[i][j]);
1144 fprintf(debug, "\n");
1146 fflush(debug);
1150 static void
1151 compute_exchange_order(int **cyclic,
1152 int **order,
1153 const int nrepl,
1154 const int maxswap)
1156 int i, j;
1158 for (j = 0; j < maxswap; j++)
1160 for (i = 0; i < nrepl; i++)
1162 if (cyclic[i][j+1] >= 0)
1164 order[cyclic[i][j+1]][j] = cyclic[i][j];
1165 order[cyclic[i][j]][j] = cyclic[i][j+1];
1168 for (i = 0; i < nrepl; i++)
1170 if (order[i][j] < 0)
1172 order[i][j] = i; /* if it's not exchanging, it should stay this round*/
1177 if (debug)
1179 fprintf(debug, "Replica Exchange Order\n");
1180 for (i = 0; i < nrepl; i++)
1182 fprintf(debug, "Replica %d:", i);
1183 for (j = 0; j < maxswap; j++)
1185 if (order[i][j] < 0)
1187 break;
1189 fprintf(debug, "%2d", order[i][j]);
1191 fprintf(debug, "\n");
1193 fflush(debug);
1197 static void
1198 prepare_to_do_exchange(struct gmx_repl_ex *re,
1199 const int replica_id,
1200 int *maxswap,
1201 gmx_bool *bThisReplicaExchanged)
1203 int i, j;
1204 /* Hold the cyclic decomposition of the (multiple) replica
1205 * exchange. */
1206 gmx_bool bAnyReplicaExchanged = FALSE;
1207 *bThisReplicaExchanged = FALSE;
1209 for (i = 0; i < re->nrepl; i++)
1211 if (re->destinations[i] != re->ind[i])
1213 /* only mark as exchanged if the index has been shuffled */
1214 bAnyReplicaExchanged = TRUE;
1215 break;
1218 if (bAnyReplicaExchanged)
1220 /* reinitialize the placeholder arrays */
1221 for (i = 0; i < re->nrepl; i++)
1223 for (j = 0; j < re->nrepl; j++)
1225 re->cyclic[i][j] = -1;
1226 re->order[i][j] = -1;
1230 /* Identify the cyclic decomposition of the permutation (very
1231 * fast if neighbor replica exchange). */
1232 cyclic_decomposition(re->destinations, re->cyclic, re->incycle, re->nrepl, maxswap);
1234 /* Now translate the decomposition into a replica exchange
1235 * order at each step. */
1236 compute_exchange_order(re->cyclic, re->order, re->nrepl, *maxswap);
1238 /* Did this replica do any exchange at any point? */
1239 for (j = 0; j < *maxswap; j++)
1241 if (replica_id != re->order[replica_id][j])
1243 *bThisReplicaExchanged = TRUE;
1244 break;
1250 gmx_bool replica_exchange(FILE *fplog, const t_commrec *cr,
1251 const gmx_multisim_t *ms, struct gmx_repl_ex *re,
1252 t_state *state, const gmx_enerdata_t *enerd,
1253 t_state *state_local, int64_t step, real time)
1255 int j;
1256 int replica_id = 0;
1257 int exchange_partner;
1258 int maxswap = 0;
1259 /* Number of rounds of exchanges needed to deal with any multiple
1260 * exchanges. */
1261 /* Where each replica ends up after the exchange attempt(s). */
1262 /* The order in which multiple exchanges will occur. */
1263 gmx_bool bThisReplicaExchanged = FALSE;
1265 if (MASTER(cr))
1267 replica_id = re->repl;
1268 test_for_replica_exchange(fplog, ms, re, enerd, det(state_local->box), step, time);
1269 prepare_to_do_exchange(re, replica_id, &maxswap, &bThisReplicaExchanged);
1271 /* Do intra-simulation broadcast so all processors belonging to
1272 * each simulation know whether they need to participate in
1273 * collecting the state. Otherwise, they might as well get on with
1274 * the next thing to do. */
1275 if (DOMAINDECOMP(cr))
1277 #if GMX_MPI
1278 MPI_Bcast(&bThisReplicaExchanged, sizeof(gmx_bool), MPI_BYTE, MASTERRANK(cr),
1279 cr->mpi_comm_mygroup);
1280 #endif
1283 if (bThisReplicaExchanged)
1285 /* Exchange the states */
1286 /* Collect the global state on the master node */
1287 if (DOMAINDECOMP(cr))
1289 dd_collect_state(cr->dd, state_local, state);
1291 else
1293 copy_state_serial(state_local, state);
1296 if (MASTER(cr))
1298 /* There will be only one swap cycle with standard replica
1299 * exchange, but there may be multiple swap cycles if we
1300 * allow multiple swaps. */
1302 for (j = 0; j < maxswap; j++)
1304 exchange_partner = re->order[replica_id][j];
1306 if (exchange_partner != replica_id)
1308 /* Exchange the global states between the master nodes */
1309 if (debug)
1311 fprintf(debug, "Exchanging %d with %d\n", replica_id, exchange_partner);
1313 exchange_state(ms, exchange_partner, state);
1316 /* For temperature-type replica exchange, we need to scale
1317 * the velocities. */
1318 if (re->type == ereTEMP || re->type == ereTL)
1320 scale_velocities(state->v,
1321 std::sqrt(re->q[ereTEMP][replica_id]/re->q[ereTEMP][re->destinations[replica_id]]));
1326 /* With domain decomposition the global state is distributed later */
1327 if (!DOMAINDECOMP(cr))
1329 /* Copy the global state to the local state data structure */
1330 copy_state_serial(state, state_local);
1334 return bThisReplicaExchanged;
1337 void print_replica_exchange_statistics(FILE *fplog, struct gmx_repl_ex *re)
1339 int i;
1341 fprintf(fplog, "\nReplica exchange statistics\n");
1343 if (re->nex == 0)
1345 fprintf(fplog, "Repl %d attempts, %d odd, %d even\n",
1346 re->nattempt[0]+re->nattempt[1], re->nattempt[1], re->nattempt[0]);
1348 fprintf(fplog, "Repl average probabilities:\n");
1349 for (i = 1; i < re->nrepl; i++)
1351 if (re->nattempt[i%2] == 0)
1353 re->prob[i] = 0;
1355 else
1357 re->prob[i] = re->prob_sum[i]/re->nattempt[i%2];
1360 print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1361 print_prob(fplog, "", re->nrepl, re->prob);
1363 fprintf(fplog, "Repl number of exchanges:\n");
1364 print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1365 print_count(fplog, "", re->nrepl, re->nexchange);
1367 fprintf(fplog, "Repl average number of exchanges:\n");
1368 for (i = 1; i < re->nrepl; i++)
1370 if (re->nattempt[i%2] == 0)
1372 re->prob[i] = 0;
1374 else
1376 re->prob[i] = (static_cast<real>(re->nexchange[i]))/re->nattempt[i%2];
1379 print_ind(fplog, "", re->nrepl, re->ind, nullptr);
1380 print_prob(fplog, "", re->nrepl, re->prob);
1382 fprintf(fplog, "\n");
1384 /* print the transition matrix */
1385 print_transition_matrix(fplog, re->nrepl, re->nmoves, re->nattempt);
1388 //! \endcond