Start modularizing force workload setup
authorSzilárd Páll <pall.szilard@gmail.com>
Tue, 19 Feb 2019 09:44:14 +0000 (10:44 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 13 Mar 2019 21:24:49 +0000 (22:24 +0100)
Extend ppForceWorkload flags with listed forces-related booleans
and add a setup function that re-initializes ppForceWorkload every
search step.

Change-Id: I71c783cf22e58e880d7d8655bc1e9b44e6428407

src/gromacs/listed_forces/listed_forces.cpp
src/gromacs/listed_forces/listed_forces.h
src/gromacs/mdlib/ppforceworkload.h
src/gromacs/mdlib/sim_util.cpp

index d99dcfe8c8fa5e0f4e3c0bc78b4dd95d29d55333..22e2c0b1a62f82475dbf3a7a79619045ae4d055f 100644 (file)
@@ -589,6 +589,28 @@ calcBondedForces(const t_idef     *idef,
     }
 }
 
+bool havePositionRestraints(const t_idef   &idef,
+                            const t_fcdata &fcd)
+{
+    return
+        ((idef.il[F_POSRES].nr > 0) ||
+         (idef.il[F_FBPOSRES].nr > 0) ||
+         fcd.orires.nr > 0 ||
+         fcd.disres.nres > 0);
+}
+
+bool haveCpuBondeds(const t_forcerec &fr)
+{
+    return fr.bondedThreading->haveBondeds;
+}
+
+bool haveCpuListedForces(const t_forcerec &fr,
+                         const t_idef     &idef,
+                         const t_fcdata   &fcd)
+{
+    return haveCpuBondeds(fr) || havePositionRestraints(idef, fcd);
+}
+
 void calc_listed(const t_commrec             *cr,
                  const gmx_multisim_t *ms,
                  struct gmx_wallcycle        *wcycle,
@@ -621,10 +643,7 @@ void calc_listed(const t_commrec             *cr,
         pbc_null = nullptr;
     }
 
-    if ((idef->il[F_POSRES].nr > 0) ||
-        (idef->il[F_FBPOSRES].nr > 0) ||
-        fcd->orires.nr > 0 ||
-        fcd->disres.nres > 0)
+    if (havePositionRestraints(*idef, *fcd))
     {
         /* TODO Use of restraints triggers further function calls
            inside the loop over calc_one_bond(), but those are too
@@ -673,7 +692,7 @@ void calc_listed(const t_commrec             *cr,
         wallcycle_sub_stop(wcycle, ewcsRESTRAINTS);
     }
 
-    if (bt->haveBondeds)
+    if (haveCpuBondeds(*fr))
     {
         wallcycle_sub_start(wcycle, ewcsLISTED);
         /* The dummy array is to have a place to store the dhdl at other values
index 165ebffc7df08a9bc23bf5c47aaf892fd2c63479..5d7999f8de475e350a59db4a94dc196a210da595 100644 (file)
@@ -157,4 +157,20 @@ do_force_listed(struct gmx_wallcycle           *wcycle,
                 int                            *global_atom_index,
                 int                             flags);
 
+/*! \brief Returns true if there are position restraints. */
+bool havePositionRestraints(const t_idef   &idef,
+                            const t_fcdata &fcd);
+
+/*! \brief Returns true if there are CPU (i.e. not GPU-offloaded) bonded interactions to compute. */
+bool haveCpuBondeds(const t_forcerec &fr);
+
+/*! \brief Returns true if there are listed interactions to compute.
+ *
+ * NOTE: the current implementation returns true if there are position restraints
+ * or any bonded interactions computed on the CPU.
+ */
+bool haveCpuListedForces(const t_forcerec &fr,
+                         const t_idef     &idef,
+                         const t_fcdata   &fcd);
+
 #endif
index 9174c8b770218f45c241e2c46c1aebb1e15c3a33..a494ce57f327c30aa16084e62bfb07856726b3d8 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019, 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.
@@ -66,6 +66,15 @@ class PpForceWorkload
     public:
         //! Whether this MD step has bonded work to run on a GPU.
         bool haveGpuBondedWork = false;
+        //! Whether this MD step has bonded work to run on he CPU.
+        bool haveCpuBondedWork = false;
+        //! Whether this MD step has listed forces bonded work to run on he CPU.
+        bool haveRestraintsWork = false;
+        //! Whether this MD step has listed forces work to run on he CPU.
+        //  Note: currently this is haveCpuBondedWork | haveRestraintsWork
+        bool haveCpuListedForceWork = false;
+        //! Whether this MD step has special forces on the CPU.
+        bool haveSpecialForces = false;
 };
 
 } // namespace gmx
