Conditionally pin GPU-related grid data
authorSzilárd Páll <pall.szilard@gmail.com>
Mon, 15 Jul 2019 11:26:36 +0000 (13:26 +0200)
committerBerk Hess <hess@kth.se>
Mon, 22 Jul 2019 08:10:11 +0000 (10:10 +0200)
Data that is transferred to the GPU when the buffer ops is offloaded is
now only pinned when the nonbonded module uses GPU offload avoidign the
runtime errors encountered when a GPU-enabled build does not detect a
GPU and therefore the CUDA runtime refuses to register the memory.

Refs #2817 #2934

Change-Id: Iabbc0d9f37fad0e88cd39a078af1346e8f713ec1

src/gromacs/nbnxm/grid.cpp
src/gromacs/nbnxm/grid.h
src/gromacs/nbnxm/gridset.cpp
src/gromacs/nbnxm/gridset.h
src/gromacs/nbnxm/nbnxm_setup.cpp
src/gromacs/nbnxm/pairsearch.cpp
src/gromacs/nbnxm/pairsearch.h

index a852b265d929800063f57123a120449ca237ac86..c69a9e8ac46128fe391ca196fb6ab67cc15af9f7 100644 (file)
@@ -101,13 +101,14 @@ static real gridAtomDensity(int        numAtoms,
     return numAtoms/(size[XX]*size[YY]*size[ZZ]);
 }
 
-void Grid::setDimensions(const int           ddZone,
-                         const int           numAtoms,
-                         const rvec          lowerCorner,
-                         const rvec          upperCorner,
-                         real                atomDensity,
-                         const real          maxAtomGroupRadius,
-                         const bool          haveFep)
+void Grid::setDimensions(const int            ddZone,
+                         const int            numAtoms,
+                         const rvec           lowerCorner,
+                         const rvec           upperCorner,
+                         real                 atomDensity,
+                         const real           maxAtomGroupRadius,
+                         const bool           haveFep,
+                         gmx::PinningPolicy   pinningPolicy)
 {
     /* For the home zone we compute the density when not set (=-1) or when =0 */
     if (ddZone == 0 && atomDensity <= 0)
@@ -182,8 +183,8 @@ void Grid::setDimensions(const int           ddZone,
     /* We need one additional cell entry for particles moved by DD */
     cxy_na_.resize(numColumns() + 1);
     cxy_ind_.resize(numColumns() + 2);
-    changePinningPolicy(&cxy_na_, gmx::PinningPolicy::PinnedIfSupported);
-    changePinningPolicy(&cxy_ind_, gmx::PinningPolicy::PinnedIfSupported);
+    changePinningPolicy(&cxy_na_, pinningPolicy);
+    changePinningPolicy(&cxy_ind_, pinningPolicy);
 
     /* Worst case scenario of 1 atom in each last cell */
     int maxNumCells;
index 488cb71c7b5e0200857f48acd9378beddc224988..a57d3ad684b692448d4459a09498908a6187d52c 100644 (file)
@@ -464,7 +464,8 @@ class Grid
                            const rvec            upperCorner,
                            real                  atomDensity,
                            real                  maxAtomGroupRadius,
-                           bool                  haveFep);
+                           bool                  haveFep,
+                           gmx::PinningPolicy    pinningPolicy);
 
         //! Sets the cell indices using indices in \p gridSetData and \p gridWork
         void setCellIndices(int                             ddZone,
index 8c82a8a9aaa0acd719405a9bcd26180fac39244a..34369acf9498be2d1dbc67a48765827d3c1a2211 100644 (file)
@@ -87,7 +87,8 @@ GridSet::GridSet(const int                 ePBC,
                  const gmx_domdec_zones_t *ddZones,
                  const PairlistType        pairlistType,
                  const bool                haveFep,
-                 const int                 numThreads) :
+                 const int                 numThreads,
+                 gmx::PinningPolicy        pinningPolicy) :
     domainSetup_(ePBC, numDDCells, ddZones),
     grids_(numDDZones(domainSetup_.haveMultipleDomainsPerDim), Grid(pairlistType, haveFep_)),
     haveFep_(haveFep),
