Added functionality that allows MdModules to disable PME-only ranks during simulation...
authorDmitry Morozov <aracsmd@gmail.com>
Fri, 29 Jan 2021 10:33:16 +0000 (10:33 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Fri, 29 Jan 2021 10:33:16 +0000 (10:33 +0000)
src/gromacs/ewald/pme.cpp
src/gromacs/ewald/pme.h
src/gromacs/ewald/tests/CMakeLists.txt
src/gromacs/ewald/tests/pme.cpp [new file with mode: 0644]
src/gromacs/mdrun/runner.cpp
src/gromacs/utility/mdmodulenotification.h

index 0daeb036852dc04eabd935dfaeca852f24fa12f3..4fb2b447dc0f75cd8911c120bf9b9013d9ef7e22 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.
@@ -1788,3 +1788,23 @@ bool gmx_pme_grid_matches(const gmx_pme_t& pme, const ivec grid_size)
 {
     return (pme.nkx == grid_size[XX] && pme.nky == grid_size[YY] && pme.nkz == grid_size[ZZ]);
 }
+
+void gmx::SeparatePmeRanksPermitted::disablePmeRanks(const std::string& reason)
+{
+    permitSeparatePmeRanks_ = false;
+
+    if (!reason.empty())
+    {
+        reasons_.push_back(reason);
+    }
+}
+
+bool gmx::SeparatePmeRanksPermitted::permitSeparatePmeRanks() const
+{
+    return permitSeparatePmeRanks_;
+}
+
+std::string gmx::SeparatePmeRanksPermitted::reasonsWhyDisabled() const
+{
+    return joinStrings(reasons_, "; ");
+}
index 29411cad0d864694c17c239677887797f93b9f36..3f2d33f6c69ff8cc969ea01b9aed3682dd3273cd 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.
@@ -50,6 +50,7 @@
 #define GMX_EWALD_PME_H
 
 #include <string>
+#include <vector>
 
 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
 #include "gromacs/gpu_utils/gpu_macros.h"
@@ -83,6 +84,29 @@ class ForceWithVirial;
 class MDLogger;
 enum class PinningPolicy : int;
 class StepWorkload;
+
+/*! \libinternal \brief Class for managing usage of separate PME-only ranks
+ *
+ * Used for checking if some parts of the code could not use PME-only ranks
+ *
+ */
+class SeparatePmeRanksPermitted
+{
+public:
+    //! Disables PME ranks permitted flag with a reason
+    void disablePmeRanks(const std::string& reason);
+    //! Return status of PME ranks usage
+    bool permitSeparatePmeRanks() const;
+    //! Returns all reasons, for not using PME ranks
+    std::string reasonsWhyDisabled() const;
+
+private:
+    //! Flag that informs whether simualtion could use dedicated PME ranks
+    bool permitSeparatePmeRanks_ = true;
+    //! Storage for all reasons, why PME ranks could not be used
+    std::vector<std::string> reasons_;
+};
+
 } // namespace gmx
 
 enum
index 61da015d2dae8308ef455d75098398bcb863373e..ee7eaa7f36433b0b09c3140cd404760373af7a86 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2016,2017,2020, by the GROMACS development team, led by
+# Copyright (c) 2016,2017,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.
@@ -39,4 +39,5 @@ gmx_add_unit_test(EwaldUnitTests ewald-test HARDWARE_DETECTION
         pmesolvetest.cpp
         pmesplinespreadtest.cpp
         pmetestcommon.cpp
+        pme.cpp
 )
