Ensure restart with update groups always works
authorMark Abraham <mark.j.abraham@gmail.com>
Wed, 18 Aug 2021 08:04:59 +0000 (08:04 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 18 Aug 2021 08:04:59 +0000 (08:04 +0000)
docs/release-notes/2021/2021.3.rst
src/gromacs/domdec/domdec.cpp
src/gromacs/domdec/domdec.h
src/gromacs/mdrun/runner.cpp

index b9ff941e9113f6073935db0d9e51e0156ba4a75a..b132c4acef3218bad83487ef7ffb3a2ad4c981e7 100644 (file)
@@ -44,6 +44,16 @@ https://gitlab.com/gromacs/gromacs/-/tree/master/python_packaging/sample_restrai
 
 :issue:`4078` and :issue:`4102`
 
+Fixed multi-rank restarts from checkpoints written by single-rank simulations
+"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+Currently a single-rank simulation never uses update groups, however a
+multi-rank run can do so. This fix ensures that the atoms within
+update groups always start in the same periodic image, which was not
+guaranteed if the checkpoint was written by a single-rank simulation.
+
+:issue:`4016`
+
 Fixes for ``gmx`` tools
 ^^^^^^^^^^^^^^^^^^^^^^^
 
index 8e27798962354c15ebac3c701defa5ce86da08fd..d75345d5927642325fd1ad0414ba74aaad046e64 100644 (file)
@@ -3221,3 +3221,42 @@ void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces)
         }
     }
 }
+
+void putUpdateGroupAtomsInSamePeriodicImage(const gmx_domdec_t&      dd,
+                                            const gmx_mtop_t&        mtop,
+                                            const matrix             box,
+                                            gmx::ArrayRef<gmx::RVec> positions)
+{
+    int atomOffset = 0;
+    for (const gmx_molblock_t& molblock : mtop.molblock)
+    {
+        const auto& updateGrouping = dd.comm->systemInfo.updateGroupingPerMoleculetype[molblock.type];
+
+        for (int mol = 0; mol < molblock.nmol; mol++)
+        {
+            for (int g = 0; g < updateGrouping.numBlocks(); g++)
+            {
+                const auto& block     = updateGrouping.block(g);
+                const int   atomBegin = atomOffset + block.begin();
+                const int   atomEnd   = atomOffset + block.end();
+                for (int a = atomBegin + 1; a < atomEnd; a++)
+                {
+                    // Make sure that atoms in the same update group
+                    // are in the same periodic image after restarts.
+                    for (int d = DIM - 1; d >= 0; d--)
+                    {
+                        while (positions[a][d] - positions[atomBegin][d] > 0.5_real * box[d][d])
+                        {
+                            positions[a] -= box[d];
+                        }
+                        while (positions[a][d] - positions[atomBegin][d] < -0.5_real * box[d][d])
+                        {
+                            positions[a] += box[d];
+                        }
+                    }
+                }
+            }
+            atomOffset += updateGrouping.fullRange().end();
+        }
+    }
+}
index 2b65e1989b822eb638f3375b86f29fe0d1f2c681..8939ac35e7278211cad903621b4194c2dc15922c 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 2005 - 2014, The GROMACS development team.
  * Copyright (c) 2015,2016,2017,2018,2019 by the GROMACS development team.
- * Copyright (c) 2020, by the GROMACS development team, led by
+ * Copyright (c) 2020,2021, 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.
@@ -345,4 +345,27 @@ void communicateGpuHaloCoordinates(const t_commrec&      cr,
  */
 void communicateGpuHaloForces(const t_commrec& cr, bool accumulateForces);
 
+/*! \brief Wraps the \c positions so that atoms from the same
+ * update group share the same periodic image wrt \c box.
+ *
+ * When DD and update groups are in use, the simulation master rank
+ * should call this to ensure that e.g. when restarting a simulation
+ * that did not use update groups that the coordinates satisfy the new
+ * requirements.
+ *
+ * This function can probably be removed when even single-rank
+ * simulations use domain decomposition, because then the choice of
+ * whether update groups are used is probably going to be the same
+ * regardless of the rank count.
+ *
+ * \param[in]    dd         The DD manager
+ * \param[in]    mtop       The system topology
+ * \param[in]    box        The global system box
+ * \param[in]    positions  The global system positions
+ */
+void putUpdateGroupAtomsInSamePeriodicImage(const gmx_domdec_t&      dd,
+                                            const gmx_mtop_t&        mtop,
+                                            const matrix             box,
+                                            gmx::ArrayRef<gmx::RVec> positions);
+
 #endif
index 2a5a7c2fe88cb93e61e4f4e8ffeb9233ea92f3c3..8d19c598c2a03a6baf8aa6a3b940af7ea1badd16 100644 (file)
@@ -1232,6 +1232,19 @@ int Mdrunner::mdrunner()
         ddBuilder.reset(nullptr);
         // Note that local state still does not exist yet.
     }
+    // Ensure that all atoms within the same update group are in the
+    // same periodic image. Otherwise, a simulation that did not use
+    // update groups (e.g. a single-rank simulation) cannot always be
+    // correctly restarted in a way that does use update groups
+    // (e.g. a multi-rank simulation).
+    if (isSimulationMasterRank)
+    {
+        const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
+        if (useUpdateGroups)
+        {
+            putUpdateGroupAtomsInSamePeriodicImage(*cr->dd, mtop, globalState->box, globalState->x);
+        }
+    }
 
     // The GPU update is decided here because we need to know whether the constraints or
     // SETTLEs can span accross the domain borders (i.e. whether or not update groups are