Limit auto assignment of bondeds to GPUs
authorBerk Hess <hess@kth.se>
Fri, 19 Oct 2018 19:20:40 +0000 (21:20 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 31 Oct 2018 20:00:46 +0000 (21:00 +0100)
When the CPU has no forces to calculatie, we should not assign
bondeds to the GPU. Here we only assign bondeds to the GPU when
the CPU does PME for electrostatics and/or LJ.

Change-Id: I82f3047311bde45c026c78e44adc4985eb07b0e7

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

index ce9b5b8eb3a1601bc7f1ddc4e95016450810f788..57031c7a529ddd0ad41ab4d067c8cfa6ac1dc3af 100644 (file)
@@ -694,7 +694,9 @@ int Mdrunner::mdrunner()
         auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
         useGpuForBonded =
             decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
-                                            bondedTarget, canUseGpuForBonded, cr->nnodes,
+                                            bondedTarget, canUseGpuForBonded,
+                                            EVDW_PME(inputrec->vdwtype),
+                                            EEL_PME_EWALD(inputrec->coulombtype),
                                             domdecOptions.numPmeRanks, gpusWereDetected);
 
         pmeRunMode   = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
index dbd9dd67460e535aa28572c5f06a6d35c61d4fbd..4488f8d5893dc47df9cbda9e294669a9b52c2ec0 100644 (file)
@@ -409,12 +409,13 @@ bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
 }
 
 bool decideWhetherToUseGpusForBonded(const bool       useGpuForNonbonded,
-                                     const bool /*useGpuForPme*/,
+                                     const bool       useGpuForPme,
                                      const bool       usingVerletScheme,
                                      const TaskTarget bondedTarget,
                                      const bool       canUseGpuForBonded,
-                                     const int        /*numRanksPerSimulation*/,
-                                     const int        /*numPmeRanksPerSimulation*/,
+                                     const bool       usingLJPme,
+                                     const bool       usingElecPmeOrEwald,
+                                     const int        numPmeRanksPerSimulation,
                                      const bool       gpusWereDetected)
 {
     if (bondedTarget == TaskTarget::Cpu)
@@ -470,8 +471,14 @@ bool decideWhetherToUseGpusForBonded(const bool       useGpuForNonbonded,
     }
 
     // If we get here, then the user permitted GPUs, which we should
-    // use for bonded interactions if any were detected.
-    return gpusWereDetected;
+    // use for bonded interactions if any were detected and the CPU
+    // is busy, for which we currently only check PME or Ewald.
+    // (It would be better to dynamically assign bondeds based on timings)
+    // Note that here we assume that the auto setting of PME ranks will not
+    // choose seperate PME ranks when nonBonded are assigned to the GPU.
+    bool usingOurCpuForPmeOrEwald = (usingLJPme || (usingElecPmeOrEwald && !useGpuForPme && numPmeRanksPerSimulation <= 0));
+
+    return gpusWereDetected && usingOurCpuForPmeOrEwald;
 }
 
 }  // namespace gmx
index 9b8279cea775c32388623fe454c51bd0987639a8..3f67d36f3dd278e09508bcab74f08474e54b900e 100644 (file)
@@ -193,8 +193,9 @@ bool decideWhetherToUseGpusForPme(bool                    useGpuForNonbonded,
  * \param[in]  usingVerletScheme         Whether the nonbondeds are using the Verlet scheme.
  * \param[in]  bondedTarget              The user's choice for mdrun -bonded for where to assign tasks.
  * \param[in]  canUseGpuForBonded        Whether the bonded interactions can run on a GPU
- * \param[in]  numRanksPerSimulation     The number of ranks in each simulation.
- * \param[in]  numPmeRanksPerSimulation  The number of PME ranks in each simulation.
+ * \param[in]  usingLJPme                Whether Vdw interactions use LJ-PME.
+ * \param[in]  usingElecPmeOrEwald       Whether a PME or Ewald type method is used for electrostatics.
+ * \param[in]  numPmeRanksPerSimulation  The number of PME ranks in each simulation, can be -1 for auto.
  * \param[in]  gpusWereDetected          Whether compatible GPUs were detected on any node.
  *
  * \returns    Whether the simulation will run bondeded tasks on GPUs.
@@ -206,7 +207,8 @@ bool decideWhetherToUseGpusForBonded(bool       useGpuForNonbonded,
                                      bool       usingVerletScheme,
                                      TaskTarget bondedTarget,
                                      bool       canUseGpuForBonded,
-                                     int        numRanksPerSimulation,
+                                     bool       usingLJPme,
+                                     bool       usingElecPmeOrEwald,
                                      int        numPmeRanksPerSimulation,
                                      bool       gpusWereDetected);