Disable PME Mixed mode with FEP
authorAndrey Alekseenko <al42and@gmail.com>
Fri, 15 Oct 2021 14:26:04 +0000 (16:26 +0200)
committerAndrey Alekseenko <al42and@gmail.com>
Mon, 18 Oct 2021 07:33:40 +0000 (09:33 +0200)
Refs #4190

docs/release-notes/2021/2021.4.rst
src/gromacs/ewald/pme.cpp
src/gromacs/ewald/pme.h
src/gromacs/mdrun/runner.cpp
src/gromacs/taskassignment/decidegpuusage.cpp
src/gromacs/taskassignment/decidegpuusage.h

index cdf6d4569aa4967fc3eeebcf0d698b5333535ef4..d162efb241e54daab60cf1261f0a7d7e9347217b 100644 (file)
@@ -25,6 +25,18 @@ a domain and its halo were more than 200000.
 
 :issue:`4167`
 
+Disabled the use of PME Mixed mode for FEP simulations
+""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+The use of Mixed mode PME (``-pme gpu -pmefft cpu``) led to incorrect
+computation of :math:`{\frac{\partial V}{\partial {\lambda}}}` in FEP
+simulations.
+
+Mixed mode is only used when explicitly requested by the user.
+
+:issue:`4190`
+
+
 Fixes for ``gmx`` tools
 ^^^^^^^^^^^^^^^^^^^^^^^
 
index 76051509ecd8d499821b82fd8ed07e7bf290a18d..b0b59e6eedbaf7fe8f37134805eab50a5280b711 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.
@@ -207,6 +207,16 @@ bool pme_gpu_supports_input(const t_inputrec& ir, std::string* error)
     return addMessageIfNotSupported(errorReasons, error);
 }
 
+bool pme_gpu_mixed_mode_supports_input(const t_inputrec& ir, std::string* error)
+{
+    std::list<std::string> errorReasons;
+    if (ir.efep != efepNO)
+    {
+        errorReasons.emplace_back("Free Energy Perturbation (in PME GPU mixed mode)");
+    }
+    return addMessageIfNotSupported(errorReasons, error);
+}
+
 /*! \brief \libinternal
  * Finds out if PME with given inputs is possible to run on GPU.
  * This function is an internal final check, validating the whole PME structure on creation,
index 29411cad0d864694c17c239677887797f93b9f36..3ac892b3ac786583d5e2a06e633f9c6119e7273f 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.
@@ -260,6 +260,17 @@ bool pme_gpu_supports_hardware(const gmx_hw_info_t& hwinfo, std::string* error);
  */
 bool pme_gpu_supports_input(const t_inputrec& ir, std::string* error);
 
+/*! \brief Checks whether the input system allows to run PME on GPU in Mixed mode.
+ * Assumes that the input system is compatible with GPU PME otherwise, that is,
+ * before calling this function one should check that \ref pme_gpu_supports_input returns \c true.
+ *
+ * \param[in]  ir     Input system.
+ * \param[out] error  If non-null, the error message if the input is not supported.
+ *
+ * \returns true if PME can run on GPU in Mixed mode with this input, false otherwise.
+ */
+bool pme_gpu_mixed_mode_supports_input(const t_inputrec& ir, std::string* error);
+
 /*! \brief
  * Returns the active PME codepath (CPU, GPU, mixed).
  * \todo This is a rather static data that should be managed by the higher level task scheduler.
index 8d19c598c2a03a6baf8aa6a3b940af7ea1badd16..423a90905cbb3eb86779f38271287c9ee0df3912 100644 (file)
@@ -815,8 +815,8 @@ int Mdrunner::mdrunner()
                     gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, GMX_THREAD_MPI),
                     hw_opt.nthreads_tmpi);
             useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi(
-                    useGpuForNonbonded, pmeTarget, numDevicesToUse, userGpuTaskAssignment, *hwinfo_,
-                    *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
+                    useGpuForNonbonded, pmeTarget, pmeFftTarget, numDevicesToUse, userGpuTaskAssignment,
+                    *hwinfo_, *inputrec, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
         }
         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
 
@@ -891,8 +891,8 @@ int Mdrunner::mdrunner()
                 nonbondedTarget, userGpuTaskAssignment, emulateGpuNonbonded, canUseGpuForNonbonded,
                 gpuAccelerationOfNonbondedIsUseful(mdlog, *inputrec, !GMX_THREAD_MPI), gpusWereDetected);
         useGpuForPme = decideWhetherToUseGpusForPme(
-                useGpuForNonbonded, pmeTarget, userGpuTaskAssignment, *hwinfo_, *inputrec,
-                cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
+                useGpuForNonbonded, pmeTarget, pmeFftTarget, userGpuTaskAssignment, *hwinfo_,
+                *inputrec, cr->sizeOfDefaultCommunicator, domdecOptions.numPmeRanks, gpusWereDetected);
         useGpuForBonded = decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme,
                                                           bondedTarget, *inputrec, mtop,
                                                           domdecOptions.numPmeRanks, gpusWereDetected);
index 72d8e2a71aa58862801788cb9435c9a8e9999a85..bc9af8d4e09ade7668c959c3345e9723dc00e5e8 100644 (file)
@@ -154,6 +154,7 @@ bool decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget        non
 
 bool decideWhetherToUseGpusForPmeWithThreadMpi(const bool              useGpuForNonbonded,
                                                const TaskTarget        pmeTarget,
+                                               const TaskTarget        pmeFftTarget,
                                                const int               numDevicesToUse,
                                                const std::vector<int>& userGpuTaskAssignment,
                                                const gmx_hw_info_t&    hardwareInfo,
@@ -165,8 +166,15 @@ bool decideWhetherToUseGpusForPmeWithThreadMpi(const bool              useGpuFor
     if ((pmeTarget == TaskTarget::Cpu) || !useGpuForNonbonded || !pme_gpu_supports_build(nullptr)
         || !pme_gpu_supports_hardware(hardwareInfo, nullptr) || !pme_gpu_supports_input(inputrec, nullptr))
     {
-        // PME can't run on a GPU. If the user required that, we issue
-        // an error later.
+        // PME can't run on a GPU. If the user required that, we issue an error later.
+        return false;
+    }
+
+    if (pmeFftTarget == TaskTarget::Cpu && !pme_gpu_mixed_mode_supports_input(inputrec, nullptr))
+    {
+        /* User requested PME FFT on CPU, but the current system is not compatible with Mixed mode,
+         * so we don't use GPUs at all.
+         * If the user had -pme gpu, we issue an error later. */
         return false;
     }
 
