Move CUDA texture setup code from NB CUDA module to cudautils.cu
authorAleksei Iupinov <a.yupinov@gmail.com>
Fri, 22 Sep 2017 12:41:02 +0000 (14:41 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 26 Sep 2017 10:59:33 +0000 (12:59 +0200)
Change-Id: I7e47a65866c29be06ce522572e90a17c775157ab

src/gromacs/gpu_utils/cuda_arch_utils.cuh
src/gromacs/gpu_utils/cudautils.cu
src/gromacs/gpu_utils/cudautils.cuh
src/gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_data_mgmt.cu
src/gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_types.h

index 7cceb1a04d96e5286fed9b796fa90ae2e3d68867..4639acdd2344db50b7858680349950a0185e0827 100644 (file)
@@ -150,6 +150,10 @@ T gmx_shfl_down_sync(const unsigned int activeMask,
 #define DISABLE_CUDA_TEXTURES 0
 #endif
 
+/*! \brief True if the use of texture fetch in the CUDA kernels is disabled. */
+static const bool c_disableCudaTextures = DISABLE_CUDA_TEXTURES;
+
+
 /* CUDA architecture technical characteristics. Needs macros because it is used
  * in the __launch_bounds__ function qualifiers and might need it in preprocessor
  * conditionals.
index 6021331ffb4740684113555b319579d683fda89a..75d87fa8ceae5f6b9f46a0d56d893ab196b6f4b8 100644 (file)
 
 #include "cudautils.cuh"
 
-#include <stdlib.h>
+#include <cassert>
+#include <cstdlib>
 
+#include "gromacs/gpu_utils/cuda_arch_utils.cuh"
 #include "gromacs/utility/smalloc.h"
 
 /*** Generic CUDA data operation wrappers ***/
@@ -246,3 +248,97 @@ void cu_realloc_buffered(void **d_dest, void *h_src,
         }
     }
 }
+
+bool use_texobj(const gmx_device_info_t *dev_info)
+{
+    assert(!c_disableCudaTextures);
+    /* Only device CC >= 3.0 (Kepler and later) support texture objects */
+    return (dev_info->prop.major >= 3);
+}
+
+/*! \brief Set up texture object for an array of type T.
+ *
+ * Set up texture object for an array of type T and bind it to the device memory
+ * \p d_ptr points to.
+ *
+ * \tparam[in] T        Raw data type
+ * \param[out] texObj   texture object to initialize
+ * \param[in]  d_ptr    pointer to device global memory to bind \p texObj to
+ * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
+ */
+template <typename T>
+static void setup1DTexture(cudaTextureObject_t &texObj,
+                           void                *d_ptr,
+                           size_t               sizeInBytes)
+{
+    assert(!c_disableCudaTextures);
+
+    cudaError_t      stat;
+    cudaResourceDesc rd;
+    cudaTextureDesc  td;
+
+    memset(&rd, 0, sizeof(rd));
+    rd.resType                = cudaResourceTypeLinear;
+    rd.res.linear.devPtr      = d_ptr;
+    rd.res.linear.desc        = cudaCreateChannelDesc<T>();
+    rd.res.linear.sizeInBytes = sizeInBytes;
+
+    memset(&td, 0, sizeof(td));
+    td.readMode                 = cudaReadModeElementType;
+    stat = cudaCreateTextureObject(&texObj, &rd, &td, NULL);
+    CU_RET_ERR(stat, "cudaCreateTextureObject failed");
+}
+
+/*! \brief Set up texture reference for an array of type T.
+ *
+ * Set up texture object for an array of type T and bind it to the device memory
+ * \p d_ptr points to.
+ *
+ * \tparam[in] T        Raw data type
+ * \param[out] texObj   texture reference to initialize
+ * \param[in]  d_ptr    pointer to device global memory to bind \p texObj to
+ * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
+ */
+template <typename T>
+static void setup1DTexture(const struct texture<T, 1, cudaReadModeElementType> *texRef,
+                           const void                                          *d_ptr,
+                           size_t                                              sizeInBytes)
+{
+    assert(!c_disableCudaTextures);
+
+    cudaError_t           stat;
+    cudaChannelFormatDesc cd;
+
+    cd   = cudaCreateChannelDesc<T>();
+    stat = cudaBindTexture(nullptr, texRef, d_ptr, &cd, sizeInBytes);
+    CU_RET_ERR(stat, "cudaBindTexture failed");
+}
+
+template <typename T>
+void initParamLookupTable(T                        * &d_ptr,
+                          cudaTextureObject_t       &texObj,
+                          const struct texture<T, 1, cudaReadModeElementType> *texRef,
+                          const T                   *h_ptr,
+                          int                        numElem,
+                          const gmx_device_info_t   *devInfo)
+{
+    const size_t sizeInBytes = numElem * sizeof(*d_ptr);
+    cudaError_t  stat        = cudaMalloc((void **)&d_ptr, sizeInBytes);
+    CU_RET_ERR(stat, "cudaMalloc failed in initParamLookupTable");
+    cu_copy_H2D(d_ptr, (void *)h_ptr, sizeInBytes);
+
+    if (!c_disableCudaTextures)
+    {
+        if (use_texobj(devInfo))
+        {
+            setup1DTexture<T>(texObj, d_ptr, sizeInBytes);
+        }
+        else
+        {
+            setup1DTexture<T>(texRef, d_ptr, sizeInBytes);
+        }
+    }
+}
+
+//! Add explicit instantiations of initParamLookupTable() here as needed
+template void initParamLookupTable<float>(float * &, cudaTextureObject_t &, const texture<float, 1, cudaReadModeElementType> *, const float *, int, const gmx_device_info_t *);
index 604e54b42cad81dd0dfad9033ce3193b475650b9..cfe5c2a381506ae3febaec07f84956bafcd1cd42 100644 (file)
@@ -165,4 +165,34 @@ float cu_event_elapsed(cudaEvent_t /*start*/, cudaEvent_t /*end*/);
 /*! Waits for event end to complete and calculates the time between start and end. */
 int cu_wait_event_time(cudaEvent_t /*end*/, cudaEvent_t /*begin*/, float * /*time*/);
 
