Respect user input for nstcomm
authorPascal Merz <pascal.merz@me.com>
Thu, 7 Jan 2021 08:39:05 +0000 (08:39 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 7 Jan 2021 08:39:05 +0000 (08:39 +0000)
This replaces do_md overwriting ir->nstcomm by a warning at grompp
time. Until now, nstcomm was set to the computed global communication
frequency if nstcomm was smaller. The user is now informed at grompp
time that their choice is not optimal, but the set value is respected.
The criterion for the warning now also encompasses any value for
nstcomm that is not a multiple of the global communication frequency.

Resolves a TODO, and refs #3854.

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

index cb052c723a681d73c0c4226888d9ac1c680928a5..9072de479c62d948a79ba1395f046d4e41f44584 100644 (file)
@@ -4,7 +4,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 by the GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,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.
@@ -79,6 +79,7 @@
 #include "gromacs/mdlib/calc_verletbuf.h"
 #include "gromacs/mdlib/compute_io.h"
 #include "gromacs/mdlib/constr.h"
+#include "gromacs/mdlib/md_support.h"
 #include "gromacs/mdlib/perf_est.h"
 #include "gromacs/mdlib/vsite.h"
 #include "gromacs/mdrun/mdmodules.h"
@@ -2480,6 +2481,25 @@ int gmx_grompp(int argc, char* argv[])
                 std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
     }
 
+    if (ir->comm_mode != ecmNO)
+    {
+        const int nstglobalcomm = computeGlobalCommunicationPeriod(ir);
+        if (ir->nstcomm % nstglobalcomm != 0)
+        {
+            warning_note(
+                    wi,
+                    gmx::formatString(
+                            "COM removal frequency is set to (%d).\n"
+                            "Other settings require a global communication frequency of %d.\n"
+                            "Note that this will require additional global communication steps,\n"
+                            "which will reduce performance when using multiple ranks.\n"
+                            "Consider setting nstcomm to a multiple of %d.",
+                            ir->nstcomm,
+                            nstglobalcomm,
+                            nstglobalcomm));
+        }
+    }
+
     if (bVerbose)
     {
         GMX_LOG(logger.info).asParagraph().appendTextFormatted("writing run input file...");
index cea710ba18dcba1deb031c85e97ffee03a5a7da3..fc82226bfb6b21312d256369804e862388041fc0 100644 (file)
@@ -4,7 +4,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 The GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,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.
@@ -469,9 +469,9 @@ static int lcd4(int i1, int i2, int i3, int i4)
     return nst;
 }
 
-int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir, const t_commrec* cr)
+int computeGlobalCommunicationPeriod(const t_inputrec* ir)
 {
-    int nstglobalcomm;
+    int nstglobalcomm = 10;
     {
         // Set up the default behaviour
         if (!(ir->nstcalcenergy > 0 || ir->nstlist > 0 || ir->etc != etcNO || ir->epc != epcNO))
@@ -479,7 +479,6 @@ int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir,
             /* The user didn't choose the period for anything
                important, so we just make sure we can send signals and
                write output suitably. */
-            nstglobalcomm = 10;
             if (ir->nstenergy > 0 && ir->nstenergy < nstglobalcomm)
             {
                 nstglobalcomm = ir->nstenergy;
@@ -502,16 +501,12 @@ int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir,
                                  ir->epc != epcNO ? ir->nstpcouple : 0);
         }
     }
+    return nstglobalcomm;
+}
 
-    // TODO change this behaviour. Instead grompp should print
-    // a (performance) note and mdrun should not change ir.
-    if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
-    {
-        GMX_LOG(mdlog.warning)
-                .asParagraph()
-                .appendTextFormatted("WARNING: Changing nstcomm from %d to %d", ir->nstcomm, nstglobalcomm);
-        ir->nstcomm = nstglobalcomm;
-    }
+int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, const t_inputrec* ir, const t_commrec* cr)
+{
+    const int nstglobalcomm = computeGlobalCommunicationPeriod(ir);
 
     if (cr->nnodes > 1)
     {
index 3544c6f14f68f1fef79e838f4f3eee72291fcb74..b6c794b8e68ff5371315f3f43a8c9bdaf2f6368a 100644 (file)
@@ -4,7 +4,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 The GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,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.
@@ -90,11 +90,16 @@ class SimulationSignaller;
  * will be computed, to check none are missing. */
 #define CGLO_CHECK_NUMBER_OF_BONDED_INTERACTIONS (1u << 12u)
 
-
 /*! \brief Return the number of steps that will take place between
  * intra-simulation communications, given the constraints of the
  * inputrec. */
-int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, t_inputrec* ir, const t_commrec* cr);
+int computeGlobalCommunicationPeriod(const t_inputrec* ir);
+
+/*! \brief Return the number of steps that will take place between
+ * intra-simulation communications, given the constraints of the
+ * inputrec, and write information to log.
+ * Calls computeGlobalCommunicationPeriod(ir) internally. */
+int computeGlobalCommunicationPeriod(const gmx::MDLogger& mdlog, const t_inputrec* ir, const t_commrec* cr);
 
 void rerun_parallel_comm(t_commrec* cr, t_trxframe* fr, gmx_bool* bLastStep);