remove InteractionLocality from PairListSet constructor
authorejjordan <ejjordan@kth.se>
Sun, 21 Feb 2021 18:47:09 +0000 (19:47 +0100)
committerArtem Zhmurov <zhmurov@gmail.com>
Mon, 22 Feb 2021 14:56:24 +0000 (14:56 +0000)
The InteractionLocality is removed as a member from PairListSet
and is no longer passed to the constructor, but instead to
PairListSets::construct. This will make it easier to remove the
explicit dependence on domdec_zones in nbnxm module.

src/gromacs/nbnxm/nbnxm_setup.cpp
src/gromacs/nbnxm/pairlist.cpp
src/gromacs/nbnxm/pairlistset.h

index 193af53cbe4609e8e4531eb102bce464ea98c45d..64501b5ab3682ed5bacf9d7bca8e4f60a771d0db 100644 (file)
@@ -297,11 +297,11 @@ PairlistSets::PairlistSets(const PairlistParams& pairlistParams,
     params_(pairlistParams),
     minimumIlistCountForGpuBalancing_(minimumIlistCountForGpuBalancing)
 {
-    localSet_ = std::make_unique<PairlistSet>(gmx::InteractionLocality::Local, params_);
+    localSet_ = std::make_unique<PairlistSet>(params_);
 
     if (haveMultipleDomains)
     {
-        nonlocalSet_ = std::make_unique<PairlistSet>(gmx::InteractionLocality::NonLocal, params_);
+        nonlocalSet_ = std::make_unique<PairlistSet>(params_);
     }
 }
 
index 8d8186beb1b7b2624b0c588ddf504196ff1539fe..9d6f7a2b73f40fb92defbfee91e07c0aabde30f1 100644 (file)
@@ -654,8 +654,7 @@ NbnxnPairlistGpu::NbnxnPairlistGpu(gmx::PinningPolicy pinningPolicy) :
 }
 
 // TODO: Move to pairlistset.cpp
-PairlistSet::PairlistSet(const InteractionLocality locality, const PairlistParams& pairlistParams) :
-    locality_(locality),
+PairlistSet::PairlistSet(const PairlistParams& pairlistParams) :
     params_(pairlistParams),
     combineLists_(sc_isGpuPairListType[pairlistParams.pairlistType]), // Currently GPU lists are always combined
     isCpuType_(!sc_isGpuPairListType[pairlistParams.pairlistType])
