Disabled replica exchange when T not in order
authorErik Lindahl <erik@kth.se>
Wed, 11 Jun 2014 21:56:42 +0000 (23:56 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 17 Jun 2014 23:12:04 +0000 (01:12 +0200)
There are issues when the replica property is not increasing
with the replica index. This patch has a partial fix. But since
there are still issues, unordered replicas now lead to a fatal error.

Fixes #1377.

Change-Id: I272d9831df85fdc18c1b2fe5724cbfc06e73444e

src/kernel/repl_ex.c

index f8046749428ca6f66bb5d7113feeab82ceabede7..0f094d4e267aa65df321b1da7e12e548db2e96c2 100644 (file)
@@ -266,6 +266,16 @@ gmx_repl_ex_t init_replica_exchange(FILE *fplog,
             {
                 if (re->q[re->type][re->ind[j]] < re->q[re->type][re->ind[i]])
                 {
+                    /* Unordered replicas are supposed to work, but there
+                     * is still an issues somewhere.
+                     * Note that at this point still re->ind[i]=i.
+                     */
+                    gmx_fatal(FARGS, "Replicas with indices %d < %d have %ss %g > %g, please order your replicas on increasing %s",
+                              i, j,
+                              erename[re->type],
+                              re->q[re->type][i], re->q[re->type][j],
+                              erename[re->type]);
+
                     k          = re->ind[i];
                     re->ind[i] = re->ind[j];
                     re->ind[j] = k;
@@ -1238,13 +1248,9 @@ compute_exchange_order(FILE     *fplog,
 
 static void
 prepare_to_do_exchange(FILE      *fplog,
-                       const int *destinations,
+                       struct gmx_repl_ex *re,
                        const int  replica_id,
-                       const int  nrepl,
                        int       *maxswap,
-                       int      **order,
-                       int      **cyclic,
-                       int       *incycle,
                        gmx_bool  *bThisReplicaExchanged)
 {
     int i, j;
@@ -1253,9 +1259,9 @@ prepare_to_do_exchange(FILE      *fplog,
     gmx_bool bAnyReplicaExchanged = FALSE;
     *bThisReplicaExchanged = FALSE;
 
-    for (i = 0; i < nrepl; i++)
+    for (i = 0; i < re->nrepl; i++)
     {
-        if (destinations[i] != i)
+        if (re->destinations[i] != re->ind[i])
         {
             /* only mark as exchanged if the index has been shuffled */
             bAnyReplicaExchanged = TRUE;
@@ -1265,27 +1271,27 @@ prepare_to_do_exchange(FILE      *fplog,
     if (bAnyReplicaExchanged)
     {
         /* reinitialize the placeholder arrays */
-        for (i = 0; i < nrepl; i++)
+        for (i = 0; i < re->nrepl; i++)
         {
-            for (j = 0; j < nrepl; j++)
+            for (j = 0; j < re->nrepl; j++)
             {
-                cyclic[i][j] = -1;
-                order[i][j]  = -1;
+                re->cyclic[i][j] = -1;
+                re->order[i][j]  = -1;
             }
         }
 
         /* Identify the cyclic decomposition of the permutation (very
          * fast if neighbor replica exchange). */
-        cyclic_decomposition(fplog, destinations, cyclic, incycle, nrepl, maxswap);
+        cyclic_decomposition(fplog, re->destinations, re->cyclic, re->incycle, re->nrepl, maxswap);
 
         /* Now translate the decomposition into a replica exchange
          * order at each step. */
-        compute_exchange_order(fplog, cyclic, order, nrepl, *maxswap);
+        compute_exchange_order(fplog, re->cyclic, re->order, re->nrepl, *maxswap);
 
         /* Did this replica do any exchange at any point? */
         for (j = 0; j < *maxswap; j++)
         {
-            if (replica_id != order[replica_id][j])
+            if (replica_id != re->order[replica_id][j])
             {
                 *bThisReplicaExchanged = TRUE;
                 break;
@@ -1312,8 +1318,7 @@ gmx_bool replica_exchange(FILE *fplog, const t_commrec *cr, struct gmx_repl_ex *
     {
         replica_id  = re->repl;
         test_for_replica_exchange(fplog, cr->ms, re, enerd, det(state_local->box), step, time);
-        prepare_to_do_exchange(fplog, re->destinations, replica_id, re->nrepl, &maxswap,
-                               re->order, re->cyclic, re->incycle, &bThisReplicaExchanged);
+        prepare_to_do_exchange(fplog, re, replica_id, &maxswap, &bThisReplicaExchanged);
     }
     /* Do intra-simulation broadcast so all processors belonging to
      * each simulation know whether they need to participate in