2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2011,2012,2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gromacs/random/random.h"
46 #include "gromacs/utility/smalloc.h"
47 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
53 #include "gromacs/random/random.h"
55 #define PROBABILITYCUTOFF 100
56 /* we don't bother evaluating if events are more rare than exp(-100) = 3.7x10^-44 */
59 ereTEMP
, ereLAMBDA
, ereENDSINGLE
, ereTL
, ereNR
61 const char *erename
[ereNR
] = { "temperature", "lambda", "end_single_marker", "temperature and lambda"};
62 /* end_single_marker merely notes the end of single variable replica exchange. All types higher than
63 it are multiple replica exchange methods */
64 /* Eventually, should add 'pressure', 'temperature and pressure', 'lambda_and_pressure', 'temperature_lambda_pressure'?;
65 Let's wait until we feel better about the pressure control methods giving exact ensembles. Right now, we assume constant pressure */
67 typedef struct gmx_repl_ex
87 /* these are helper arrays for replica exchange; allocated here so they
88 don't have to be allocated each time */
96 /* helper arrays to hold the quantities that are exchanged */
105 static gmx_bool
repl_quantity(const gmx_multisim_t
*ms
,
106 struct gmx_repl_ex
*re
, int ere
, real q
)
112 snew(qall
, ms
->nsim
);
114 gmx_sum_sim(ms
->nsim
, qall
, ms
);
117 for (s
= 1; s
< ms
->nsim
; s
++)
119 if (qall
[s
] != qall
[0])
127 /* Set the replica exchange type and quantities */
130 snew(re
->q
[ere
], re
->nrepl
);
131 for (s
= 0; s
< ms
->nsim
; s
++)
133 re
->q
[ere
][s
] = qall
[s
];
140 gmx_repl_ex_t
init_replica_exchange(FILE *fplog
,
141 const gmx_multisim_t
*ms
,
142 const t_state
*state
,
143 const t_inputrec
*ir
,
144 int nst
, int nex
, int init_seed
)
148 struct gmx_repl_ex
*re
;
150 gmx_bool bLambda
= FALSE
;
152 fprintf(fplog
, "\nInitializing Replica Exchange\n");
154 if (ms
== NULL
|| ms
->nsim
== 1)
156 gmx_fatal(FARGS
, "Nothing to exchange with only one replica, maybe you forgot to set the -multi option of mdrun?");
158 if (!EI_DYNAMICS(ir
->eI
))
160 gmx_fatal(FARGS
, "Replica exchange is only supported by dynamical simulations");
161 /* Note that PAR(cr) is defined by cr->nnodes > 1, which is
162 * distinct from MULTISIM(cr). A multi-simulation only runs
163 * with real MPI parallelism, but this does not imply PAR(cr)
166 * Since we are using a dynamical integrator, the only
167 * decomposition is DD, so PAR(cr) and DOMAINDECOMP(cr) are
168 * synonymous. The only way for cr->nnodes > 1 to be true is
169 * if we are using DD. */
175 re
->nrepl
= ms
->nsim
;
176 snew(re
->q
, ereENDSINGLE
);
178 fprintf(fplog
, "Repl There are %d replicas:\n", re
->nrepl
);
180 check_multi_int(fplog
, ms
, state
->natoms
, "the number of atoms", FALSE
);
181 check_multi_int(fplog
, ms
, ir
->eI
, "the integrator", FALSE
);
182 check_multi_int64(fplog
, ms
, ir
->init_step
+ir
->nsteps
, "init_step+nsteps", FALSE
);
183 check_multi_int64(fplog
, ms
, (ir
->init_step
+nst
-1)/nst
,
184 "first exchange step: init_step/-replex", FALSE
);
185 check_multi_int(fplog
, ms
, ir
->etc
, "the temperature coupling", FALSE
);
186 check_multi_int(fplog
, ms
, ir
->opts
.ngtc
,
187 "the number of temperature coupling groups", FALSE
);
188 check_multi_int(fplog
, ms
, ir
->epc
, "the pressure coupling", FALSE
);
189 check_multi_int(fplog
, ms
, ir
->efep
, "free energy", FALSE
);
190 check_multi_int(fplog
, ms
, ir
->fepvals
->n_lambda
, "number of lambda states", FALSE
);
192 re
->temp
= ir
->opts
.ref_t
[0];
193 for (i
= 1; (i
< ir
->opts
.ngtc
); i
++)
195 if (ir
->opts
.ref_t
[i
] != re
->temp
)
197 fprintf(fplog
, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
198 fprintf(stderr
, "\nWARNING: The temperatures of the different temperature coupling groups are not identical\n\n");
203 bTemp
= repl_quantity(ms
, re
, ereTEMP
, re
->temp
);
204 if (ir
->efep
!= efepNO
)
206 bLambda
= repl_quantity(ms
, re
, ereLAMBDA
, (real
)ir
->fepvals
->init_fep_state
);
208 if (re
->type
== -1) /* nothing was assigned */
210 gmx_fatal(FARGS
, "The properties of the %d systems are all the same, there is nothing to exchange", re
->nrepl
);
212 if (bLambda
&& bTemp
)
219 please_cite(fplog
, "Sugita1999a");
220 if (ir
->epc
!= epcNO
)
223 fprintf(fplog
, "Repl Using Constant Pressure REMD.\n");
224 please_cite(fplog
, "Okabe2001a");
226 if (ir
->etc
== etcBERENDSEN
)
228 gmx_fatal(FARGS
, "REMD with the %s thermostat does not produce correct potential energy distributions, consider using the %s thermostat instead",
229 ETCOUPLTYPE(ir
->etc
), ETCOUPLTYPE(etcVRESCALE
));
234 if (ir
->fepvals
->delta_lambda
!= 0) /* check this? */
236 gmx_fatal(FARGS
, "delta_lambda is not zero");
241 snew(re
->pres
, re
->nrepl
);
242 if (ir
->epct
== epctSURFACETENSION
)
244 pres
= ir
->ref_p
[ZZ
][ZZ
];
250 for (i
= 0; i
< DIM
; i
++)
252 if (ir
->compress
[i
][i
] != 0)
254 pres
+= ir
->ref_p
[i
][i
];
260 re
->pres
[re
->repl
] = pres
;
261 gmx_sum_sim(re
->nrepl
, re
->pres
, ms
);
264 /* Make an index for increasing replica order */
265 /* only makes sense if one or the other is varying, not both!
266 if both are varying, we trust the order the person gave. */
267 snew(re
->ind
, re
->nrepl
);
268 for (i
= 0; i
< re
->nrepl
; i
++)
273 if (re
->type
< ereENDSINGLE
)
276 for (i
= 0; i
< re
->nrepl
; i
++)
278 for (j
= i
+1; j
< re
->nrepl
; j
++)
280 if (re
->q
[re
->type
][re
->ind
[j
]] < re
->q
[re
->type
][re
->ind
[i
]])
282 /* Unordered replicas are supposed to work, but there
283 * is still an issues somewhere.
284 * Note that at this point still re->ind[i]=i.
286 gmx_fatal(FARGS
, "Replicas with indices %d < %d have %ss %g > %g, please order your replicas on increasing %s",
289 re
->q
[re
->type
][i
], re
->q
[re
->type
][j
],
293 re
->ind
[i
] = re
->ind
[j
];
296 else if (re
->q
[re
->type
][re
->ind
[j
]] == re
->q
[re
->type
][re
->ind
[i
]])
298 gmx_fatal(FARGS
, "Two replicas have identical %ss", erename
[re
->type
]);
304 /* keep track of all the swaps, starting with the initial placement. */
305 snew(re
->allswaps
, re
->nrepl
);
306 for (i
= 0; i
< re
->nrepl
; i
++)
308 re
->allswaps
[i
] = re
->ind
[i
];
314 fprintf(fplog
, "\nReplica exchange in temperature\n");
315 for (i
= 0; i
< re
->nrepl
; i
++)
317 fprintf(fplog
, " %5.1f", re
->q
[re
->type
][re
->ind
[i
]]);
319 fprintf(fplog
, "\n");
322 fprintf(fplog
, "\nReplica exchange in lambda\n");
323 for (i
= 0; i
< re
->nrepl
; i
++)
325 fprintf(fplog
, " %3d", (int)re
->q
[re
->type
][re
->ind
[i
]]);
327 fprintf(fplog
, "\n");
330 fprintf(fplog
, "\nReplica exchange in temperature and lambda state\n");
331 for (i
= 0; i
< re
->nrepl
; i
++)
333 fprintf(fplog
, " %5.1f", re
->q
[ereTEMP
][re
->ind
[i
]]);
335 fprintf(fplog
, "\n");
336 for (i
= 0; i
< re
->nrepl
; i
++)
338 fprintf(fplog
, " %5d", (int)re
->q
[ereLAMBDA
][re
->ind
[i
]]);
340 fprintf(fplog
, "\n");
343 gmx_incons("Unknown replica exchange quantity");
347 fprintf(fplog
, "\nRepl p");
348 for (i
= 0; i
< re
->nrepl
; i
++)
350 fprintf(fplog
, " %5.2f", re
->pres
[re
->ind
[i
]]);
353 for (i
= 0; i
< re
->nrepl
; i
++)
355 if ((i
> 0) && (re
->pres
[re
->ind
[i
]] < re
->pres
[re
->ind
[i
-1]]))
357 fprintf(fplog
, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
358 fprintf(stderr
, "\nWARNING: The reference pressures decrease with increasing temperatures\n\n");
367 re
->seed
= (int)gmx_rng_make_seed();
373 gmx_sumi_sim(1, &(re
->seed
), ms
);
377 re
->seed
= init_seed
;
379 fprintf(fplog
, "\nReplica exchange interval: %d\n", re
->nst
);
380 fprintf(fplog
, "\nReplica random seed: %d\n", re
->seed
);
381 re
->rng
= gmx_rng_init(re
->seed
);
386 snew(re
->prob_sum
, re
->nrepl
);
387 snew(re
->nexchange
, re
->nrepl
);
388 snew(re
->nmoves
, re
->nrepl
);
389 for (i
= 0; i
< re
->nrepl
; i
++)
391 snew(re
->nmoves
[i
], re
->nrepl
);
393 fprintf(fplog
, "Replica exchange information below: x=exchange, pr=probability\n");
395 /* generate space for the helper functions so we don't have to snew each time */
397 snew(re
->destinations
, re
->nrepl
);
398 snew(re
->incycle
, re
->nrepl
);
399 snew(re
->tmpswap
, re
->nrepl
);
400 snew(re
->cyclic
, re
->nrepl
);
401 snew(re
->order
, re
->nrepl
);
402 for (i
= 0; i
< re
->nrepl
; i
++)
404 snew(re
->cyclic
[i
], re
->nrepl
);
405 snew(re
->order
[i
], re
->nrepl
);
407 /* allocate space for the functions storing the data for the replicas */
408 /* not all of these arrays needed in all cases, but they don't take
409 up much space, since the max size is nrepl**2 */
410 snew(re
->prob
, re
->nrepl
);
411 snew(re
->bEx
, re
->nrepl
);
412 snew(re
->beta
, re
->nrepl
);
413 snew(re
->Vol
, re
->nrepl
);
414 snew(re
->Epot
, re
->nrepl
);
415 snew(re
->de
, re
->nrepl
);
416 for (i
= 0; i
< re
->nrepl
; i
++)
418 snew(re
->de
[i
], re
->nrepl
);
424 static void exchange_reals(const gmx_multisim_t gmx_unused
*ms
, int gmx_unused b
, real
*v
, int n
)
434 MPI_Sendrecv(v, n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
435 buf,n*sizeof(real),MPI_BYTE,MSRANK(ms,b),0,
436 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
441 MPI_Isend(v
, n
*sizeof(real
), MPI_BYTE
, MSRANK(ms
, b
), 0,
442 ms
->mpi_comm_masters
, &mpi_req
);
443 MPI_Recv(buf
, n
*sizeof(real
), MPI_BYTE
, MSRANK(ms
, b
), 0,
444 ms
->mpi_comm_masters
, MPI_STATUS_IGNORE
);
445 MPI_Wait(&mpi_req
, MPI_STATUS_IGNORE
);
448 for (i
= 0; i
< n
; i
++)
457 static void exchange_ints(const gmx_multisim_t gmx_unused
*ms
, int gmx_unused b
, int *v
, int n
)
467 MPI_Sendrecv(v, n*sizeof(int),MPI_BYTE,MSRANK(ms,b),0,
468 buf,n*sizeof(int),MPI_BYTE,MSRANK(ms,b),0,
469 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
474 MPI_Isend(v
, n
*sizeof(int), MPI_BYTE
, MSRANK(ms
, b
), 0,
475 ms
->mpi_comm_masters
, &mpi_req
);
476 MPI_Recv(buf
, n
*sizeof(int), MPI_BYTE
, MSRANK(ms
, b
), 0,
477 ms
->mpi_comm_masters
, MPI_STATUS_IGNORE
);
478 MPI_Wait(&mpi_req
, MPI_STATUS_IGNORE
);
481 for (i
= 0; i
< n
; i
++)
489 static void exchange_doubles(const gmx_multisim_t gmx_unused
*ms
, int gmx_unused b
, double *v
, int n
)
499 MPI_Sendrecv(v, n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
500 buf,n*sizeof(double),MPI_BYTE,MSRANK(ms,b),0,
501 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
506 MPI_Isend(v
, n
*sizeof(double), MPI_BYTE
, MSRANK(ms
, b
), 0,
507 ms
->mpi_comm_masters
, &mpi_req
);
508 MPI_Recv(buf
, n
*sizeof(double), MPI_BYTE
, MSRANK(ms
, b
), 0,
509 ms
->mpi_comm_masters
, MPI_STATUS_IGNORE
);
510 MPI_Wait(&mpi_req
, MPI_STATUS_IGNORE
);
513 for (i
= 0; i
< n
; i
++)
521 static void exchange_rvecs(const gmx_multisim_t gmx_unused
*ms
, int gmx_unused b
, rvec
*v
, int n
)
531 MPI_Sendrecv(v[0], n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
532 buf[0],n*sizeof(rvec),MPI_BYTE,MSRANK(ms,b),0,
533 ms->mpi_comm_masters,MPI_STATUS_IGNORE);
538 MPI_Isend(v
[0], n
*sizeof(rvec
), MPI_BYTE
, MSRANK(ms
, b
), 0,
539 ms
->mpi_comm_masters
, &mpi_req
);
540 MPI_Recv(buf
[0], n
*sizeof(rvec
), MPI_BYTE
, MSRANK(ms
, b
), 0,
541 ms
->mpi_comm_masters
, MPI_STATUS_IGNORE
);
542 MPI_Wait(&mpi_req
, MPI_STATUS_IGNORE
);
545 for (i
= 0; i
< n
; i
++)
547 copy_rvec(buf
[i
], v
[i
]);
553 static void exchange_state(const gmx_multisim_t
*ms
, int b
, t_state
*state
)
555 /* When t_state changes, this code should be updated. */
557 ngtc
= state
->ngtc
* state
->nhchainlength
;
558 nnhpres
= state
->nnhpres
* state
->nhchainlength
;
559 exchange_rvecs(ms
, b
, state
->box
, DIM
);
560 exchange_rvecs(ms
, b
, state
->box_rel
, DIM
);
561 exchange_rvecs(ms
, b
, state
->boxv
, DIM
);
562 exchange_reals(ms
, b
, &(state
->veta
), 1);
563 exchange_reals(ms
, b
, &(state
->vol0
), 1);
564 exchange_rvecs(ms
, b
, state
->svir_prev
, DIM
);
565 exchange_rvecs(ms
, b
, state
->fvir_prev
, DIM
);
566 exchange_rvecs(ms
, b
, state
->pres_prev
, DIM
);
567 exchange_doubles(ms
, b
, state
->nosehoover_xi
, ngtc
);
568 exchange_doubles(ms
, b
, state
->nosehoover_vxi
, ngtc
);
569 exchange_doubles(ms
, b
, state
->nhpres_xi
, nnhpres
);
570 exchange_doubles(ms
, b
, state
->nhpres_vxi
, nnhpres
);
571 exchange_doubles(ms
, b
, state
->therm_integral
, state
->ngtc
);
572 exchange_rvecs(ms
, b
, state
->x
, state
->natoms
);
573 exchange_rvecs(ms
, b
, state
->v
, state
->natoms
);
574 exchange_rvecs(ms
, b
, state
->sd_X
, state
->natoms
);
577 static void copy_rvecs(rvec
*s
, rvec
*d
, int n
)
583 for (i
= 0; i
< n
; i
++)
585 copy_rvec(s
[i
], d
[i
]);
590 static void copy_doubles(const double *s
, double *d
, int n
)
596 for (i
= 0; i
< n
; i
++)
603 static void copy_reals(const real
*s
, real
*d
, int n
)
609 for (i
= 0; i
< n
; i
++)
616 static void copy_ints(const int *s
, int *d
, int n
)
622 for (i
= 0; i
< n
; i
++)
629 #define scopy_rvecs(v, n) copy_rvecs(state->v, state_local->v, n);
630 #define scopy_doubles(v, n) copy_doubles(state->v, state_local->v, n);
631 #define scopy_reals(v, n) copy_reals(state->v, state_local->v, n);
632 #define scopy_ints(v, n) copy_ints(state->v, state_local->v, n);
634 static void copy_state_nonatomdata(t_state
*state
, t_state
*state_local
)
636 /* When t_state changes, this code should be updated. */
638 ngtc
= state
->ngtc
* state
->nhchainlength
;
639 nnhpres
= state
->nnhpres
* state
->nhchainlength
;
640 scopy_rvecs(box
, DIM
);
641 scopy_rvecs(box_rel
, DIM
);
642 scopy_rvecs(boxv
, DIM
);
643 state_local
->veta
= state
->veta
;
644 state_local
->vol0
= state
->vol0
;
645 scopy_rvecs(svir_prev
, DIM
);
646 scopy_rvecs(fvir_prev
, DIM
);
647 scopy_rvecs(pres_prev
, DIM
);
648 scopy_doubles(nosehoover_xi
, ngtc
);
649 scopy_doubles(nosehoover_vxi
, ngtc
);
650 scopy_doubles(nhpres_xi
, nnhpres
);
651 scopy_doubles(nhpres_vxi
, nnhpres
);
652 scopy_doubles(therm_integral
, state
->ngtc
);
653 scopy_rvecs(x
, state
->natoms
);
654 scopy_rvecs(v
, state
->natoms
);
655 scopy_rvecs(sd_X
, state
->natoms
);
656 copy_ints(&(state
->fep_state
), &(state_local
->fep_state
), 1);
657 scopy_reals(lambda
, efptNR
);
660 static void scale_velocities(t_state
*state
, real fac
)
666 for (i
= 0; i
< state
->natoms
; i
++)
668 svmul(fac
, state
->v
[i
], state
->v
[i
]);
673 static void print_transition_matrix(FILE *fplog
, int n
, int **nmoves
, int *nattempt
)
678 ntot
= nattempt
[0] + nattempt
[1];
679 fprintf(fplog
, "\n");
680 fprintf(fplog
, "Repl");
681 for (i
= 0; i
< n
; i
++)
683 fprintf(fplog
, " "); /* put the title closer to the center */
685 fprintf(fplog
, "Empirical Transition Matrix\n");
687 fprintf(fplog
, "Repl");
688 for (i
= 0; i
< n
; i
++)
690 fprintf(fplog
, "%8d", (i
+1));
692 fprintf(fplog
, "\n");
694 for (i
= 0; i
< n
; i
++)
696 fprintf(fplog
, "Repl");
697 for (j
= 0; j
< n
; j
++)
700 if (nmoves
[i
][j
] > 0)
702 Tprint
= nmoves
[i
][j
]/(2.0*ntot
);
704 fprintf(fplog
, "%8.4f", Tprint
);
706 fprintf(fplog
, "%3d\n", i
);
710 static void print_ind(FILE *fplog
, const char *leg
, int n
, int *ind
, gmx_bool
*bEx
)
714 fprintf(fplog
, "Repl %2s %2d", leg
, ind
[0]);
715 for (i
= 1; i
< n
; i
++)
717 fprintf(fplog
, " %c %2d", (bEx
!= 0 && bEx
[i
]) ? 'x' : ' ', ind
[i
]);
719 fprintf(fplog
, "\n");
722 static void print_allswitchind(FILE *fplog
, int n
, int *pind
, int *allswaps
, int *tmpswap
)
726 for (i
= 0; i
< n
; i
++)
728 tmpswap
[i
] = allswaps
[i
];
730 for (i
= 0; i
< n
; i
++)
732 allswaps
[i
] = tmpswap
[pind
[i
]];
735 fprintf(fplog
, "\nAccepted Exchanges: ");
736 for (i
= 0; i
< n
; i
++)
738 fprintf(fplog
, "%d ", pind
[i
]);
740 fprintf(fplog
, "\n");
742 /* the "Order After Exchange" is the state label corresponding to the configuration that
743 started in state listed in order, i.e.
748 configuration starting in simulation 3 is now in simulation 0,
749 configuration starting in simulation 0 is now in simulation 1,
750 configuration starting in simulation 1 is now in simulation 2,
751 configuration starting in simulation 2 is now in simulation 3
753 fprintf(fplog
, "Order After Exchange: ");
754 for (i
= 0; i
< n
; i
++)
756 fprintf(fplog
, "%d ", allswaps
[i
]);
758 fprintf(fplog
, "\n\n");
761 static void print_prob(FILE *fplog
, const char *leg
, int n
, real
*prob
)
766 fprintf(fplog
, "Repl %2s ", leg
);
767 for (i
= 1; i
< n
; i
++)
771 sprintf(buf
, "%4.2f", prob
[i
]);
772 fprintf(fplog
, " %3s", buf
[0] == '1' ? "1.0" : buf
+1);
779 fprintf(fplog
, "\n");
782 static void print_count(FILE *fplog
, const char *leg
, int n
, int *count
)
786 fprintf(fplog
, "Repl %2s ", leg
);
787 for (i
= 1; i
< n
; i
++)
789 fprintf(fplog
, " %4d", count
[i
]);
791 fprintf(fplog
, "\n");
794 static real
calc_delta(FILE *fplog
, gmx_bool bPrint
, struct gmx_repl_ex
*re
, int a
, int b
, int ap
, int bp
)
797 real ediff
, dpV
, delta
= 0;
798 real
*Epot
= re
->Epot
;
801 real
*beta
= re
->beta
;
803 /* Two cases; we are permuted and not. In all cases, setting ap = a and bp = b will reduce
804 to the non permuted case */
810 * Okabe et. al. Chem. Phys. Lett. 335 (2001) 435-439
812 ediff
= Epot
[b
] - Epot
[a
];
813 delta
= -(beta
[bp
] - beta
[ap
])*ediff
;
816 /* two cases: when we are permuted, and not. */
818 ediff = E_new - E_old
819 = [H_b(x_a) + H_a(x_b)] - [H_b(x_b) + H_a(x_a)]
820 = [H_b(x_a) - H_a(x_a)] + [H_a(x_b) - H_b(x_b)]
821 = de[b][a] + de[a][b] */
824 ediff = E_new - E_old
825 = [H_bp(x_a) + H_ap(x_b)] - [H_bp(x_b) + H_ap(x_a)]
826 = [H_bp(x_a) - H_ap(x_a)] + [H_ap(x_b) - H_bp(x_b)]
827 = [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)]
828 = [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)]
829 = (de[bp][a] - de[ap][a]) + (de[ap][b] - de[bp][b]) */
830 /* but, in the current code implementation, we flip configurations, not indices . . .
831 So let's examine that.
832 = [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)]
833 = [H_b(x_ap) - H_a(x_ap)] + [H_a(x_bp) - H_b(x_pb)]
834 = (de[b][ap] - de[a][ap]) + (de[a][bp] - de[b][bp]
835 So, if we exchange b<=> bp and a<=> ap, we return to the same result.
836 So the simple solution is to flip the
837 position of perturbed and original indices in the tests.
840 ediff
= (de
[bp
][a
] - de
[ap
][a
]) + (de
[ap
][b
] - de
[bp
][b
]);
841 delta
= ediff
*beta
[a
]; /* assume all same temperature in this case */
845 /* delta = reduced E_new - reduced E_old
846 = [beta_b H_b(x_a) + beta_a H_a(x_b)] - [beta_b H_b(x_b) + beta_a H_a(x_a)]
847 = [beta_b H_b(x_a) - beta_a H_a(x_a)] + [beta_a H_a(x_b) - beta_b H_b(x_b)]
848 = [beta_b dH_b(x_a) + beta_b H_a(x_a) - beta_a H_a(x_a)] +
849 [beta_a dH_a(x_b) + beta_a H_b(x_b) - beta_b H_b(x_b)]
850 = [beta_b dH_b(x_a) + [beta_a dH_a(x_b) +
851 beta_b (H_a(x_a) - H_b(x_b)]) - beta_a (H_a(x_a) - H_b(x_b))
852 = beta_b dH_b(x_a) + beta_a dH_a(x_b) - (beta_b - beta_a)(H_b(x_b) - H_a(x_a) */
853 /* delta = beta[b]*de[b][a] + beta[a]*de[a][b] - (beta[b] - beta[a])*(Epot[b] - Epot[a]; */
854 /* permuted (big breath!) */
855 /* delta = reduced E_new - reduced E_old
856 = [beta_bp H_bp(x_a) + beta_ap H_ap(x_b)] - [beta_bp H_bp(x_b) + beta_ap H_ap(x_a)]
857 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
858 = [beta_bp H_bp(x_a) - beta_ap H_ap(x_a)] + [beta_ap H_ap(x_b) - beta_bp H_bp(x_b)]
859 - beta_pb H_a(x_a) + beta_ap H_a(x_a) + beta_pb H_a(x_a) - beta_ap H_a(x_a)
860 - beta_ap H_b(x_b) + beta_bp H_b(x_b) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
861 = [(beta_bp H_bp(x_a) - beta_bp H_a(x_a)) - (beta_ap H_ap(x_a) - beta_ap H_a(x_a))] +
862 [(beta_ap H_ap(x_b) - beta_ap H_b(x_b)) - (beta_bp H_bp(x_b) - beta_bp H_b(x_b))]
863 + beta_pb H_a(x_a) - beta_ap H_a(x_a) + beta_ap H_b(x_b) - beta_bp H_b(x_b)
864 = [beta_bp (H_bp(x_a) - H_a(x_a)) - beta_ap (H_ap(x_a) - H_a(x_a))] +
865 [beta_ap (H_ap(x_b) - H_b(x_b)) - beta_bp (H_bp(x_b) - H_b(x_b))]
866 + beta_pb (H_a(x_a) - H_b(x_b)) - beta_ap (H_a(x_a) - H_b(x_b))
867 = ([beta_bp de[bp][a] - beta_ap de[ap][a]) + beta_ap de[ap][b] - beta_bp de[bp][b])
868 + (beta_pb-beta_ap)(H_a(x_a) - H_b(x_b)) */
869 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
]);
872 gmx_incons("Unknown replica exchange quantity");
876 fprintf(fplog
, "Repl %d <-> %d dE_term = %10.3e (kT)\n", a
, b
, delta
);
880 /* revist the calculation for 5.0. Might be some improvements. */
881 dpV
= (beta
[ap
]*re
->pres
[ap
]-beta
[bp
]*re
->pres
[bp
])*(Vol
[b
]-Vol
[a
])/PRESFAC
;
884 fprintf(fplog
, " dpV = %10.3e d = %10.3e\n", dpV
, delta
+ dpV
);
892 test_for_replica_exchange(FILE *fplog
,
893 const gmx_multisim_t
*ms
,
894 struct gmx_repl_ex
*re
,
895 gmx_enerdata_t
*enerd
,
900 int m
, i
, j
, a
, b
, ap
, bp
, i0
, i1
, tmp
;
902 gmx_bool bPrint
, bMultiEx
;
903 gmx_bool
*bEx
= re
->bEx
;
904 real
*prob
= re
->prob
;
905 int *pind
= re
->destinations
; /* permuted index */
906 gmx_bool bEpot
= FALSE
;
907 gmx_bool bDLambda
= FALSE
;
908 gmx_bool bVol
= FALSE
;
910 bMultiEx
= (re
->nex
> 1); /* multiple exchanges at each state */
911 fprintf(fplog
, "Replica exchange at step %" GMX_PRId64
" time %.5f\n", step
, time
);
915 for (i
= 0; i
< re
->nrepl
; i
++)
920 re
->Vol
[re
->repl
] = vol
;
922 if ((re
->type
== ereTEMP
|| re
->type
== ereTL
))
924 for (i
= 0; i
< re
->nrepl
; i
++)
929 re
->Epot
[re
->repl
] = enerd
->term
[F_EPOT
];
930 /* temperatures of different states*/
931 for (i
= 0; i
< re
->nrepl
; i
++)
933 re
->beta
[i
] = 1.0/(re
->q
[ereTEMP
][i
]*BOLTZ
);
938 for (i
= 0; i
< re
->nrepl
; i
++)
940 re
->beta
[i
] = 1.0/(re
->temp
*BOLTZ
); /* we have a single temperature */
943 if (re
->type
== ereLAMBDA
|| re
->type
== ereTL
)
946 /* lambda differences. */
947 /* de[i][j] is the energy of the jth simulation in the ith Hamiltonian
948 minus the energy of the jth simulation in the jth Hamiltonian */
949 for (i
= 0; i
< re
->nrepl
; i
++)
951 for (j
= 0; j
< re
->nrepl
; j
++)
956 for (i
= 0; i
< re
->nrepl
; i
++)
958 re
->de
[i
][re
->repl
] = (enerd
->enerpart_lambda
[(int)re
->q
[ereLAMBDA
][i
]+1]-enerd
->enerpart_lambda
[0]);
962 /* now actually do the communication */
965 gmx_sum_sim(re
->nrepl
, re
->Vol
, ms
);
969 gmx_sum_sim(re
->nrepl
, re
->Epot
, ms
);
973 for (i
= 0; i
< re
->nrepl
; i
++)
975 gmx_sum_sim(re
->nrepl
, re
->de
[i
], ms
);
979 /* make a duplicate set of indices for shuffling */
980 for (i
= 0; i
< re
->nrepl
; i
++)
982 pind
[i
] = re
->ind
[i
];
987 /* multiple random switch exchange */
989 for (i
= 0; i
< re
->nex
+ nself
; i
++)
993 gmx_rng_cycle_2uniform(step
, i
*2, re
->seed
, RND_SEED_REPLEX
, rnd
);
994 /* randomly select a pair */
995 /* in theory, could reduce this by identifying only which switches had a nonneglibible
996 probability of occurring (log p > -100) and only operate on those switches */
997 /* find out which state it is from, and what label that state currently has. Likely
998 more work that useful. */
999 i0
= (int)(re
->nrepl
*rnd
[0]);
1000 i1
= (int)(re
->nrepl
*rnd
[1]);
1004 continue; /* self-exchange, back up and do it again */
1007 a
= re
->ind
[i0
]; /* what are the indices of these states? */
1012 bPrint
= FALSE
; /* too noisy */
1013 /* calculate the energy difference */
1014 /* if the code changes to flip the STATES, rather than the configurations,
1015 use the commented version of the code */
1016 /* delta = calc_delta(fplog,bPrint,re,a,b,ap,bp); */
1017 delta
= calc_delta(fplog
, bPrint
, re
, ap
, bp
, a
, b
);
1019 /* we actually only use the first space in the prob and bEx array,
1020 since there are actually many switches between pairs. */
1030 if (delta
> PROBABILITYCUTOFF
)
1036 prob
[0] = exp(-delta
);
1038 /* roll a number to determine if accepted */
1039 gmx_rng_cycle_2uniform(step
, i
*2+1, re
->seed
, RND_SEED_REPLEX
, rnd
);
1040 bEx
[0] = rnd
[0] < prob
[0];
1042 re
->prob_sum
[0] += prob
[0];
1046 /* swap the states */
1048 pind
[i0
] = pind
[i1
];
1052 re
->nattempt
[0]++; /* keep track of total permutation trials here */
1053 print_allswitchind(fplog
, re
->nrepl
, pind
, re
->allswaps
, re
->tmpswap
);
1057 /* standard nearest neighbor replica exchange */
1059 m
= (step
/ re
->nst
) % 2;
1060 for (i
= 1; i
< re
->nrepl
; i
++)
1065 bPrint
= (re
->repl
== a
|| re
->repl
== b
);
1068 delta
= calc_delta(fplog
, bPrint
, re
, a
, b
, a
, b
);
1079 if (delta
> PROBABILITYCUTOFF
)
1085 prob
[i
] = exp(-delta
);
1087 /* roll a number to determine if accepted */
1088 gmx_rng_cycle_2uniform(step
, i
, re
->seed
, RND_SEED_REPLEX
, rnd
);
1089 bEx
[i
] = rnd
[0] < prob
[i
];
1091 re
->prob_sum
[i
] += prob
[i
];
1095 /* swap these two */
1097 pind
[i
-1] = pind
[i
];
1099 re
->nexchange
[i
]++; /* statistics for back compatibility */
1108 /* print some statistics */
1109 print_ind(fplog
, "ex", re
->nrepl
, re
->ind
, bEx
);
1110 print_prob(fplog
, "pr", re
->nrepl
, prob
);
1111 fprintf(fplog
, "\n");
1115 /* record which moves were made and accepted */
1116 for (i
= 0; i
< re
->nrepl
; i
++)
1118 re
->nmoves
[re
->ind
[i
]][pind
[i
]] += 1;
1119 re
->nmoves
[pind
[i
]][re
->ind
[i
]] += 1;
1121 fflush(fplog
); /* make sure we can see what the last exchange was */
1124 static void write_debug_x(t_state
*state
)
1130 for (i
= 0; i
< state
->natoms
; i
+= 10)
1132 fprintf(debug
, "dx %5d %10.5f %10.5f %10.5f\n", i
, state
->x
[i
][XX
], state
->x
[i
][YY
], state
->x
[i
][ZZ
]);
1138 cyclic_decomposition(const int *destinations
,
1147 for (i
= 0; i
< nrepl
; i
++)
1151 for (i
= 0; i
< nrepl
; i
++) /* one cycle for each replica */
1162 for (j
= 0; j
< nrepl
; j
++) /* potentially all cycles are part, but we will break first */
1164 p
= destinations
[p
]; /* start permuting */
1172 break; /* we've reached the original element, the cycle is complete, and we marked the end. */
1176 cyclic
[i
][c
] = p
; /* each permutation gives a new member of the cycle */
1182 *nswap
= maxlen
- 1;
1186 for (i
= 0; i
< nrepl
; i
++)
1188 fprintf(debug
, "Cycle %d:", i
);
1189 for (j
= 0; j
< nrepl
; j
++)
1191 if (cyclic
[i
][j
] < 0)
1195 fprintf(debug
, "%2d", cyclic
[i
][j
]);
1197 fprintf(debug
, "\n");
1204 compute_exchange_order(FILE *fplog
,
1212 for (j
= 0; j
< maxswap
; j
++)
1214 for (i
= 0; i
< nrepl
; i
++)
1216 if (cyclic
[i
][j
+1] >= 0)
1218 order
[cyclic
[i
][j
+1]][j
] = cyclic
[i
][j
];
1219 order
[cyclic
[i
][j
]][j
] = cyclic
[i
][j
+1];
1222 for (i
= 0; i
< nrepl
; i
++)
1224 if (order
[i
][j
] < 0)
1226 order
[i
][j
] = i
; /* if it's not exchanging, it should stay this round*/
1233 fprintf(fplog
, "Replica Exchange Order\n");
1234 for (i
= 0; i
< nrepl
; i
++)
1236 fprintf(fplog
, "Replica %d:", i
);
1237 for (j
= 0; j
< maxswap
; j
++)
1239 if (order
[i
][j
] < 0)
1243 fprintf(debug
, "%2d", order
[i
][j
]);
1245 fprintf(fplog
, "\n");
1252 prepare_to_do_exchange(FILE *fplog
,
1253 struct gmx_repl_ex
*re
,
1254 const int replica_id
,
1256 gmx_bool
*bThisReplicaExchanged
)
1259 /* Hold the cyclic decomposition of the (multiple) replica
1261 gmx_bool bAnyReplicaExchanged
= FALSE
;
1262 *bThisReplicaExchanged
= FALSE
;
1264 for (i
= 0; i
< re
->nrepl
; i
++)
1266 if (re
->destinations
[i
] != re
->ind
[i
])
1268 /* only mark as exchanged if the index has been shuffled */
1269 bAnyReplicaExchanged
= TRUE
;
1273 if (bAnyReplicaExchanged
)
1275 /* reinitialize the placeholder arrays */
1276 for (i
= 0; i
< re
->nrepl
; i
++)
1278 for (j
= 0; j
< re
->nrepl
; j
++)
1280 re
->cyclic
[i
][j
] = -1;
1281 re
->order
[i
][j
] = -1;
1285 /* Identify the cyclic decomposition of the permutation (very
1286 * fast if neighbor replica exchange). */
1287 cyclic_decomposition(re
->destinations
, re
->cyclic
, re
->incycle
, re
->nrepl
, maxswap
);
1289 /* Now translate the decomposition into a replica exchange
1290 * order at each step. */
1291 compute_exchange_order(fplog
, re
->cyclic
, re
->order
, re
->nrepl
, *maxswap
);
1293 /* Did this replica do any exchange at any point? */
1294 for (j
= 0; j
< *maxswap
; j
++)
1296 if (replica_id
!= re
->order
[replica_id
][j
])
1298 *bThisReplicaExchanged
= TRUE
;
1305 gmx_bool
replica_exchange(FILE *fplog
, const t_commrec
*cr
, struct gmx_repl_ex
*re
,
1306 t_state
*state
, gmx_enerdata_t
*enerd
,
1307 t_state
*state_local
, gmx_int64_t step
, real time
)
1311 int exchange_partner
;
1313 /* Number of rounds of exchanges needed to deal with any multiple
1315 /* Where each replica ends up after the exchange attempt(s). */
1316 /* The order in which multiple exchanges will occur. */
1317 gmx_bool bThisReplicaExchanged
= FALSE
;
1321 replica_id
= re
->repl
;
1322 test_for_replica_exchange(fplog
, cr
->ms
, re
, enerd
, det(state_local
->box
), step
, time
);
1323 prepare_to_do_exchange(fplog
, re
, replica_id
, &maxswap
, &bThisReplicaExchanged
);
1325 /* Do intra-simulation broadcast so all processors belonging to
1326 * each simulation know whether they need to participate in
1327 * collecting the state. Otherwise, they might as well get on with
1328 * the next thing to do. */
1329 if (DOMAINDECOMP(cr
))
1332 MPI_Bcast(&bThisReplicaExchanged
, sizeof(gmx_bool
), MPI_BYTE
, MASTERRANK(cr
),
1333 cr
->mpi_comm_mygroup
);
1337 if (bThisReplicaExchanged
)
1339 /* Exchange the states */
1340 /* Collect the global state on the master node */
1341 if (DOMAINDECOMP(cr
))
1343 dd_collect_state(cr
->dd
, state_local
, state
);
1347 copy_state_nonatomdata(state_local
, state
);
1352 /* There will be only one swap cycle with standard replica
1353 * exchange, but there may be multiple swap cycles if we
1354 * allow multiple swaps. */
1356 for (j
= 0; j
< maxswap
; j
++)
1358 exchange_partner
= re
->order
[replica_id
][j
];
1360 if (exchange_partner
!= replica_id
)
1362 /* Exchange the global states between the master nodes */
1365 fprintf(debug
, "Exchanging %d with %d\n", replica_id
, exchange_partner
);
1367 exchange_state(cr
->ms
, exchange_partner
, state
);
1370 /* For temperature-type replica exchange, we need to scale
1371 * the velocities. */
1372 if (re
->type
== ereTEMP
|| re
->type
== ereTL
)
1374 scale_velocities(state
, sqrt(re
->q
[ereTEMP
][replica_id
]/re
->q
[ereTEMP
][re
->destinations
[replica_id
]]));
1379 /* With domain decomposition the global state is distributed later */
1380 if (!DOMAINDECOMP(cr
))
1382 /* Copy the global state to the local state data structure */
1383 copy_state_nonatomdata(state
, state_local
);
1387 return bThisReplicaExchanged
;
1390 void print_replica_exchange_statistics(FILE *fplog
, struct gmx_repl_ex
*re
)
1394 fprintf(fplog
, "\nReplica exchange statistics\n");
1398 fprintf(fplog
, "Repl %d attempts, %d odd, %d even\n",
1399 re
->nattempt
[0]+re
->nattempt
[1], re
->nattempt
[1], re
->nattempt
[0]);
1401 fprintf(fplog
, "Repl average probabilities:\n");
1402 for (i
= 1; i
< re
->nrepl
; i
++)
1404 if (re
->nattempt
[i
%2] == 0)
1410 re
->prob
[i
] = re
->prob_sum
[i
]/re
->nattempt
[i
%2];
1413 print_ind(fplog
, "", re
->nrepl
, re
->ind
, NULL
);
1414 print_prob(fplog
, "", re
->nrepl
, re
->prob
);
1416 fprintf(fplog
, "Repl number of exchanges:\n");
1417 print_ind(fplog
, "", re
->nrepl
, re
->ind
, NULL
);
1418 print_count(fplog
, "", re
->nrepl
, re
->nexchange
);
1420 fprintf(fplog
, "Repl average number of exchanges:\n");
1421 for (i
= 1; i
< re
->nrepl
; i
++)
1423 if (re
->nattempt
[i
%2] == 0)
1429 re
->prob
[i
] = ((real
)re
->nexchange
[i
])/re
->nattempt
[i
%2];
1432 print_ind(fplog
, "", re
->nrepl
, re
->ind
, NULL
);
1433 print_prob(fplog
, "", re
->nrepl
, re
->prob
);
1435 fprintf(fplog
, "\n");
1437 /* print the transition matrix */
1438 print_transition_matrix(fplog
, re
->nrepl
, re
->nmoves
, re
->nattempt
);