@@ -3868,7 +3867,8 @@ static Range<int> getJZoneRange(const gmx_domdec_zones_t* ddZones,
 //! Prepares CPU lists produced by the search for dynamic pruning
 static void prepareListsForDynamicPruning(gmx::ArrayRef<NbnxnPairlistCpu> lists);
 
-void PairlistSet::constructPairlists(const Nbnxm::GridSet&         gridSet,
+void PairlistSet::constructPairlists(gmx::InteractionLocality      locality,
+                                     const Nbnxm::GridSet&         gridSet,
                                      gmx::ArrayRef<PairsearchWork> searchWork,
                                      nbnxn_atomdata_t*             nbat,
                                      const ListOfLists<int>&       exclusions,
@@ -3887,7 +3887,7 @@ void PairlistSet::constructPairlists(const Nbnxm::GridSet&         gridSet,
 
     nbat->bUseBufferFlags = (nbat->out.size() > 1);
     /* We should re-init the flags before making the first list */
-    if (nbat->bUseBufferFlags && locality_ == InteractionLocality::Local)
+    if (nbat->bUseBufferFlags && locality == InteractionLocality::Local)
     {
         resizeAndZeroBufferFlags(&nbat->buffer_flags, nbat->numAtoms());
     }
@@ -3897,7 +3897,7 @@ void PairlistSet::constructPairlists(const Nbnxm::GridSet&         gridSet,
     if (!isCpuType_ && minimumIlistCountForGpuBalancing > 0)
     {
         get_nsubpair_target(
-                gridSet, locality_, rlist, minimumIlistCountForGpuBalancing, &nsubpair_target, &nsubpair_tot_est);
+                gridSet, locality, rlist, minimumIlistCountForGpuBalancing, &nsubpair_target, &nsubpair_tot_est);
     }
 
     /* Clear all pair-lists */
@@ -3919,16 +3919,16 @@ void PairlistSet::constructPairlists(const Nbnxm::GridSet&         gridSet,
     }
 
     const gmx_domdec_zones_t* ddZones = gridSet.domainSetup().zones;
-    GMX_ASSERT(locality_ == InteractionLocality::Local || ddZones != nullptr,
+    GMX_ASSERT(locality == InteractionLocality::Local || ddZones != nullptr,
                "Nonlocal interaction locality with null ddZones.");
 
-    const auto iZoneRange = getIZoneRange(gridSet.domainSetup(), locality_);
+    const auto iZoneRange = getIZoneRange(gridSet.domainSetup(), locality);
 
     for (const int iZone : iZoneRange)
     {
         const Grid& iGrid = gridSet.grids()[iZone];
 
-        const auto jZoneRange = getJZoneRange(ddZones, locality_, iZone);
+        const auto jZoneRange = getJZoneRange(ddZones, locality, iZone);
 
         for (int jZone : jZoneRange)
         {
@@ -3947,7 +3947,7 @@ void PairlistSet::constructPairlists(const Nbnxm::GridSet&         gridSet,
             /* With GPU: generate progressively smaller lists for
              * load balancing for local only or non-local with 2 zones.
              */
-            const bool progBal = (locality_ == InteractionLocality::Local || ddZones->n <= 2);
+            const bool progBal = (locality == InteractionLocality::Local || ddZones->n <= 2);
 
 #pragma omp parallel for num_threads(numLists) schedule(static)
             for (int th = 0; th < numLists; th++)
@@ -4180,7 +4180,8 @@ void PairlistSets::construct(const InteractionLocality iLocality,
             "exclusions should either be empty or the number of lists should match the number of "
             "local i-atoms");
 
-    pairlistSet(iLocality).constructPairlists(gridSet,
+    pairlistSet(iLocality).constructPairlists(iLocality,
+                                              gridSet,
                                               pairSearch->work(),
                                               nbat,
                                               exclusions,
index 3728dc9863a83f750efff75d0c76e8bdb51817ad..c7f80e9f4859401e4c4a0ea5dcad06da78fc1997 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,2020,2021, 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.
@@ -83,12 +83,13 @@ class PairlistSet
 {
 public:
     //! Constructor: initializes the pairlist set as empty
-    PairlistSet(gmx::InteractionLocality locality, const PairlistParams& listParams);
+    PairlistSet(const PairlistParams& listParams);
 
     ~PairlistSet();
 
     //! Constructs the pairlists in the set using the coordinates in \p nbat
-    void constructPairlists(const Nbnxm::GridSet&         gridSet,
+    void constructPairlists(gmx::InteractionLocality      locality,
+                            const Nbnxm::GridSet&         gridSet,
                             gmx::ArrayRef<PairsearchWork> searchWork,
                             nbnxn_atomdata_t*             nbat,
                             const gmx::ListOfLists<int>&  exclusions,
@@ -99,9 +100,6 @@ public:
     //! Dispatch the kernel for dynamic pairlist pruning
     void dispatchPruneKernel(const nbnxn_atomdata_t* nbat, const rvec* shift_vec);
 
-    //! Returns the locality
-    gmx::InteractionLocality locality() const { return locality_; }
-
     //! Returns the lists of CPU pairlists
     gmx::ArrayRef<const NbnxnPairlistCpu> cpuLists() const { return cpuLists_; }
 
@@ -122,8 +120,6 @@ public:
     gmx::ArrayRef<const std::unique_ptr<t_nblist>> fepLists() const { return fepLists_; }
 
 private:
-    //! The locality of the pairlist set
-    gmx::InteractionLocality locality_;
     //! List of pairlists in CPU layout
     std::vector<NbnxnPairlistCpu> cpuLists_;
     //! List of working list for rebalancing CPU lists