Fix mdrun -nb GPU on non-gpu builds
authorPaul Bauer <paul.bauer.q@gmail.com>
Fri, 21 Dec 2018 13:03:03 +0000 (14:03 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Sat, 22 Dec 2018 23:20:17 +0000 (00:20 +0100)
Added function to check if build supports GPU for nonbonded.

Fixes #2812

Change-Id: I7f907d402dcd7c9f8964d42412125c23c0615909

src/gromacs/gpu_utils/gpu_utils.cpp
src/gromacs/gpu_utils/gpu_utils.h
src/gromacs/mdrun/runner.cpp
src/gromacs/taskassignment/decidegpuusage.cpp
src/gromacs/taskassignment/decidegpuusage.h

index 753c5d88484c8d3fe15e184db592688d8d4683d3..905de06b728ce30fa7a2f511ed92f5fab6bb7dbc 100644 (file)
@@ -46,7 +46,9 @@
 #include <cassert>
 
 #include "gromacs/hardware/gpu_hw_info.h"
+#include "gromacs/utility/arrayref.h"
 #include "gromacs/utility/smalloc.h"
+#include "gromacs/utility/stringutil.h"
 
 #if !GMX_GPU
 /*! \brief Set allocation functions used by the GPU host
@@ -95,3 +97,33 @@ const char *getGpuCompatibilityDescription(const gmx_gpu_info_t &gpu_info,
             gpu_detect_res_str[egpuNonexistent] :
             gpu_detect_res_str[gpu_info_get_stat(gpu_info, index)]);
 }
+/*! \brief Help build a descriptive message in \c error if there are
+ * \c errorReasons why nonbondeds on a GPU are not supported.
+ *
+ * \returns Whether the lack of errorReasons indicate there is support. */
+static bool
+addMessageIfNotSupported(gmx::ArrayRef <const std::string> errorReasons,
+                         std::string                      *error)
+{
+    bool isSupported = errorReasons.empty();
+    if (!isSupported && error)
+    {
+        *error  = "Nonbonded interactions cannot run on GPUs: ";
+        *error += joinStrings(errorReasons, "; ") + ".";
+    }
+    return isSupported;
+}
+
+bool buildSupportsNonbondedOnGpu(std::string *error)
+{
+    std::vector<std::string> errorReasons;
+    if (GMX_DOUBLE)
+    {
+        errorReasons.emplace_back("double precision");
+    }
+    if (GMX_GPU == GMX_GPU_NONE)
+    {
+        errorReasons.emplace_back("non-GPU build of GROMACS");
+    }
+    return addMessageIfNotSupported(errorReasons, error);
+}
index 69e8ee410cb2ca205c01049f021f68b00410a10f..88058dba56b38a065d5edd1e9aae7172561136f1 100644 (file)
@@ -249,6 +249,13 @@ void gpu_set_host_malloc_and_free(bool               bUseGpuKernels,
 //! Get status of device with specified index
 int gpu_info_get_stat(const gmx_gpu_info_t &info, int index);
 
+/*! \brief Check if GROMACS has been built with GPU support.
+ *
+ * \param[in] error Pointer to error string or nullptr.
+ * \todo Move this to NB module once it exists.
+ */
+bool buildSupportsNonbondedOnGpu(std::string *error);
+
 /*! \brief Starts the GPU profiler if mdrun is being profiled.
  *
  *  When a profiler run is in progress (based on the presence of the NVPROF_ID
index 060c83f33a562d917c606e4fe8c26529d489320c..110c9e3c3720b31a89c1e9317dd674f2383000cd 100644 (file)
@@ -614,9 +614,10 @@ int Mdrunner::mdrunner()
             // If the user specified the number of ranks, then we must
             // respect that, but in default mode, we need to allow for
             // the number of GPUs to choose the number of ranks.
-
+            auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
             useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
                     (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
+                    canUseGpuForNonbonded,
                     inputrec->cutoff_scheme == ecutsVERLET,
                     gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
                     hw_opt.nthreads_tmpi);
@@ -682,10 +683,13 @@ int Mdrunner::mdrunner()
         // different nodes, which is the user's responsibilty to
         // handle. If unsuitable, we will notice that during task
         // assignment.
-        bool gpusWereDetected  = hwinfo->ngpu_compatible_tot > 0;
-        bool usingVerletScheme = inputrec->cutoff_scheme == ecutsVERLET;
+        bool gpusWereDetected      = hwinfo->ngpu_compatible_tot > 0;
+        bool usingVerletScheme     = inputrec->cutoff_scheme == ecutsVERLET;
+        auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
         useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
-                                                                emulateGpuNonbonded, usingVerletScheme,
+                                                                emulateGpuNonbonded,
+                                                                canUseGpuForNonbonded,
+                                                                usingVerletScheme,
                                                                 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
                                                                 gpusWereDetected);
         auto canUseGpuForPme   = pme_gpu_supports_build(*hwinfo, nullptr) && pme_gpu_supports_input(*inputrec, mtop, nullptr);
index 4488f8d5893dc47df9cbda9e294669a9b52c2ec0..8e74f9b0872cebae8854770b08175478032226d3 100644 (file)
@@ -105,6 +105,7 @@ decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget          nonbon
                                                 const std::vector<int>   &gpuIdsToUse,
                                                 const std::vector<int>   &userGpuTaskAssignment,
                                                 const EmulateGpuNonbonded emulateGpuNonbonded,
+                                                const bool                buildSupportsNonbondedOnGpu,
                                                 const bool                usingVerletScheme,
                                                 const bool                nonbondedOnGpuIsUseful,
                                                 const int                 numRanksPerSimulation)
@@ -113,7 +114,8 @@ decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget          nonbon
     if (nonbondedTarget == TaskTarget::Cpu ||
         emulateGpuNonbonded == EmulateGpuNonbonded::Yes ||
         !usingVerletScheme ||
-        !nonbondedOnGpuIsUseful)
+        !nonbondedOnGpuIsUseful ||
+        !buildSupportsNonbondedOnGpu)
     {
         // If the user required NB on GPUs, we issue an error later.
         return false;
@@ -235,6 +237,7 @@ decideWhetherToUseGpusForPmeWithThreadMpi(const bool              useGpuForNonbo
 bool decideWhetherToUseGpusForNonbonded(const TaskTarget           nonbondedTarget,
                                         const std::vector<int>    &userGpuTaskAssignment,
                                         const EmulateGpuNonbonded  emulateGpuNonbonded,
+                                        const bool                 buildSupportsNonbondedOnGpu,
                                         const bool                 usingVerletScheme,
                                         const bool                 nonbondedOnGpuIsUseful,
                                         const bool                 gpusWereDetected)
@@ -251,6 +254,15 @@ bool decideWhetherToUseGpusForNonbonded(const TaskTarget           nonbondedTarg
         return false;
     }
 
+    if (!buildSupportsNonbondedOnGpu && nonbondedTarget == TaskTarget::Gpu)
+    {
+        GMX_THROW(InconsistentInputError
+                      ("Nonbonded interactions on the GPU were requested with -nb gpu, "
+                      "but the GROMACS binary has been built without GPU support. "
+                      "Either run without selecting GPU options, or recompile GROMACS "
+                      "with GPU support enabled"));
+    }
+
     // TODO refactor all these TaskTarget::Gpu checks into one place?
     // e.g. use a subfunction that handles only the cases where
     // TaskTargets are not Cpu?
index 3f67d36f3dd278e09508bcab74f08474e54b900e..34e61ffa9284bafb42da9d761200548990e33906 100644 (file)
@@ -68,11 +68,12 @@ enum class TaskTarget : int
  * user. So we need to consider this before any automated choice of
  * the number of thread-MPI ranks.
  *
- * \param[in]  nonbondedTarget           The user's choice for mdrun -nb for where to assign short-ranged nonbonded interaction tasks.
- * \param[in]  gpuIdsToUse               The compatible GPUs that the user permitted us to use.
- * \param[in]  userGpuTaskAssignment     The user-specified assignment of GPU tasks to device IDs.
- * \param[in]  emulateGpuNonbonded       Whether we will emulate GPU calculation of nonbonded interactions.
- * \param[in]  usingVerletScheme         Whether the nonbondeds are using the Verlet scheme.
+ * \param[in]  nonbondedTarget             The user's choice for mdrun -nb for where to assign short-ranged nonbonded interaction tasks.
+ * \param[in]  gpuIdsToUse                 The compatible GPUs that the user permitted us to use.
+ * \param[in]  userGpuTaskAssignment       The user-specified assignment of GPU tasks to device IDs.
+ * \param[in]  emulateGpuNonbonded         Whether we will emulate GPU calculation of nonbonded interactions.
+ * \param[in]  buildSupportsNonbondedOnGpu Whether GROMACS was built with GPU support.
+ * \param[in]  usingVerletScheme           Whether the nonbondeds are using the Verlet scheme.
  * \param[in]  nonbondedOnGpuIsUseful    Whether computing nonbonded interactions on a GPU is useful for this calculation.
  * \param[in]  numRanksPerSimulation     The number of ranks in each simulation.
  *
@@ -80,13 +81,14 @@ enum class TaskTarget : int
  *
  * \throws     std::bad_alloc          If out of memory
  *             InconsistentInputError  If the user requirements are inconsistent. */
-bool decideWhetherToUseGpusForNonbondedWithThreadMpi(TaskTarget                nonbondedTarget,
-                                                     const std::vector<int>   &gpuIdsToUse,
-                                                     const std::vector<int>   &userGpuTaskAssignment,
-                                                     EmulateGpuNonbonded       emulateGpuNonbonded,
-                                                     bool                      usingVerletScheme,
-                                                     bool                      nonbondedOnGpuIsUseful,
-                                                     int                       numRanksPerSimulation);
+bool decideWhetherToUseGpusForNonbondedWithThreadMpi(TaskTarget              nonbondedTarget,
+                                                     const std::vector<int> &gpuIdsToUse,
+                                                     const std::vector<int> &userGpuTaskAssignment,
+                                                     EmulateGpuNonbonded     emulateGpuNonbonded,
+                                                     bool                    buildSupportsNonbondedOnGpu,
+                                                     bool                    usingVerletScheme,
+                                                     bool                    nonbondedOnGpuIsUseful,
+                                                     int                     numRanksPerSimulation);
 
 /*! \brief Decide whether this thread-MPI simulation will run
  * PME tasks on GPUs.
@@ -132,23 +134,25 @@ bool decideWhetherToUseGpusForPmeWithThreadMpi(bool                    useGpuFor
  * decision is made in this routine, along with many more
  * consistency checks.
  *
- * \param[in]  nonbondedTarget           The user's choice for mdrun -nb for where to assign short-ranged nonbonded interaction tasks.
- * \param[in]  userGpuTaskAssignment     The user-specified assignment of GPU tasks to device IDs.
- * \param[in]  emulateGpuNonbonded       Whether we will emulate GPU calculation of nonbonded interactions.
- * \param[in]  usingVerletScheme         Whether the nonbondeds are using the Verlet scheme.
- * \param[in]  nonbondedOnGpuIsUseful    Whether computing nonbonded interactions on a GPU is useful for this calculation.
- * \param[in]  gpusWereDetected          Whether compatible GPUs were detected on any node.
+ * \param[in]  nonbondedTarget             The user's choice for mdrun -nb for where to assign short-ranged nonbonded interaction tasks.
+ * \param[in]  userGpuTaskAssignment       The user-specified assignment of GPU tasks to device IDs.
+ * \param[in]  emulateGpuNonbonded         Whether we will emulate GPU calculation of nonbonded interactions.
+ * \param[in]  buildSupportsNonbondedOnGpu Whether GROMACS was build with GPU support.
+ * \param[in]  usingVerletScheme           Whether the nonbondeds are using the Verlet scheme.
+ * \param[in]  nonbondedOnGpuIsUseful      Whether computing nonbonded interactions on a GPU is useful for this calculation.
+ * \param[in]  gpusWereDetected            Whether compatible GPUs were detected on any node.
  *
  * \returns    Whether the simulation will run nonbonded and PME tasks, respectively, on GPUs.
  *
  * \throws     std::bad_alloc          If out of memory
  *             InconsistentInputError  If the user requirements are inconsistent. */
-bool decideWhetherToUseGpusForNonbonded(TaskTarget                 nonbondedTarget,
-                                        const std::vector<int>    &userGpuTaskAssignment,
-                                        EmulateGpuNonbonded        emulateGpuNonbonded,
-                                        bool                       usingVerletScheme,
-                                        bool                       nonbondedOnGpuIsUseful,
-                                        bool                       gpusWereDetected);
+bool decideWhetherToUseGpusForNonbonded(TaskTarget              nonbondedTarget,
+                                        const std::vector<int> &userGpuTaskAssignment,
+                                        EmulateGpuNonbonded     emulateGpuNonbonded,
+                                        bool                    buildSupportsNonbondedOnGpu,
+                                        bool                    usingVerletScheme,
+                                        bool                    nonbondedOnGpuIsUseful,
+                                        bool                    gpusWereDetected);
 
 /*! \brief Decide whether the simulation will try to run tasks of
  * different types on GPUs.