Add simulation workload flag for dipole computation
authorSzilárd Páll <pall.szilard@gmail.com>
Mon, 3 Feb 2020 10:38:50 +0000 (11:38 +0100)
committerArtem Zhmurov <zhmurov@gmail.com>
Mon, 17 Feb 2020 21:21:19 +0000 (22:21 +0100)
Change-Id: I4e7828d4dc73fb9caea100e0303779853be80c70

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

index 61878b8f5af58b95b8d7c6b1c05183d001ced540..54dae71d61b7e961c67cc5a065d38d2cff54ffed 100644 (file)
@@ -955,16 +955,13 @@ void do_force(FILE*                               fplog,
 
     clear_mat(vir_force);
 
-    if (stepWork.stateChanged)
+    if (stepWork.stateChanged && simulationWork.computeMuTot)
     {
-        if (inputrecNeedMutot(inputrec))
-        {
-            /* Calculate total (local) dipole moment in a temporary common array.
-             * This makes it possible to sum them over nodes faster.
-             */
-            calc_mu(start, homenr, x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB,
-                    mdatoms->nChargePerturbed, mu, mu + DIM);
-        }
+        /* Calculate total (local) dipole moment in a temporary common array.
+         * This makes it possible to sum them over nodes faster.
+         */
+        calc_mu(start, homenr, x.unpaddedArrayRef(), mdatoms->chargeA, mdatoms->chargeB,
+                mdatoms->nChargePerturbed, mu, mu + DIM);
     }
 
     if (fr->pbcType != PbcType::No)
@@ -1355,7 +1352,7 @@ void do_force(FILE*                               fplog,
         wallcycle_stop(wcycle, ewcLAUNCH_GPU);
     }
 
-    if (stepWork.stateChanged && inputrecNeedMutot(inputrec))
+    if (stepWork.stateChanged && simulationWork.computeMuTot)
     {
         if (PAR(cr))
         {
index 41f609d1cc3659226b6da206881e5fa6e577f6a9..f68714ce70b21e40f0aa19d13d2b6cc91d3d76af 100644 (file)
@@ -1557,10 +1557,10 @@ int Mdrunner::mdrunner()
         // make it work.
         MdrunScheduleWorkload runScheduleWork;
         // Also populates the simulation constant workload description.
-        runScheduleWork.simulationWork = createSimulationWorkload(
-                useGpuForNonbonded, pmeRunMode, useGpuForBonded, useGpuForUpdate,
-                devFlags.enableGpuBufferOps, devFlags.enableGpuHaloExchange,
-                devFlags.enableGpuPmePPComm, haveEwaldSurfaceContribution(*inputrec));
+        runScheduleWork.simulationWork =
+                createSimulationWorkload(*inputrec, useGpuForNonbonded, pmeRunMode, useGpuForBonded,
+                                         useGpuForUpdate, devFlags.enableGpuBufferOps,
+                                         devFlags.enableGpuHaloExchange, devFlags.enableGpuPmePPComm);
 
         std::unique_ptr<gmx::StatePropagatorDataGpu> stateGpu;
         if (gpusWereDetected
index 9bf230a5938222e00dff0449441eb7351a502b9d..5dfe43c333c55c3519a5469f808af0cc2c38e88f 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018,2019, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,2020, 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.
@@ -142,6 +142,8 @@ public:
 class SimulationWorkload
 {
 public:
+    //! Whether total dipole needs to be computed
+    bool computeMuTot = false;
     //! If we have calculation of short range nonbondeds on CPU
     bool useCpuNonbonded = false;
     //! If we have calculation of short range nonbondeds on GPU
index 16849d4e5dd611856ca9e537b5b6a7663bcfd79f..b9ba44e94caeacc89e45f68272b4fc44002e2026 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2019, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020, 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.
 namespace gmx
 {
 
-SimulationWorkload createSimulationWorkload(bool       useGpuForNonbonded,
-                                            PmeRunMode pmeRunMode,
-                                            bool       useGpuForBonded,
-                                            bool       useGpuForUpdate,
-                                            bool       useGpuForBufferOps,
-                                            bool       useGpuHaloExchange,
-                                            bool       useGpuPmePpComm,
-                                            bool       haveEwaldSurfaceContribution)
+SimulationWorkload createSimulationWorkload(const t_inputrec& inputrec,
+                                            bool              useGpuForNonbonded,
+                                            PmeRunMode        pmeRunMode,
+                                            bool              useGpuForBonded,
+                                            bool              useGpuForUpdate,
+                                            bool              useGpuForBufferOps,
+                                            bool              useGpuHaloExchange,
+                                            bool              useGpuPmePpComm)
 {
     SimulationWorkload simulationWorkload;
+    simulationWorkload.computeMuTot    = inputrecNeedMutot(&inputrec);
     simulationWorkload.useCpuNonbonded = !useGpuForNonbonded;
     simulationWorkload.useGpuNonbonded = useGpuForNonbonded;
     simulationWorkload.useCpuPme       = (pmeRunMode == PmeRunMode::CPU);
@@ -71,7 +72,7 @@ SimulationWorkload createSimulationWorkload(bool       useGpuForNonbonded,
     simulationWorkload.useGpuHaloExchange       = useGpuHaloExchange;
     simulationWorkload.useGpuPmePpCommunication = useGpuPmePpComm && (pmeRunMode == PmeRunMode::GPU);
     simulationWorkload.useGpuDirectCommunication    = useGpuHaloExchange || useGpuPmePpComm;
-    simulationWorkload.haveEwaldSurfaceContribution = haveEwaldSurfaceContribution;
+    simulationWorkload.haveEwaldSurfaceContribution = haveEwaldSurfaceContribution(inputrec);
 
     return simulationWorkload;
 }
index 538fa823ed2bc4fbc1b9e248dd897416b450e9e6..b8d81c02b5107b2ca5f6671a74ee231077d983f7 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2019, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020, 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.
@@ -56,6 +56,7 @@ namespace gmx
  * Build datastructure that contains decisions whether to run different workload
  * task on GPUs.
  *
+ * \param[in] inputrec           The input record
  * \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.
@@ -65,17 +66,16 @@ namespace gmx
  * \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,
-                                            bool       haveEwaldSurfaceContribution);
+SimulationWorkload createSimulationWorkload(const t_inputrec& inputrec,
+                                            bool              useGpuForNonbonded,
+                                            PmeRunMode        pmeRunMode,
+                                            bool              useGpuForBonded,
+                                            bool              useGpuForUpdate,
+                                            bool              useGpuForBufferOps,
+                                            bool              useGpuHaloExchange,
+                                            bool              useGpuPmePpComm);
 
 } // namespace gmx