@@ -327,6 +335,7 @@ bool decideWhetherToUseGpusForNonbonded(const TaskTarget          nonbondedTarge
 
 bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
                                   const TaskTarget        pmeTarget,
+                                  const TaskTarget        pmeFftTarget,
                                   const std::vector<int>& userGpuTaskAssignment,
                                   const gmx_hw_info_t&    hardwareInfo,
                                   const t_inputrec&       inputrec,
@@ -374,6 +383,16 @@ bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
         }
         return false;
     }
+    if (pmeFftTarget == TaskTarget::Cpu && !pme_gpu_mixed_mode_supports_input(inputrec, &message))
+    {
+        /* User requested PME FFT on CPU, but the current system is not compatible with Mixed mode,
+         * so we don't use GPUs at all. */
+        if (pmeTarget == TaskTarget::Gpu)
+        {
+            GMX_THROW(NotImplementedError("Cannot compute PME interactions in Mixed mode, because " + message));
+        }
+        return false;
+    }
 
     if (pmeTarget == TaskTarget::Cpu)
     {
index 765636c448dec7cb1f69c683270c0bfc1b3322f6..8f16f25ae80b198aabf6cb2116ff741bc976ca32 100644 (file)
@@ -136,6 +136,7 @@ bool decideWhetherToUseGpusForNonbondedWithThreadMpi(TaskTarget              non
  * \param[in]  useGpuForNonbonded        Whether GPUs will be used for nonbonded interactions.
  * \param[in]  pmeTarget                 The user's choice for mdrun -pme for where to assign
  *                                       long-ranged PME nonbonded interaction tasks.
+ * \param[in]  pmeFftTarget              The user's choice for mdrun -pmefft for where to run FFT.
  * \param[in]  numDevicesToUse           The number of compatible GPUs that the user permitted us to use.
  * \param[in]  userGpuTaskAssignment     The user-specified assignment of GPU tasks to device IDs.
  * \param[in]  hardwareInfo              Hardware information
@@ -149,6 +150,7 @@ bool decideWhetherToUseGpusForNonbondedWithThreadMpi(TaskTarget              non
  *             InconsistentInputError  If the user requirements are inconsistent. */
 bool decideWhetherToUseGpusForPmeWithThreadMpi(bool                    useGpuForNonbonded,
                                                TaskTarget              pmeTarget,
+                                               TaskTarget              pmeFftTarget,
                                                int                     numDevicesToUse,
                                                const std::vector<int>& userGpuTaskAssignment,
                                                const gmx_hw_info_t&    hardwareInfo,
@@ -208,6 +210,7 @@ bool decideWhetherToUseGpusForNonbonded(TaskTarget              nonbondedTarget,
  *
  * \param[in]  useGpuForNonbonded        Whether GPUs will be used for nonbonded interactions.
  * \param[in]  pmeTarget                 The user's choice for mdrun -pme for where to assign long-ranged PME nonbonded interaction tasks.
+ * \param[in]  pmeFftTarget              The user's choice for mdrun -pmefft for where to do FFT for PME.
  * \param[in]  userGpuTaskAssignment     The user-specified assignment of GPU tasks to device IDs.
  * \param[in]  hardwareInfo              Hardware information
  * \param[in]  inputrec                  The user input
@@ -221,6 +224,7 @@ bool decideWhetherToUseGpusForNonbonded(TaskTarget              nonbondedTarget,
  *             InconsistentInputError  If the user requirements are inconsistent. */
 bool decideWhetherToUseGpusForPme(bool                    useGpuForNonbonded,
                                   TaskTarget              pmeTarget,
+                                  TaskTarget              pmeFftTarget,
                                   const std::vector<int>& userGpuTaskAssignment,
                                   const gmx_hw_info_t&    hardwareInfo,
                                   const t_inputrec&       inputrec,