Add haveEwaldSurfaceContribution to SimulationWorkload
authorBerk Hess <hess@kth.se>
Thu, 31 Oct 2019 19:20:01 +0000 (20:20 +0100)
committerPaul Bauer <paul.bauer.q@gmail.com>
Fri, 1 Nov 2019 12:38:27 +0000 (13:38 +0100)
When update is offloaded, the coordinates are needed on host if
an Ewald surface contribution is to be computed. This flag is set
to reflect that and used when the decision of whether the
coordinates are to be copied on this step is made.

Change-Id: Iba40387f065e6fae2551540136e925f7c764d6aa

src/gromacs/mdlib/force.cpp
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdrun/runner.cpp
src/gromacs/mdtypes/inputrec.cpp
src/gromacs/mdtypes/inputrec.h
src/gromacs/mdtypes/simulation_workload.h
src/gromacs/taskassignment/decidesimulationworkload.cpp
src/gromacs/taskassignment/decidesimulationworkload.h

index 99f0b97a21b65741d084e31dbac61c9831365b2d..a82e5a98af5b7a42eb819b31d42fb0447e646ae2 100644 (file)
@@ -195,16 +195,14 @@ do_force_lowlevel(t_forcerec                               *fr,
         thisRankHasDuty(cr, DUTY_PME) &&
         (pme_run_mode(fr->pmedata) == PmeRunMode::CPU);
 
-    const bool haveEwaldSurfaceTerms =
-        EEL_PME_EWALD(fr->ic->eeltype) &&
-        (ir->ewald_geometry != eewg3D || ir->epsilon_surface != 0);
+    const bool haveEwaldSurfaceTerm = haveEwaldSurfaceContribution(*ir);
 
     /* Do long-range electrostatics and/or LJ-PME
      * and compute PME surface terms when necessary.
      */
     if (computePmeOnCpu ||
         fr->ic->eeltype == eelEWALD ||
-        haveEwaldSurfaceTerms)
+        haveEwaldSurfaceTerm)
     {
         int  status            = 0;
         real Vlr_q             = 0, Vlr_lj = 0;
@@ -217,8 +215,8 @@ do_force_lowlevel(t_forcerec                               *fr,
 
         if (EEL_PME_EWALD(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype))
         {
-            /* Calculate Ewald surface terms, when necessary */
-            if (haveEwaldSurfaceTerms)
+            /* Calculate the Ewald surface force and energy contributions, when necessary */
+            if (haveEwaldSurfaceTerm)
             {
                 wallcycle_sub_start(wcycle, ewcsEWALD_CORRECTION);
 
index 6a2cfbe34e7e6ac1fa75abc0f7a32206a40be7f9..f5a73578ba7d711fbd31801ea756b71373224771 100644 (file)
@@ -799,10 +799,9 @@ setupDomainLifetimeWorkload(const t_inputrec         &inputrec,
     // We assume we have local force work if there are CPU
     // force tasks including PME or nonbondeds.
     domainWork.haveCpuLocalForceWork  = domainWork.haveSpecialForces || domainWork.haveCpuListedForceWork || domainWork.haveFreeEnergyWork ||
-        simulationWork.useCpuNonbonded || simulationWork.useCpuPme ||
-        (EEL_PME_EWALD(inputrec.coulombtype) && (inputrec.ewald_geometry == eewg3DC ||
-                                                 inputrec.epsilon_surface != 0)) ||
+        simulationWork.useCpuNonbonded || simulationWork.useCpuPme || simulationWork.haveEwaldSurfaceContribution ||
         inputrec.nwall > 0;
+
     return domainWork;
 }
 
index 1c2fb3e42d3044b315e693b900918260f66cd9ac..dc647f63613073cace817c2fcf31f44d78181707 100644 (file)
@@ -1583,7 +1583,8 @@ int Mdrunner::mdrunner()
                                                                   useGpuForUpdate,
                                                                   devFlags.enableGpuBufferOps,
                                                                   devFlags.enableGpuHaloExchange,
-                                                                  devFlags.enableGpuPmePPComm);
+                                                                  devFlags.enableGpuPmePPComm,
+                                                                  haveEwaldSurfaceContribution(*inputrec));
 
         std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
         if (gpusWereDetected && ((useGpuForPme && thisRankHasDuty(cr, DUTY_PME)) || runScheduleWork.simulationWork.useGpuBufferOps))
index 7445dfc4cc7f0b9074924d2e33e86362b0cdf32f..354d39f8ea08000a5a678e86597799a393651423 100644 (file)
@@ -1549,3 +1549,9 @@ real maxReferenceTemperature(const t_inputrec &ir)
 
     return maxTemperature;
 }
+
+bool haveEwaldSurfaceContribution(const t_inputrec &ir)
+{
+    return EEL_PME_EWALD(ir.coulombtype) && (ir.ewald_geometry == eewg3DC ||
+                                             ir.epsilon_surface != 0);
+}
index 6cff9d6df6458df0b88706a8f322f497cc41c549..86dc181cd197ebd2ac14b7c95043d37956e84b8d 100644 (file)
@@ -693,4 +693,8 @@ int ndof_com(const t_inputrec *ir);
  */
 real maxReferenceTemperature(const t_inputrec &ir);
 
+/*! \brief Returns whether there is an Ewald surface contribution
+ */
+bool haveEwaldSurfaceContribution(const t_inputrec &ir);
+
 #endif /* GMX_MDTYPES_INPUTREC_H */
index 6594081c509036399646ee972944d7f6278be54d..6100e70e6d473b6445845127878fa184e2c57aa8 100644 (file)
@@ -158,6 +158,8 @@ class SimulationWorkload
         bool useGpuPmePpCommunication  = false;
         //! If direct GPU-GPU communication is enabled.
         bool useGpuDirectCommunication = false;
+        //! If there is an Ewald surface (dipole) term to compute
+        bool haveEwaldSurfaceContribution = false;
 };
 
 class MdrunScheduleWorkload