diff --git a/src/gromacs/ewald/tests/pme.cpp b/src/gromacs/ewald/tests/pme.cpp
new file mode 100644 (file)
index 0000000..5452cb4
--- /dev/null
@@ -0,0 +1,122 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief
+ * Tests for functionality of the pme related classes:
+ * class SeparatePmeRanksPermitted
+ *
+ * \author Dmitry Morozov <dmitry.morozov@jyu.fi>
+ * \ingroup module_ewald
+ */
+#include "gmxpre.h"
+
+#include <gtest/gtest.h>
+
+#include "gromacs/utility/real.h"
+#include "gromacs/utility/smalloc.h"
+#include "gromacs/utility/stringcompare.h"
+
+#include "testutils/testasserts.h"
+#include "testutils/testmatchers.h"
+
+#include "pmetestcommon.h"
+
+namespace gmx
+{
+
+namespace test
+{
+
+class SeparatePmeRanksPermittedTest : public ::testing::Test
+{
+public:
+    void disableFirstReason() { separatePmeRanksPermitted_.disablePmeRanks("First reason"); }
+
+    void disableSecondReason() { separatePmeRanksPermitted_.disablePmeRanks("Second reason"); }
+
+    void disableEmptyReason() { separatePmeRanksPermitted_.disablePmeRanks(""); }
+
+protected:
+    SeparatePmeRanksPermitted separatePmeRanksPermitted_;
+};
+
+TEST_F(SeparatePmeRanksPermittedTest, ZeroPmeDisableReasons)
+{
+    // Expect that SeparatePmeRanksPermitted is enabled by default
+    EXPECT_TRUE(separatePmeRanksPermitted_.permitSeparatePmeRanks());
+}
+
+TEST_F(SeparatePmeRanksPermittedTest, CanBeDisabled)
+{
+    // Test if disablePmeRanks works
+    EXPECT_NO_THROW(disableFirstReason(););
+}
+
+TEST_F(SeparatePmeRanksPermittedTest, OneDisableReasonFlag)
+{
+    disableFirstReason();
+
+    // Expect that SeparatePmeRanksPermitted is disabled now
+    EXPECT_FALSE(separatePmeRanksPermitted_.permitSeparatePmeRanks());
+}
+
+TEST_F(SeparatePmeRanksPermittedTest, OneDisableReasonText)
+{
+    disableFirstReason();
+
+    // Expect that reasonsWhyDisabled works with one reason
+    EXPECT_TRUE(separatePmeRanksPermitted_.reasonsWhyDisabled() == "First reason");
+}
+
+TEST_F(SeparatePmeRanksPermittedTest, TwoDisableReasonText)
+{
+    disableFirstReason();
+    disableSecondReason();
+
+    // Expect that reasonsWhyDisabled works with two reasons
+    EXPECT_TRUE(separatePmeRanksPermitted_.reasonsWhyDisabled() == "First reason; Second reason");
+}
+
+TEST_F(SeparatePmeRanksPermittedTest, EmptyDisableReasonText)
+{
+    disableEmptyReason();
+
+    // Expect that reasonsWhyDisabled works with empty reason
+    EXPECT_TRUE(separatePmeRanksPermitted_.reasonsWhyDisabled().empty());
+}
+
+} // namespace test
+
+} // namespace gmx
index 943913668c7dbe98a36ea5220df26080b047c622..8ce13dfff4d52f3934b0041b9545605b740c2053 100644 (file)
@@ -64,6 +64,7 @@
 #include "gromacs/domdec/localatomsetmanager.h"
 #include "gromacs/domdec/partition.h"
 #include "gromacs/ewald/ewald_utils.h"
+#include "gromacs/ewald/pme.h"
 #include "gromacs/ewald/pme_gpu_program.h"
 #include "gromacs/ewald/pme_only.h"
 #include "gromacs/ewald/pme_pp_comm_gpu.h"
@@ -1061,40 +1062,59 @@ int Mdrunner::mdrunner()
                   "these are not compatible with mdrun -rerun");
     }
 
+    /* Object for collecting reasons for not using PME-only ranks */
+    SeparatePmeRanksPermitted separatePmeRanksPermitted;
+
+    /* Permit MDModules to notify whether they want to use PME-only ranks */
+    mdModulesNotifier.notify(&separatePmeRanksPermitted);
+
+    /* If simulation is not using PME then disable PME-only ranks */
     if (!(EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype)))
     {
-        if (domdecOptions.numPmeRanks > 0)
-        {
-            gmx_fatal_collective(FARGS,
-                                 cr->mpiDefaultCommunicator,
-                                 MASTER(cr),
-                                 "PME-only ranks are requested, but the system does not use PME "
-                                 "for electrostatics or LJ");
-        }
-
-        domdecOptions.numPmeRanks = 0;
+        separatePmeRanksPermitted.disablePmeRanks(
+                "PME-only ranks are requested, but the system does not use PME "
+                "for electrostatics or LJ");
     }
 
