Add freeDeviceBuffer GPU utility function
authorAleksei Iupinov <a.yupinov@gmail.com>
Mon, 12 Feb 2018 09:16:15 +0000 (10:16 +0100)
committerSzilárd Páll <pall.szilard@gmail.com>
Thu, 15 Feb 2018 18:00:22 +0000 (19:00 +0100)
This allows to make the GPU pairlist deletion code look same
for CUDA and OpenCL. More cleanup TODOs are left.

Change-Id: I4f7775ae2bda65b69696f7aacdb1d3fea9c62ac5

src/gromacs/gpu_utils/cudautils.cuh
src/gromacs/gpu_utils/oclutils.h
src/gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.cu
src/gromacs/mdlib/nbnxn_ocl/nbnxn_ocl_data_mgmt.cpp

index 5d985c81fce2214c2d32411f441f6194f815df0b..d74ecd4d8dea253fe7c4663428294a5890a396bc 100644 (file)
@@ -297,4 +297,19 @@ static inline bool haveStreamTasksCompleted(cudaStream_t s)
     return true;
 }
 
+/*! \brief Free a device-side buffer.
+ * TODO: fully replace cu_free_buffered with this.
+ *
+ * \param[in] buffer  Pointer to the buffer to free.
+ */
+template <typename DeviceBuffer>
+void freeDeviceBuffer(DeviceBuffer *buffer)
+{
+    GMX_ASSERT(buffer, "needs a buffer pointer");
+    if (*buffer)
+    {
+        GMX_RELEASE_ASSERT(cudaFree(*buffer) == cudaSuccess, "cudaFree failed");
+    }
+}
+
 #endif
index 2b2c82a88aaee9c0587c313105f3d9fc247586b0..e4bc91b320430a01089f9eb854320126d8a92ee4 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2014,2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2014,2015,2016,2017,2018, 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.
@@ -177,4 +177,19 @@ static inline bool haveStreamTasksCompleted(cl_command_queue gmx_unused s)
     return false;
 }
 
+/*! \brief Free a device-side buffer.
+ * TODO: fully replace free_ocl_buffer and ocl_free_buffered with this.
+ *
+ * \param[in] buffer  Pointer to the buffer to free.
+ */
+template <typename DeviceBuffer>
+void freeDeviceBuffer(DeviceBuffer *buffer)
+{
+    GMX_ASSERT(buffer, "needs a buffer pointer");
+    if (*buffer)
+    {
+        GMX_RELEASE_ASSERT(clReleaseMemObject(*buffer) == CL_SUCCESS, "clReleaseMemObject failed");
+    }
+}
+
 #endif
index cc8ba63b0ebcb713eeb3bf843714f3734c56a00e..ea5375667d0c3783b65f48804d36937538572fcb 100644 (file)
@@ -734,7 +734,6 @@ void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
     cudaError_t      stat;
     cu_atomdata_t   *atdat;
     cu_nbparam_t    *nbparam;
-    cu_plist_t      *plist, *plist_nl;
 
     if (nb == NULL)
     {
@@ -743,8 +742,6 @@ void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
 
     atdat       = nb->atdat;
     nbparam     = nb->nbparam;
-    plist       = nb->plist[eintLocal];
-    plist_nl    = nb->plist[eintNonlocal];
 
     nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
 
@@ -792,25 +789,25 @@ void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
     cu_free_buffered(atdat->lj_comb);
 
-    cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
-    cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
-    cu_free_buffered(plist->imask, &plist->nimask, &plist->imask_nalloc);
-    cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
+    /* Free plist */
+    auto *plist = nb->plist[eintLocal];
+    freeDeviceBuffer(&plist->sci);
+    freeDeviceBuffer(&plist->cj4);
+    freeDeviceBuffer(&plist->imask);
+    freeDeviceBuffer(&plist->excl);
+    sfree(plist);
     if (nb->bUseTwoStreams)
     {
-        cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
-        cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
-        cu_free_buffered(plist_nl->imask, &plist_nl->nimask, &plist_nl->imask_nalloc);
-        cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
+        auto *plist_nl = nb->plist[eintNonlocal];
+        freeDeviceBuffer(&plist_nl->sci);
+        freeDeviceBuffer(&plist_nl->cj4);
+        freeDeviceBuffer(&plist_nl->imask);
+        freeDeviceBuffer(&plist_nl->excl);
+        sfree(plist_nl);
     }
 
     sfree(atdat);
     sfree(nbparam);
-    sfree(plist);
-    if (nb->bUseTwoStreams)
-    {
-        sfree(plist_nl);
-    }
     sfree(nb->timings);
     sfree(nb);
 
index e65160374f680152d246561ebd528bdba1b81db8..189c06944f064df9659a349f20d12e2088594198 100644 (file)
@@ -1180,18 +1180,20 @@ void nbnxn_gpu_free(gmx_nbnxn_ocl_t *nb)
     sfree(nb->nbparam);
 
     /* Free plist */
-    free_ocl_buffer(&(nb->plist[eintLocal]->sci));
-    free_ocl_buffer(&(nb->plist[eintLocal]->cj4));
-    free_ocl_buffer(&(nb->plist[eintLocal]->imask));
-    free_ocl_buffer(&(nb->plist[eintLocal]->excl));
-    sfree(nb->plist[eintLocal]);
+    auto *plist = nb->plist[eintLocal];
+    freeDeviceBuffer(&plist->sci);
+    freeDeviceBuffer(&plist->cj4);
+    freeDeviceBuffer(&plist->imask);
+    freeDeviceBuffer(&plist->excl);
+    sfree(plist);
     if (nb->bUseTwoStreams)
     {
-        free_ocl_buffer(&(nb->plist[eintNonlocal]->sci));
-        free_ocl_buffer(&(nb->plist[eintNonlocal]->cj4));
-        free_ocl_buffer(&(nb->plist[eintNonlocal]->imask));
-        free_ocl_buffer(&(nb->plist[eintNonlocal]->excl));
-        sfree(nb->plist[eintNonlocal]);
+        auto *plist_nl = nb->plist[eintNonlocal];
+        freeDeviceBuffer(&plist_nl->sci);
+        freeDeviceBuffer(&plist_nl->cj4);
+        freeDeviceBuffer(&plist_nl->imask);
+        freeDeviceBuffer(&plist_nl->excl);
+        sfree(plist_nl);
     }
 
     /* Free nbst */