index 0457c20f0c89f7105765d9c3aafb329b7d6254ef..324ee3f79bdcf9f943166928eddc49efaf133793 100644 (file)
@@ -65,6 +65,7 @@
 #include "gromacs/listed_forces/bonded.h"
 #include "gromacs/listed_forces/disre.h"
 #include "gromacs/listed_forces/gpubonded.h"
+#include "gromacs/listed_forces/listed_forces.h"
 #include "gromacs/listed_forces/manage_threading.h"
 #include "gromacs/listed_forces/orires.h"
 #include "gromacs/math/arrayrefwithpadding.h"
@@ -466,6 +467,26 @@ static void checkPotentialEnergyValidity(int64_t               step,
     }
 }
 
+/*! \brief Return true if there are special forces computed this step.
+ *
+ * The conditionals exactly correspond to those in computeSpecialForces().
+ */
+static bool
+haveSpecialForces(const t_inputrec              *inputrec,
+                  ForceProviders                *forceProviders,
+                  int                            forceFlags,
+                  const gmx_edsam               *ed)
+{
+    const bool computeForces = (forceFlags & GMX_FORCE_FORCES) != 0;
+
+    return
+        ((computeForces && forceProviders->hasForceProvider()) ||         // forceProviders
+         (inputrec->bPull && pull_have_potential(inputrec->pull_work)) || // pull
+         inputrec->bRot ||                                                // enforced rotation
+         (ed != nullptr) ||                                               // flooding
+         (inputrec->bIMD && computeForces));                              // IMD
+}
+
 /*! \brief Compute forces and/or energies for special algorithms
  *
  * The intention is to collect all calls to algorithms that compute
@@ -758,6 +779,30 @@ setupForceOutputs(const t_forcerec                    *fr,
 }
 
 
+/*! \brief Set up flags that indicate what type of work is there to compute.
+ *
+ * Currently we only update it at search steps,
+ * but some properties may change more frequently (e.g. virial/non-virial step),
+ * so when including those either the frequency of update (per-step) or the scope
+ * of a flag will change (i.e. a set of flags for nstlist steps).
+ *
+ */
+static void
+setupForceWorkload(gmx::PpForceWorkload *forceWork,
+                   const t_inputrec     *inputrec,
+                   const t_forcerec     *fr,
+                   const gmx_edsam      *ed,
+                   const t_idef         &idef,
+                   const t_fcdata       *fcd,
+                   const int             forceFlags
+                   )
+{
+    forceWork->haveSpecialForces      = haveSpecialForces(inputrec, fr->forceProviders, forceFlags, ed);
+    forceWork->haveCpuBondedWork      = haveCpuBondeds(*fr);
+    forceWork->haveRestraintsWork     = havePositionRestraints(idef, *fcd);
+    forceWork->haveCpuListedForceWork = haveCpuListedForces(*fr, idef, *fcd);
+}
+
 static void do_force_cutsVERLET(FILE *fplog,
                                 const t_commrec *cr,
                                 const gmx_multisim_t *ms,
@@ -812,6 +857,17 @@ static void do_force_cutsVERLET(FILE *fplog,
         ((flags & GMX_FORCE_ENERGY) ? GMX_PME_CALC_ENER_VIR : 0) |
         ((flags & GMX_FORCE_FORCES) ? GMX_PME_CALC_F : 0);
 
+    if (bNS)
+    {
+        setupForceWorkload(ppForceWorkload,
+                           inputrec,
+                           fr,
+                           ed,
+                           top->idef,
+                           fcd,
+                           flags);
+    }
+
     /* At a search step we need to start the first balancing region
      * somewhere early inside the step after communication during domain
      * decomposition (and not during the previous step as usual).