@@ -96,8 +97,8 @@ GridSet::GridSet(const int                 ePBC,
     gridWork_(numThreads)
 {
     clear_mat(box_);
-    changePinningPolicy(&gridSetData_.cells, gmx::PinningPolicy::PinnedIfSupported);
-    changePinningPolicy(&gridSetData_.atomIndices, gmx::PinningPolicy::PinnedIfSupported);
+    changePinningPolicy(&gridSetData_.cells, pinningPolicy);
+    changePinningPolicy(&gridSetData_.atomIndices, pinningPolicy);
 }
 
 void GridSet::setLocalAtomOrder()
@@ -181,11 +182,14 @@ void GridSet::putOnGrid(const matrix                    box,
     /* We always use the home zone (grid[0]) for setting the cell size,
      * since determining densities for non-local zones is difficult.
      */
+    // grid data used in GPU transfers inherits the gridset pinnin policy
+    auto pinPolicy = gridSetData_.cells.get_allocator().pinningPolicy();
     grid.setDimensions(ddZone, n - numAtomsMoved,
                        lowerCorner, upperCorner,
                        atomDensity,
                        maxAtomGroupRadius,
-                       haveFep_);
+                       haveFep_,
+                       pinPolicy);
 
     for (GridWork &work : gridWork_)
     {
index 71f1cd101bfd59728dcb676914c79b6a1231bb3c..436438196e9a9967b199d0a5acf71ca7d3551ad8 100644 (file)
@@ -103,7 +103,8 @@ class GridSet
                 const gmx_domdec_zones_t *ddZones,
                 PairlistType              pairlistType,
                 bool                      haveFep,
-                int                       numThreads);
+                int                       numThreads,
+                gmx::PinningPolicy        pinningPolicy);
 
         //! Puts the atoms in \p ddZone on the grid and copies the coordinates to \p nbat
         void putOnGrid(const matrix                    box,
index 794f6f98d42b06d9d1ed71deb465bf837d258d18..18c8182dd5c944865856034b83581eaa51bcb7db 100644 (file)
@@ -427,10 +427,11 @@ init_nb_verlet(const gmx::MDLogger     &mdlog,
         enbnxninitcombrule = enbnxninitcombruleNONE;
     }
 
-    auto nbat =
-        std::make_unique<nbnxn_atomdata_t>(useGpu ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
+    auto pinPolicy = (useGpu ? gmx::PinningPolicy::PinnedIfSupported : gmx::PinningPolicy::CannotBePinned);
 
-    int mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
+    auto nbat      = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
+
+    int  mimimumNumEnergyGroupNonbonded = ir->opts.ngener;
     if (ir->opts.ngener - ir->nwall == 1)
     {
         /* We have only one non-wall energy group, we do not need energy group
@@ -474,7 +475,8 @@ init_nb_verlet(const gmx::MDLogger     &mdlog,
                                      DOMAINDECOMP(cr) ? domdec_zones(cr->dd) : nullptr,
                                      pairlistParams.pairlistType,
                                      bFEP_NonBonded,
-                                     gmx_omp_nthreads_get(emntPairsearch));
+                                     gmx_omp_nthreads_get(emntPairsearch),
+                                     pinPolicy);
 
     return std::make_unique<nonbonded_verlet_t>(std::move(pairlistSets),
                                                 std::move(pairSearch),
index d43e39be253cfd4fe459da597a58c66b778b6644..021f6cb25a4f1044df4167de6198b058a0ddc267 100644 (file)
@@ -113,8 +113,9 @@ PairSearch::PairSearch(const int                 ePBC,
                        const gmx_domdec_zones_t *ddZones,
                        const PairlistType        pairlistType,
                        const bool                haveFep,
-                       const int                 maxNumThreads) :
-    gridSet_(ePBC, numDDCells, ddZones, pairlistType, haveFep, maxNumThreads),
+                       const int                 maxNumThreads,
+                       gmx::PinningPolicy        pinningPolicy) :
+    gridSet_(ePBC, numDDCells, ddZones, pairlistType, haveFep, maxNumThreads, pinningPolicy),
     work_(maxNumThreads)
 {
     cycleCounting_.recordCycles_ = (getenv("GMX_NBNXN_CYCLE") != nullptr);
index b4ec14b72c1775615ac7d1e32768cc05a29a0e9b..10ab4d15f9bac7289f7e48361c28ff323d685ac1 100644 (file)
@@ -211,7 +211,8 @@ class PairSearch
                    const gmx_domdec_zones_t *zones,
                    PairlistType              pairlistType,
                    bool                      haveFep,
-                   int                       maxNumthreads);
+                   int                       maxNumthreads,
+                   gmx::PinningPolicy        pinningPolicy);
 
         //! Sets the order of the local atoms to the order grid atom ordering
         void setLocalAtomOrder()