Fix reading of checkpoints with Parrinello-Rahman (modular simulator)
authorPascal Merz <pascal.merz@me.com>
Thu, 13 Feb 2020 02:32:36 +0000 (19:32 -0700)
committerPascal Merz <pascal.merz@me.com>
Fri, 14 Feb 2020 02:14:41 +0000 (19:14 -0700)
Using modular simulator, simulations using Parrinello-Rahman barostat
could not be read from checkpoint, throwing an error in the checkpoint
loading routine. While the legacy implementation of the P-R barostat
required the pressure at the previous step to be checkpointed, the
modular implementation does not require this. load_checkpoint is,
however, expecting this field to be present and throws an error.

This change fixes this by setting the globalState flags in dependence
of whether the modular simulator will be used, avoiding read_checkpoint
to expect this entry.

Note that tests ensuring this bug not to reappear are introduced in the
child change I3bcd0729.

Refs #3377 (fixes point 3)

Change-Id: If8afd294b8c79ceef66e71293d9d93cf2f7d0df8

docs/release-notes/2020/2020.1.rst
src/gromacs/mdlib/md_support.cpp
src/gromacs/mdlib/md_support.h
src/gromacs/mdrun/runner.cpp

index 366ac52109e2874d34012b9df433afc61f86c360..0c7c3b7f364c8a9de5b79c0b549f64c6cb10d0be 100644 (file)
@@ -64,6 +64,14 @@ mdrun was deadlocking instead of exiting with a failure.
 
 :issue:`3373`
 
+Fix checkpoint restarts using Parrinello-Rahman and md-vv
+"""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+Checkpoints using Parrinello-Rahman and md-vv (only implemented in
+the new modular simulator approach) could not be read.
+
+:issue:`3377`
+
 Fixes for ``gmx`` tools
 ^^^^^^^^^^^^^^^^^^^^^^^
 
index 3f6829005d0262905fd63e38f700fc3735f4483f..8f06530c9d29edd44385741fb779df9dc38683fc 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,2017,2018,2019, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019,2020, 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.
@@ -499,7 +499,7 @@ void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep)
 }
 
 // TODO Most of this logic seems to belong in the respective modules
-void set_state_entries(t_state* state, const t_inputrec* ir)
+void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimulator)
 {
     /* The entries in the state in the tpx file might not correspond
      * with what is needed, so we correct this here.
@@ -529,7 +529,10 @@ void set_state_entries(t_state* state, const t_inputrec* ir)
         if ((ir->epc == epcPARRINELLORAHMAN) || (ir->epc == epcMTTK))
         {
             state->flags |= (1 << estBOXV);
-            state->flags |= (1 << estPRES_PREV);
+            if (!useModularSimulator)
+            {
+                state->flags |= (1 << estPRES_PREV);
+            }
         }
         if (inputrecNptTrotter(ir) || (inputrecNphTrotter(ir)))
         {
index 7023f4ec5128ab3ffa99a2e0893e9d2fad6934cd..2918ed87f7747d82eef7070bac411baf9c10b560 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,2017,2018,2019, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019,2020, 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.
@@ -104,7 +104,7 @@ bool multisim_int_all_are_equal(const gmx_multisim_t* ms, int64_t value);
 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep);
 
 //! \brief Allocate and initialize node-local state entries
-void set_state_entries(t_state* state, const t_inputrec* ir);
+void set_state_entries(t_state* state, const t_inputrec* ir, bool useModularSimulator);
 
 /* Set the lambda values in the global state from a frame read with rerun */
 void setCurrentLambdasRerun(int64_t           step,
index 78e6e4a29e8739697c3c76b73a868adba7e1e71d..f47caae0905ed13fa1a7355fd3db169081b018ea 100644 (file)
@@ -944,7 +944,7 @@ int Mdrunner::mdrunner()
         }
 
         /* now make sure the state is initialized and propagated */
-        set_state_entries(globalState.get(), inputrec);
+        set_state_entries(globalState.get(), inputrec, useModularSimulator);
     }
 
     /* NM and TPI parallelize over force/energy calculations, not atoms,