+/*! \brief Return whether texture objects are used on this device.
+ *
+ * \todo This should be static in cudautils.cu, as soon as texture destruction code is moved there as well
+ *
+ * \param[in]   pointer to the GPU device info structure to inspect for texture objects support
+ * \return      true if texture objects are used on this device
+ */
+bool use_texobj(const gmx_device_info_t *dev_info);
+
+/*! \brief Initialize parameter lookup table.
+ *
+ * Initializes device memory, copies data from host and binds
+ * a texture to allocated device memory to be used for parameter lookup.
+ *
+ * \tparam[in] T         Raw data type
+ * \param[out] d_ptr     device pointer to the memory to be allocated
+ * \param[out] texObj    texture object to be initialized
+ * \param[out] texRef    texture reference to be initialized
+ * \param[in]  h_ptr     pointer to the host memory to be uploaded to the device
+ * \param[in]  numElem   number of elements in the h_ptr
+ * \param[in]  devInfo   pointer to the info struct of the device in use
+ */
+template <typename T>
+void initParamLookupTable(T                        * &d_ptr,
+                          cudaTextureObject_t       &texObj,
+                          const struct texture<T, 1, cudaReadModeElementType> *texRef,
+                          const T                   *h_ptr,
+                          int                        numElem,
+                          const gmx_device_info_t   *devInfo);
+
 #endif
index 79ef8e842e53f188aab4b5ecab761094e24907a6..24f316d5498c1d2162f8e47f57eb1d89f95a94bf 100644 (file)
@@ -86,19 +86,6 @@ static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
                                           const gmx_device_info_t *dev_info);
 
