Merge branch release-2018
[alexxy/gromacs.git] / src / gromacs / gpu_utils / cudautils.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #ifndef GMX_GPU_UTILS_CUDAUTILS_CUH
36 #define GMX_GPU_UTILS_CUDAUTILS_CUH
37
38 #include "config.h"
39
40 #include <stdio.h>
41 #if HAVE_NVML
42 #include <nvml.h>
43 #endif /* HAVE_NVML */
44
45 #include <array>
46 #include <string>
47
48 #include "gromacs/gpu_utils/gputraits.cuh"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/utility/exceptions.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/stringutil.h"
55
56 namespace gmx
57 {
58 namespace
59 {
60
61 /*! \brief Helper function to ensure no pending error silently
62  * disrupts error handling.
63  *
64  * Asserts in a debug build if an unhandled error is present. Issues a
65  * warning at run time otherwise.
66  *
67  * \todo This is similar to CU_CHECK_PREV_ERR, which should be
68  * consolidated.
69  */
70 static inline void ensureNoPendingCudaError(const char *errorMessage)
71 {
72     // Ensure there is no pending error that would otherwise affect
73     // the behaviour of future error handling.
74     cudaError_t stat = cudaGetLastError();
75     if (stat == cudaSuccess)
76     {
77         return;
78     }
79
80     // If we would find an error in a release build, we do not know
81     // what is appropriate to do about it, so assert only for debug
82     // builds.
83     auto fullMessage = formatString("%s An unhandled error from a previous CUDA operation was detected. %s: %s",
84                                     errorMessage, cudaGetErrorName(stat), cudaGetErrorString(stat));
85     GMX_ASSERT(stat == cudaSuccess, fullMessage.c_str());
86     // TODO When we evolve a better logging framework, use that
87     // for release-build error reporting.
88     gmx_warning(fullMessage.c_str());
89 }
90
91 }   // namespace
92 }   // namespace
93
94 enum class GpuApiCallBehavior;
95
96 /* TODO error checking needs to be rewritten. We have 2 types of error checks needed
97    based on where they occur in the code:
98    - non performance-critical: these errors are unsafe to be ignored and must be
99      _always_ checked for, e.g. initializations
100    - performance critical: handling errors might hurt performance so care need to be taken
101      when/if we should check for them at all, e.g. in cu_upload_X. However, we should be
102      able to turn the check for these errors on!
103
104    Probably we'll need two sets of the macros below...
105
106  */
107 #define CHECK_CUDA_ERRORS
108
109 #ifdef CHECK_CUDA_ERRORS
110
111 /*! Check for CUDA error on the return status of a CUDA RT API call. */
112 #define CU_RET_ERR(status, msg) \
113     do { \
114         if (status != cudaSuccess) \
115         { \
116             gmx_fatal(FARGS, "%s: %s\n", msg, cudaGetErrorString(status)); \
117         } \
118     } while (0)
119
120 /*! Check for any previously occurred uncaught CUDA error. */
121 #define CU_CHECK_PREV_ERR() \
122     do { \
123         cudaError_t _CU_CHECK_PREV_ERR_status = cudaGetLastError(); \
124         if (_CU_CHECK_PREV_ERR_status != cudaSuccess) { \
125             gmx_warning("Just caught a previously occurred CUDA error (%s), will try to continue.", cudaGetErrorString(_CU_CHECK_PREV_ERR_status)); \
126         } \
127     } while (0)
128
129 #else /* CHECK_CUDA_ERRORS */
130
131 #define CU_RET_ERR(status, msg) do { } while (0)
132 #define CU_CHECK_PREV_ERR()     do { } while (0)
133 #define HANDLE_NVML_RET_ERR(status, msg) do { } while (0)
134
135 #endif /* CHECK_CUDA_ERRORS */
136
137 /*! \brief CUDA device information.
138  *
139  * The CUDA device information is queried and set at detection and contains
140  * both information about the device/hardware returned by the runtime as well
141  * as additional data like support status.
142  *
143  * \todo extract an object to manage NVML details
144  */
145 struct gmx_device_info_t
146 {
147     int                 id;                      /* id of the CUDA device */
148     cudaDeviceProp      prop;                    /* CUDA device properties */
149     int                 stat;                    /* result of the device check */
150     unsigned int        nvml_orig_app_sm_clock;  /* The original SM clock before we changed it */
151     unsigned int        nvml_orig_app_mem_clock; /* The original memory clock before we changed it */
152     gmx_bool            nvml_app_clocks_changed; /* If application clocks have been changed */
153     unsigned int        nvml_set_app_sm_clock;   /* The SM clock we set */
154     unsigned int        nvml_set_app_mem_clock;  /* The memory clock we set */
155 #if HAVE_NVML
156     nvmlDevice_t        nvml_device_id;          /* NVML device id */
157     // TODO This can become a bool with a more useful name
158     nvmlEnableState_t   nvml_is_restricted;      /* Status of application clocks permission */
159 #endif                                           /* HAVE_NVML */
160 };
161
162 /*! Launches synchronous or asynchronous device to host memory copy.
163  *
164  *  The copy is launched in stream s or if not specified, in stream 0.
165  */
166 int cu_copy_D2H(void *h_dest, void *d_src, size_t bytes, GpuApiCallBehavior transferKind, cudaStream_t s /*= 0*/);
167
168 /*! Launches synchronous host to device memory copy in stream 0. */
169 int cu_copy_D2H_sync(void * /*h_dest*/, void * /*d_src*/, size_t /*bytes*/);
170
171 /*! Launches asynchronous host to device memory copy in stream s. */
172 int cu_copy_D2H_async(void * /*h_dest*/, void * /*d_src*/, size_t /*bytes*/, cudaStream_t /*s = 0*/);
173
174 /*! Launches synchronous or asynchronous host to device memory copy.
175  *
176  *  The copy is launched in stream s or if not specified, in stream 0.
177  */
178 int cu_copy_H2D(void *d_dest, void *h_src, size_t bytes, GpuApiCallBehavior transferKind, cudaStream_t /*s = 0*/);
179
180 /*! Launches synchronous host to device memory copy. */
181 int cu_copy_H2D_sync(void * /*d_dest*/, void * /*h_src*/, size_t /*bytes*/);
182
183 /*! Launches asynchronous host to device memory copy in stream s. */
184 int cu_copy_H2D_async(void * /*d_dest*/, void * /*h_src*/, size_t /*bytes*/, cudaStream_t /*s = 0*/);
185
186 // TODO: the 2 functions below are pretty much a constructor/destructor of a simple
187 // GPU table object. There is also almost self-contained fetchFromParamLookupTable()
188 // in cuda_kernel_utils.cuh. They could all live in a separate class/struct file.
189
190 /*! \brief Initialize parameter lookup table.
191  *
192  * Initializes device memory, copies data from host and binds
193  * a texture to allocated device memory to be used for parameter lookup.
194  *
195  * \tparam[in] T         Raw data type
196  * \param[out] d_ptr     device pointer to the memory to be allocated
197  * \param[out] texObj    texture object to be initialized
198  * \param[in]  h_ptr     pointer to the host memory to be uploaded to the device
199  * \param[in]  numElem   number of elements in the h_ptr
200  * \param[in]  devInfo   pointer to the info struct of the device in use
201  */
202 template <typename T>
203 void initParamLookupTable(T                        * &d_ptr,
204                           cudaTextureObject_t        &texObj,
205                           const T                    *h_ptr,
206                           int                         numElem,
207                           const gmx_device_info_t    *devInfo);
208
209 /*! \brief Destroy parameter lookup table.
210  *
211  * Unbinds texture object, deallocates device memory.
212  *
213  * \tparam[in] T         Raw data type
214  * \param[in]  d_ptr     Device pointer to the memory to be deallocated
215  * \param[in]  texObj    Texture object to be deinitialized
216  * \param[in]  devInfo   Pointer to the info struct of the device in use
217  */
218 template <typename T>
219 void destroyParamLookupTable(T                       *d_ptr,
220                              cudaTextureObject_t      texObj,
221                              const gmx_device_info_t *devInfo);
222
223 /*! \brief Add a triplets stored in a float3 to an rvec variable.
224  *
225  * \param[out]  a Rvec to increment
226  * \param[in]   b Float triplet to increment with.
227  */
228 static inline void rvec_inc(rvec a, const float3 b)
229 {
230     rvec tmp = {b.x, b.y, b.z};
231     rvec_inc(a, tmp);
232 }
233
234 /*! \brief Wait for all taks in stream \p s to complete.
235  *
236  * \param[in] s stream to synchronize with
237  */
238 static inline void gpuStreamSynchronize(cudaStream_t s)
239 {
240     cudaError_t stat = cudaStreamSynchronize(s);
241     CU_RET_ERR(stat, "cudaStreamSynchronize failed");
242 }
243
244 /*! \brief  Returns true if all tasks in \p s have completed.
245  *
246  * \param[in] s stream to check
247  *
248  *  \returns     True if all tasks enqueued in the stream \p s (at the time of this call) have completed.
249  */
250 static inline bool haveStreamTasksCompleted(cudaStream_t s)
251 {
252     cudaError_t stat = cudaStreamQuery(s);
253
254     if (stat == cudaErrorNotReady)
255     {
256         // work is still in progress in the stream
257         return false;
258     }
259
260     GMX_ASSERT(stat !=  cudaErrorInvalidResourceHandle, "Stream idnetifier not valid");
261
262     // cudaSuccess and cudaErrorNotReady are the expected return values
263     CU_RET_ERR(stat, "Unexpected cudaStreamQuery failure");
264
265     GMX_ASSERT(stat == cudaSuccess, "Values other than cudaSuccess should have been explicitly handled");
266
267     return true;
268 }
269
270 /* Kernel launch helpers */
271
272 /*! \brief
273  * A function for setting up a single CUDA kernel argument.
274  * This is the tail of the compile-time recursive function below.
275  * It has to be seen by the compiler first.
276  *
277  * \tparam        totalArgsCount  Number of the kernel arguments
278  * \tparam        KernelPtr       Kernel function handle type
279  * \param[in]     argIndex        Index of the current argument
280  */
281 template <size_t totalArgsCount, typename KernelPtr>
282 void prepareGpuKernelArgument(KernelPtr /*kernel*/,
283                               std::array<void *, totalArgsCount> */* kernelArgsPtr */,
284                               size_t gmx_used_in_debug argIndex)
285 {
286     GMX_ASSERT(argIndex == totalArgsCount, "Tail expansion");
287 }
288
289 /*! \brief
290  * Compile-time recursive function for setting up a single CUDA kernel argument.
291  * This function copies a kernel argument pointer \p argPtr into \p kernelArgsPtr,
292  * and calls itself on the next argument, eventually calling the tail function above.
293  *
294  * \tparam        CurrentArg      Type of the current argument
295  * \tparam        RemainingArgs   Types of remaining arguments after the current one
296  * \tparam        totalArgsCount  Number of the kernel arguments
297  * \tparam        KernelPtr       Kernel function handle type
298  * \param[in]     kernel          Kernel function handle
299  * \param[in,out] kernelArgsPtr   Pointer to the argument array to be filled in
300  * \param[in]     argIndex        Index of the current argument
301  * \param[in]     argPtr          Pointer to the current argument
302  * \param[in]     otherArgsPtrs   Pack of pointers to arguments remaining to process after the current one
303  */
304 template <typename CurrentArg, typename ... RemainingArgs, size_t totalArgsCount, typename KernelPtr>
305 void prepareGpuKernelArgument(KernelPtr kernel,
306                               std::array<void *, totalArgsCount> *kernelArgsPtr,
307                               size_t argIndex,
308                               const CurrentArg *argPtr,
309                               const RemainingArgs *... otherArgsPtrs)
310 {
311     (*kernelArgsPtr)[argIndex] = (void *)argPtr;
312     prepareGpuKernelArgument(kernel, kernelArgsPtr, argIndex + 1, otherArgsPtrs ...);
313 }
314
315 /*! \brief
316  * A wrapper function for setting up all the CUDA kernel arguments.
317  * Calls the recursive functions above.
318  *
319  * \tparam    Args            Types of all the kernel arguments
320  * \param[in] kernel          Kernel function handle
321  * \param[in] argsPtrs        Pointers to all the kernel arguments
322  * \returns A prepared parameter pack to be used with launchGpuKernel() as the last argument.
323  */
324 template <typename ... Args>
325 std::array<void *, sizeof ... (Args)> prepareGpuKernelArguments(void                     (*kernel)(Args...),
326                                                                 const KernelLaunchConfig & /*config */,
327                                                                 const Args *...          argsPtrs)
328 {
329     std::array<void *, sizeof ... (Args)> kernelArgs;
330     prepareGpuKernelArgument(kernel, &kernelArgs, 0, argsPtrs ...);
331     return kernelArgs;
332 }
333
334 /*! \brief Launches the CUDA kernel and handles the errors.
335  *
336  * \tparam    Args            Types of all the kernel arguments
337  * \param[in] kernel          Kernel function handle
338  * \param[in] config          Kernel configuration for launching
339  * \param[in] kernelName      Human readable kernel description, for error handling only
340  * \param[in] kernelArgs      Array of the pointers to the kernel arguments, prepared by prepareGpuKernelArguments()
341  * \throws gmx::InternalError on kernel launch failure
342  */
343 template <typename... Args>
344 void launchGpuKernel(void                                       (*kernel)(Args...),
345                      const KernelLaunchConfig                    &config,
346                      CommandEvent                                */*timingEvent */,
347                      const char                                  *kernelName,
348                      const std::array<void *, sizeof ... (Args)> &kernelArgs)
349 {
350     dim3 blockSize(config.blockSize[0], config.blockSize[1], config.blockSize[2]);
351     dim3 gridSize(config.gridSize[0], config.gridSize[1], config.gridSize[2]);
352     cudaLaunchKernel((void *)kernel, gridSize, blockSize, const_cast<void **>(kernelArgs.data()), config.sharedMemorySize, config.stream);
353
354     cudaError_t status = cudaGetLastError();
355     if (cudaSuccess != status)
356     {
357         const std::string errorMessage = "GPU kernel (" +  std::string(kernelName) +
358             ") failed to launch: " + std::string(cudaGetErrorString(status));
359         GMX_THROW(gmx::InternalError(errorMessage));
360     }
361 }
362
363 #endif