index ea82f95f42bccd1eef0cef622a3d343a0204ec28..933aeaef6757bce127cb44d19e50019de6aea467 100644 (file)
 namespace gmx
 {
 
-SimulationWorkload createSimulationWorkload(bool       useGpuForNonbonded,
-                                            PmeRunMode pmeRunMode,
-                                            bool       useGpuForBonded,
-                                            bool       useGpuForUpdate,
-                                            bool       useGpuForBufferOps,
-                                            bool       useGpuHaloExchange,
-                                            bool       useGpuPmePpComm)
+SimulationWorkload createSimulationWorkload(bool              useGpuForNonbonded,
+                                            PmeRunMode        pmeRunMode,
+                                            bool              useGpuForBonded,
+                                            bool              useGpuForUpdate,
+                                            bool              useGpuForBufferOps,
+                                            bool              useGpuHaloExchange,
+                                            bool              useGpuPmePpComm,
+                                            bool              haveEwaldSurfaceContribution)
 {
     SimulationWorkload simulationWorkload;
-    simulationWorkload.useCpuNonbonded           = !useGpuForNonbonded;
-    simulationWorkload.useGpuNonbonded           = useGpuForNonbonded;
-    simulationWorkload.useCpuPme                 = (pmeRunMode == PmeRunMode::CPU);
-    simulationWorkload.useGpuPme                 = (pmeRunMode == PmeRunMode::GPU || pmeRunMode == PmeRunMode::Mixed);
-    simulationWorkload.useGpuPmeFft              = (pmeRunMode == PmeRunMode::Mixed);
-    simulationWorkload.useGpuBonded              = useGpuForBonded;
-    simulationWorkload.useGpuUpdate              = useGpuForUpdate;
-    simulationWorkload.useGpuBufferOps           = useGpuForBufferOps || useGpuForUpdate;
-    simulationWorkload.useGpuHaloExchange        = useGpuHaloExchange;
-    simulationWorkload.useGpuPmePpCommunication  = useGpuPmePpComm;
-    simulationWorkload.useGpuDirectCommunication = useGpuHaloExchange || useGpuPmePpComm;
+    simulationWorkload.useCpuNonbonded              = !useGpuForNonbonded;
+    simulationWorkload.useGpuNonbonded              = useGpuForNonbonded;
+    simulationWorkload.useCpuPme                    = (pmeRunMode == PmeRunMode::CPU);
+    simulationWorkload.useGpuPme                    = (pmeRunMode == PmeRunMode::GPU || pmeRunMode == PmeRunMode::Mixed);
+    simulationWorkload.useGpuPmeFft                 = (pmeRunMode == PmeRunMode::Mixed);
+    simulationWorkload.useGpuBonded                 = useGpuForBonded;
+    simulationWorkload.useGpuUpdate                 = useGpuForUpdate;
+    simulationWorkload.useGpuBufferOps              = useGpuForBufferOps || useGpuForUpdate;
+    simulationWorkload.useGpuHaloExchange           = useGpuHaloExchange;
+    simulationWorkload.useGpuPmePpCommunication     = useGpuPmePpComm;
+    simulationWorkload.useGpuDirectCommunication    = useGpuHaloExchange || useGpuPmePpComm;
+    simulationWorkload.haveEwaldSurfaceContribution = haveEwaldSurfaceContribution;
 
     return simulationWorkload;
 }
index b389da53c168a36d98396f004d532f196f800a17..ef53314b392f226fc44cd5de9cc50ef6b6db4686 100644 (file)
@@ -44,6 +44,7 @@
 
 #include <vector>
 
+#include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/simulation_workload.h"
 
 enum class PmeRunMode;
@@ -55,25 +56,26 @@ namespace gmx
  * Build datastructure that contains decisions whether to run different workload
  * task on GPUs.
  *
- * \param[in] useGpuForNonbonded If we have short-range nonbonded interactions
+ * \param[in] useGpuForNonbonded Whether we have short-range nonbonded interactions
  *                               calculations on GPU(s).
  * \param[in] pmeRunMode         Run mode indicating what resource is PME execured on.
- * \param[in] useGpuForBonded    If bonded interactions are calculated on GPU(s).
- * \param[in] useGpuForUpdate    If coordinate update and constraint solving is performed on
+ * \param[in] useGpuForBonded    Whether bonded interactions are calculated on GPU(s).
+ * \param[in] useGpuForUpdate    Whether coordinate update and constraint solving is performed on
  *                               GPU(s).
- * \param[in] useGpuForBufferOps If buffer ops / reduction are calculated on GPU(s).
- * \param[in] useGpuHaloExchange If GPU direct communication is used in halo exchange.
- * \param[in] useGpuPmePpComm    If GPu direct communication is used in PME-PP communication.
+ * \param[in] useGpuForBufferOps Whether buffer ops / reduction are calculated on GPU(s).
+ * \param[in] useGpuHaloExchange Whether GPU direct communication is used in halo exchange.
+ * \param[in] useGpuPmePpComm    Whether GPU direct communication is used in PME-PP communication.
+ * \param[in] haveEwaldSurfaceContribution Whether there is an Ewald surface contribution
  * \returns Simulation lifetime constant workload description.
  */
-SimulationWorkload createSimulationWorkload(bool       useGpuForNonbonded,
-                                            PmeRunMode pmeRunMode,
-                                            bool       useGpuForBonded,
-                                            bool       useGpuForUpdate,
-                                            bool       useGpuForBufferOps,
-                                            bool       useGpuHaloExchange,
-                                            bool       useGpuPmePpComm);
-
+SimulationWorkload createSimulationWorkload(bool              useGpuForNonbonded,
+                                            PmeRunMode        pmeRunMode,
+                                            bool              useGpuForBonded,
+                                            bool              useGpuForUpdate,
+                                            bool              useGpuForBufferOps,
+                                            bool              useGpuHaloExchange,
+                                            bool              useGpuPmePpComm,
+                                            bool              haveEwaldSurfaceContribution);
 
 }  // namespace gmx