+    /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
+     * improve performance with many threads per GPU, since our OpenMP
+     * scaling is bad, but it's difficult to automate the setup.
+     */
     if (useGpuForNonbonded && domdecOptions.numPmeRanks < 0)
     {
-        /* With NB GPUs we don't automatically use PME-only CPU ranks. PME ranks can
-         * improve performance with many threads per GPU, since our OpenMP
-         * scaling is bad, but it's difficult to automate the setup.
-         */
+        separatePmeRanksPermitted.disablePmeRanks(
+                "PME-only CPU ranks are not automatically used when "
+                "non-bonded interactions are computed on GPUs");
+    }
+
+    /* If GPU is used for PME then only 1 PME rank is permitted */
+    if (useGpuForPme && (domdecOptions.numPmeRanks < 0 || domdecOptions.numPmeRanks > 1))
+    {
+        separatePmeRanksPermitted.disablePmeRanks(
+                "PME GPU decomposition is not supported. Only one separate PME-only GPU rank "
+                "can be used.");
+    }
+
+    /* Disable PME-only ranks if some parts of the code requested so and it's up to GROMACS to decide */
+    if (!separatePmeRanksPermitted.permitSeparatePmeRanks() && domdecOptions.numPmeRanks < 0)
+    {
         domdecOptions.numPmeRanks = 0;
+        GMX_LOG(mdlog.info)
+                .asParagraph()
+                .appendText("Simulation will not use PME-only ranks because: "
+                            + separatePmeRanksPermitted.reasonsWhyDisabled());
     }
-    if (useGpuForPme)
+
+    /* If some parts of the code could not use PME-only ranks and
+     * user explicitly used mdrun -npme option then throw an error */
+    if (!separatePmeRanksPermitted.permitSeparatePmeRanks() && domdecOptions.numPmeRanks > 0)
     {
-        if (domdecOptions.numPmeRanks < 0)
-        {
-            domdecOptions.numPmeRanks = 0;
-            // TODO possibly print a note that one can opt-in for a separate PME GPU rank?
-        }
-        else
-        {
-            GMX_RELEASE_ASSERT(domdecOptions.numPmeRanks <= 1,
-                               "PME GPU decomposition is not supported");
-        }
+        gmx_fatal_collective(FARGS,
+                             cr->mpiDefaultCommunicator,
+                             MASTER(cr),
+                             "Requested -npme %d option is not viable because: %s",
+                             domdecOptions.numPmeRanks,
+                             separatePmeRanksPermitted.reasonsWhyDisabled().c_str());
     }
 
     /* NMR restraints must be initialized before load_checkpoint,
@@ -1458,10 +1478,10 @@ int Mdrunner::mdrunner()
                           hw_opt.nthreads_omp_pme,
                           !thisRankHasDuty(cr, DUTY_PP));
 
-    // Enable FP exception detection, but not in
-    // Release mode and not for compilers with known buggy FP
-    // exception support (clang with any optimization) or suspected
-    // buggy FP exception support (gcc 7.* with optimization).
+// Enable FP exception detection, but not in
+// Release mode and not for compilers with known buggy FP
+// exception support (clang with any optimization) or suspected
+// buggy FP exception support (gcc 7.* with optimization).
 #if !defined NDEBUG                                                                         \
         && !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
              && defined __OPTIMIZE__)
index b5df1e334096f7213d7e146870077140e390b936..8854152fe68347d92ae830fd434bec3f2561f09c 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 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.
@@ -59,6 +59,7 @@ class KeyValueTreeObject;
 class KeyValueTreeObjectBuilder;
 class LocalAtomSetManager;
 class IndexGroupsAndNames;
+class SeparatePmeRanksPermitted;
 struct MdModulesCheckpointReadingDataOnMaster;
 struct MdModulesCheckpointReadingBroadcast;
 struct MdModulesWriteCheckpointData;
@@ -200,6 +201,8 @@ struct MdModulesNotifier
      * MdModulesEnergyOutputToDensityFittingRequestChecker* enables modules to
      *                      report if they want to write their energy output
      *                      to the density fitting field in the energy files
+     * SeparatePmeRanksPermitted* enables modules to report if they want
+     *                      to disable dedicated PME ranks
      * const PbcType& provides modules with the periodic boundary condition type
      *                that is used during the simulation
      * const SimulationTimeStep& provides modules with the simulation time-step
@@ -211,6 +214,7 @@ struct MdModulesNotifier
     registerMdModuleNotification<const KeyValueTreeObject&,
                                  LocalAtomSetManager*,
                                  MdModulesEnergyOutputToDensityFittingRequestChecker*,
+                                 SeparatePmeRanksPermitted*,
                                  const PbcType&,
                                  const SimulationTimeStep&,
                                  const t_commrec&>::type simulationSetupNotifications_;