-
-/*! \brief Return whether texture objects are used on this device.
- *
- * \param[in]   pointer to the GPU device info structure to inspect for texture objects support
- * \return      true if texture objects are used on this device
- */
-static bool use_texobj(const gmx_device_info_t *dev_info)
-{
-    assert(!c_disableCudaTextures);
-    /* Only device CC >= 3.0 (Kepler and later) support texture objects */
-    return (dev_info->prop.major >= 3);
-}
-
 /*! \brief Return whether combination rules are used.
  *
  * \param[in]   pointer to nonbonded paramter struct
@@ -110,104 +97,6 @@ static inline bool useLjCombRule(const cu_nbparam_t  *nbparam)
             nbparam->vdwtype == evdwCuCUTCOMBLB);
 }
 
-/*! \brief Set up texture object for an array of type T.
- *
- * Set up texture object for an array of type T and bind it to the device memory
- * \p d_ptr points to.
- *
- * \tparam[in] T        Raw data type
- * \param[out] texObj   texture object to initialize
- * \param[in]  d_ptr    pointer to device global memory to bind \p texObj to
- * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
- */
-template <typename T>
-static void setup1DTexture(cudaTextureObject_t &texObj,
-                           void                *d_ptr,
-                           size_t               sizeInBytes)
-{
-    assert(!c_disableCudaTextures);
-
-    cudaError_t      stat;
-    cudaResourceDesc rd;
-    cudaTextureDesc  td;
-
-    memset(&rd, 0, sizeof(rd));
-    rd.resType                = cudaResourceTypeLinear;
-    rd.res.linear.devPtr      = d_ptr;
-    rd.res.linear.desc        = cudaCreateChannelDesc<T>();
-    rd.res.linear.sizeInBytes = sizeInBytes;
-
-    memset(&td, 0, sizeof(td));
-    td.readMode                 = cudaReadModeElementType;
-    stat = cudaCreateTextureObject(&texObj, &rd, &td, NULL);
-    CU_RET_ERR(stat, "cudaCreateTextureObject failed");
-}
-
-/*! \brief Set up texture reference for an array of type T.
- *
- * Set up texture object for an array of type T and bind it to the device memory
- * \p d_ptr points to.
- *
- * \tparam[in] T        Raw data type
- * \param[out] texObj   texture reference to initialize
- * \param[in]  d_ptr    pointer to device global memory to bind \p texObj to
- * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
- */
-template <typename T>
-static void setup1DTexture(const struct texture<T, 1, cudaReadModeElementType> *texRef,
-                           const void                                          *d_ptr,
-                           size_t                                              sizeInBytes)
-{
-    assert(!c_disableCudaTextures);
-
-    cudaError_t           stat;
-    cudaChannelFormatDesc cd;
-
-    cd   = cudaCreateChannelDesc<T>();
-    stat = cudaBindTexture(nullptr, texRef, d_ptr, &cd, sizeInBytes);
-    CU_RET_ERR(stat, "cudaBindTexture failed");
-}
-
-/*! \brief Initialize parameter lookup table.
- *
- * Initializes device memory, copies data from host and binds
- * a texture to allocated device memory to be used for LJ/Ewald/... parameter
- * lookup.
- *
- * \tparam[in] T         Raw data type
- * \param[out] d_ptr     device pointer to the memory to be allocated
- * \param[out] texObj    texture object to be initialized
- * \param[out] texRef    texture reference to be initialized
- * \param[in]  h_ptr     pointer to the host memory to be uploaded to the device
- * \param[in]  numElem   number of elements in the h_ptr
- * \param[in]  devInfo   pointer to the info struct of the device in use
- */
-template <typename T>
-static void initParamLookupTable(T                        * &d_ptr,
-                                 cudaTextureObject_t       &texObj,
-                                 const struct texture<T, 1, cudaReadModeElementType> *texRef,
-                                 const T                   *h_ptr,
-                                 int                        numElem,
-                                 const gmx_device_info_t   *devInfo)
-{
-    const size_t sizeInBytes = numElem * sizeof(*d_ptr);
-    cudaError_t  stat        = cudaMalloc((void **)&d_ptr, sizeInBytes);
-    CU_RET_ERR(stat, "cudaMalloc failed in initParamLookupTable");
-    cu_copy_H2D(d_ptr, (void *)h_ptr, sizeInBytes);
-
-    if (!c_disableCudaTextures)
-    {
-        if (use_texobj(devInfo))
-        {
-            setup1DTexture<T>(texObj, d_ptr, sizeInBytes);
-        }
-        else
-        {
-            setup1DTexture<T>(texRef, d_ptr, sizeInBytes);
-        }
-    }
-}
-
 /*! \brief Initialized the Ewald Coulomb correction GPU table.
 
     Tabulates the Ewald Coulomb force and initializes the size/scale
index 609123f346f784ebb27d7eaae10b0808b87846df..bbd283d075429cec6f75d595ae12f5e8631e2262 100644 (file)
@@ -74,10 +74,6 @@ static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
 /*! \brief cluster size = number of atoms per cluster. */
 static const int c_clSize          = c_nbnxnGpuClusterSize;
 
-/*! \brief True if the use of texture fetch in the CUDA kernels is disabled. */
-static const bool c_disableCudaTextures = DISABLE_CUDA_TEXTURES;
-
-
 #ifdef __cplusplus
 extern "C" {
 #endif