Fix multi-domain reruns
authorMark Abraham <mark.j.abraham@gmail.com>
Fri, 27 Jan 2017 15:58:09 +0000 (16:58 +0100)
committerSzilárd Páll <pall.szilard@gmail.com>
Mon, 30 Jan 2017 18:26:50 +0000 (19:26 +0100)
3da127 broke multi-domain rerun by failing to consider the
logic within rerun_parallel_comm.

Fixes #2105

Change-Id: I3b38a248a4fba5b3d46f0275add75948ef7a9c53

src/gromacs/mdlib/md_support.cpp
src/gromacs/mdlib/md_support.h

index 4ad4c7732b722e5dd0a4f687aaec1b603a9ca050..2663e635d1bf4c6bc0f7490ee3b7db1302e9069f 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -626,11 +626,11 @@ int check_nstglobalcomm(FILE *fplog, t_commrec *cr,
 }
 
 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
-                         gmx_bool *bNotLastFrame)
+                         gmx_bool *bLastStep)
 {
     rvec    *xp, *vp;
 
-    if (MASTER(cr) && !*bNotLastFrame)
+    if (MASTER(cr) && *bLastStep)
     {
         fr->natoms = -1;
     }
@@ -640,7 +640,7 @@ void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
     fr->x = xp;
     fr->v = vp;
 
-    *bNotLastFrame = (fr->natoms >= 0);
+    *bLastStep = (fr->natoms < 0);
 
 }
 
index ce594fc874068fdfc21f2aded79ad6a49aae78b9..b2b8b25b8ce5b6e39f3f0ba71c09f4282c1a9429 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -110,7 +110,7 @@ bool multisim_int_all_are_equal(const gmx_multisim_t *ms,
                                 gmx_int64_t           value);
 
 void rerun_parallel_comm(t_commrec *cr, t_trxframe *fr,
-                         gmx_bool *bNotLastFrame);
+                         gmx_bool *bLastStep);
 
 /* get the conserved energy associated with the ensemble type*/
 real compute_conserved_from_auxiliary(t_inputrec *ir, t_state *state,