Apply clang-format to source tree
[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,2019, 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 <stdio.h>
39
40 #include <array>
41 #include <string>
42
43 #include "gromacs/gpu_utils/gputraits.cuh"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/utility/exceptions.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/gmxassert.h"
49 #include "gromacs/utility/stringutil.h"
50
51 namespace gmx
52 {
53 namespace
54 {
55
56 /*! \brief Helper function to ensure no pending error silently
57  * disrupts error handling.
58  *
59  * Asserts in a debug build if an unhandled error is present. Issues a
60  * warning at run time otherwise.
61  *
62  * \todo This is similar to CU_CHECK_PREV_ERR, which should be
63  * consolidated.
64  */
65 static inline void ensureNoPendingCudaError(const char* errorMessage)
66 {
67     // Ensure there is no pending error that would otherwise affect
68     // the behaviour of future error handling.
69     cudaError_t stat = cudaGetLastError();
70     if (stat == cudaSuccess)
71     {
72         return;
73     }
74
75     // If we would find an error in a release build, we do not know
76     // what is appropriate to do about it, so assert only for debug
77     // builds.
78     auto fullMessage = formatString(
79             "%s An unhandled error from a previous CUDA operation was detected. %s: %s",
80             errorMessage, cudaGetErrorName(stat), cudaGetErrorString(stat));
81     GMX_ASSERT(stat == cudaSuccess, fullMessage.c_str());
82     // TODO When we evolve a better logging framework, use that
83     // for release-build error reporting.
84     gmx_warning("%s", fullMessage.c_str());
85 }
86
87 } // namespace
88 } // namespace gmx
89
90 enum class GpuApiCallBehavior;
91
92 /* TODO error checking needs to be rewritten. We have 2 types of error checks needed
93    based on where they occur in the code:
94    - non performance-critical: these errors are unsafe to be ignored and must be
95      _always_ checked for, e.g. initializations
96    - performance critical: handling errors might hurt performance so care need to be taken
97      when/if we should check for them at all, e.g. in cu_upload_X. However, we should be
98      able to turn the check for these errors on!
99
100    Probably we'll need two sets of the macros below...
101
102  */
103 #define CHECK_CUDA_ERRORS
104
105 #ifdef CHECK_CUDA_ERRORS
106
107 /*! Check for CUDA error on the return status of a CUDA RT API call. */
108 #    define CU_RET_ERR(status, msg)                                            \
109         do                                                                     \
110         {                                                                      \
111             if (status != cudaSuccess)                                         \
112             {                                                                  \
113                 gmx_fatal(FARGS, "%s: %s\n", msg, cudaGetErrorString(status)); \
114             }                                                                  \
115         } while (0)
116
117 /*! Check for any previously occurred uncaught CUDA error. */
118 #    define CU_CHECK_PREV_ERR()                                                           \
119         do                                                                                \
120         {                                                                                 \
121             cudaError_t _CU_CHECK_PREV_ERR_status = cudaGetLastError();                   \
122             if (_CU_CHECK_PREV_ERR_status != cudaSuccess)                                 \
123             {                                                                             \
124                 gmx_warning(                                                              \
125                         "Just caught a previously occurred CUDA error (%s), will try to " \
126                         "continue.",                                                      \
127                         cudaGetErrorString(_CU_CHECK_PREV_ERR_status));                   \
128             }                                                                             \
129         } while (0)
130
131 #else /* CHECK_CUDA_ERRORS */
132
133 #    define CU_RET_ERR(status, msg) \
134         do                          \
135         {                           \
136         } while (0)
137 #    define CU_CHECK_PREV_ERR() \
138         do                      \
139         {                       \
140         } while (0)
141
142 #endif /* CHECK_CUDA_ERRORS */
143
144 /*! \brief CUDA device information.
145  *
146  * The CUDA device information is queried and set at detection and contains
147  * both information about the device/hardware returned by the runtime as well
148  * as additional data like support status.
149  */
150 struct gmx_device_info_t
151 {
152     int            id;   /* id of the CUDA device */
153     cudaDeviceProp prop; /* CUDA device properties */
154     int            stat; /* result of the device check */
155 };
156
157 /*! Launches synchronous or asynchronous device to host memory copy.
158  *
159  *  The copy is launched in stream s or if not specified, in stream 0.
160  */
161 int cu_copy_D2H(void* h_dest, void* d_src, size_t bytes, GpuApiCallBehavior transferKind, cudaStream_t /*s = nullptr*/);
162
163 /*! Launches synchronous host to device memory copy in stream 0. */
164 int cu_copy_D2H_sync(void* /*h_dest*/, void* /*d_src*/, size_t /*bytes*/);
165
166 /*! Launches asynchronous host to device memory copy in stream s. */
167 int cu_copy_D2H_async(void* /*h_dest*/, void* /*d_src*/, size_t /*bytes*/, cudaStream_t /*s = nullptr*/);
168
169 /*! Launches synchronous or asynchronous host to device memory copy.
170  *
171  *  The copy is launched in stream s or if not specified, in stream 0.
172  */
173 int cu_copy_H2D(void* d_dest, const void* h_src, size_t bytes, GpuApiCallBehavior transferKind, cudaStream_t /*s = nullptr*/);
174
175 /*! Launches synchronous host to device memory copy. */
176 int cu_copy_H2D_sync(void* /*d_dest*/, const void* /*h_src*/, size_t /*bytes*/);
177
178 /*! Launches asynchronous host to device memory copy in stream s. */
179 int cu_copy_H2D_async(void* /*d_dest*/, const void* /*h_src*/, size_t /*bytes*/, cudaStream_t /*s = nullptr*/);
180
181 // TODO: the 2 functions below are pretty much a constructor/destructor of a simple
182 // GPU table object. There is also almost self-contained fetchFromParamLookupTable()
183 // in cuda_kernel_utils.cuh. They could all live in a separate class/struct file.
184
185 /*! \brief Initialize parameter lookup table.
186  *
187  * Initializes device memory, copies data from host and binds
188  * a texture to allocated device memory to be used for parameter lookup.
189  *
190  * \tparam[in] T         Raw data type
191  * \param[out] d_ptr     device pointer to the memory to be allocated
192  * \param[out] texObj    texture object to be initialized
193  * \param[in]  h_ptr     pointer to the host memory to be uploaded to the device
194  * \param[in]  numElem   number of elements in the h_ptr
195  */
196 template<typename T>
197 void initParamLookupTable(T*& d_ptr, cudaTextureObject_t& texObj, const T* h_ptr, int numElem);
198
199 // Add extern declarations so each translation unit understands that
200 // there will be a definition provided.
201 extern template void initParamLookupTable<int>(int*&, cudaTextureObject_t&, const int*, int);
202 extern template void initParamLookupTable<float>(float*&, cudaTextureObject_t&, const float*, int);
203
204 /*! \brief Destroy parameter lookup table.
205  *
206  * Unbinds texture object, deallocates device memory.
207  *
208  * \tparam[in] T         Raw data type
209  * \param[in]  d_ptr     Device pointer to the memory to be deallocated
210  * \param[in]  texObj    Texture object to be deinitialized
211  */
212 template<typename T>
213 void destroyParamLookupTable(T* d_ptr, cudaTextureObject_t texObj);
214
215 // Add extern declarations so each translation unit understands that
216 // there will be a definition provided.
217 extern template void destroyParamLookupTable<int>(int*, cudaTextureObject_t);
218 extern template void destroyParamLookupTable<float>(float*, cudaTextureObject_t);
219
220 /*! \brief Add a triplets stored in a float3 to an rvec variable.
221  *
222  * \param[out]  a Rvec to increment
223  * \param[in]   b Float triplet to increment with.
224  */
225 static inline void rvec_inc(rvec a, const float3 b)
226 {
227     rvec tmp = { b.x, b.y, b.z };
228     rvec_inc(a, tmp);
229 }
230
231 /*! \brief Wait for all taks in stream \p s to complete.
232  *
233  * \param[in] s stream to synchronize with
234  */
235 static inline void gpuStreamSynchronize(cudaStream_t s)
236 {
237     cudaError_t stat = cudaStreamSynchronize(s);
238     CU_RET_ERR(stat, "cudaStreamSynchronize failed");
239 }
240
241 /*! \brief  Returns true if all tasks in \p s have completed.
242  *
243  * \param[in] s stream to check
244  *
245  *  \returns     True if all tasks enqueued in the stream \p s (at the time of this call) have completed.
246  */
247 static inline bool haveStreamTasksCompleted(cudaStream_t s)
248 {
249     cudaError_t stat = cudaStreamQuery(s);
250
251     if (stat == cudaErrorNotReady)
252     {
253         // work is still in progress in the stream
254         return false;
255     }
256
257     GMX_ASSERT(stat != cudaErrorInvalidResourceHandle, "Stream idnetifier not valid");
258
259     // cudaSuccess and cudaErrorNotReady are the expected return values
260     CU_RET_ERR(stat, "Unexpected cudaStreamQuery failure");
261
262     GMX_ASSERT(stat == cudaSuccess,
263                "Values other than cudaSuccess should have been explicitly handled");
264
265     return true;
266 }
267
268 /* Kernel launch helpers */
269
270 /*! \brief
271  * A function for setting up a single CUDA kernel argument.
272  * This is the tail of the compile-time recursive function below.
273  * It has to be seen by the compiler first.
274  *
275  * \tparam        totalArgsCount  Number of the kernel arguments
276  * \tparam        KernelPtr       Kernel function handle type
277  * \param[in]     argIndex        Index of the current argument
278  */
279 template<size_t totalArgsCount, typename KernelPtr>
280 void prepareGpuKernelArgument(KernelPtr /*kernel*/,
281                               std::array<void*, totalArgsCount>* /* kernelArgsPtr */,
282                               size_t gmx_used_in_debug argIndex)
283 {
284     GMX_ASSERT(argIndex == totalArgsCount, "Tail expansion");
285 }
286
287 /*! \brief
288  * Compile-time recursive function for setting up a single CUDA kernel argument.
289  * This function copies a kernel argument pointer \p argPtr into \p kernelArgsPtr,
290  * and calls itself on the next argument, eventually calling the tail function above.
291  *
292  * \tparam        CurrentArg      Type of the current argument
293  * \tparam        RemainingArgs   Types of remaining arguments after the current one
294  * \tparam        totalArgsCount  Number of the kernel arguments
295  * \tparam        KernelPtr       Kernel function handle type
296  * \param[in]     kernel          Kernel function handle
297  * \param[in,out] kernelArgsPtr   Pointer to the argument array to be filled in
298  * \param[in]     argIndex        Index of the current argument
299  * \param[in]     argPtr          Pointer to the current argument
300  * \param[in]     otherArgsPtrs   Pack of pointers to arguments remaining to process after the current one
301  */
302 template<typename CurrentArg, typename... RemainingArgs, size_t totalArgsCount, typename KernelPtr>
303 void prepareGpuKernelArgument(KernelPtr                          kernel,
304                               std::array<void*, totalArgsCount>* kernelArgsPtr,
305                               size_t                             argIndex,
306                               const CurrentArg*                  argPtr,
307                               const RemainingArgs*... otherArgsPtrs)
308 {
309     (*kernelArgsPtr)[argIndex] = (void*)argPtr;
310     prepareGpuKernelArgument(kernel, kernelArgsPtr, argIndex + 1, otherArgsPtrs...);
311 }
312
313 /*! \brief
314  * A wrapper function for setting up all the CUDA kernel arguments.
315  * Calls the recursive functions above.
316  *
317  * \tparam    KernelPtr       Kernel function handle type
318  * \tparam    Args            Types of all the kernel arguments
319  * \param[in] kernel          Kernel function handle
320  * \param[in] argsPtrs        Pointers to all the kernel arguments
321  * \returns A prepared parameter pack to be used with launchGpuKernel() as the last argument.
322  */
323 template<typename KernelPtr, typename... Args>
324 std::array<void*, sizeof...(Args)> prepareGpuKernelArguments(KernelPtr kernel,
325                                                              const KernelLaunchConfig& /*config */,
326                                                              const Args*... argsPtrs)
327 {
328     std::array<void*, sizeof...(Args)> kernelArgs;
329     prepareGpuKernelArgument(kernel, &kernelArgs, 0, argsPtrs...);
330     return kernelArgs;
331 }
332
333 /*! \brief Launches the CUDA kernel and handles the errors.
334  *
335  * \tparam    Args            Types of all the kernel arguments
336  * \param[in] kernel          Kernel function handle
337  * \param[in] config          Kernel configuration for launching
338  * \param[in] kernelName      Human readable kernel description, for error handling only
339  * \param[in] kernelArgs      Array of the pointers to the kernel arguments, prepared by
340  * prepareGpuKernelArguments() \throws gmx::InternalError on kernel launch failure
341  */
342 template<typename... Args>
343 void launchGpuKernel(void (*kernel)(Args...),
344                      const KernelLaunchConfig& config,
345                      CommandEvent* /*timingEvent */,
346                      const char*                               kernelName,
347                      const std::array<void*, sizeof...(Args)>& kernelArgs)
348 {
349     dim3 blockSize(config.blockSize[0], config.blockSize[1], config.blockSize[2]);
350     dim3 gridSize(config.gridSize[0], config.gridSize[1], config.gridSize[2]);
351     cudaLaunchKernel((void*)kernel, gridSize, blockSize, const_cast<void**>(kernelArgs.data()),
352                      config.sharedMemorySize, config.stream);
353
354     cudaError_t status = cudaGetLastError();
355     if (cudaSuccess != status)
356     {
357         const std::string errorMessage =
358                 "GPU kernel (" + std::string(kernelName)
359                 + ") failed to launch: " + std::string(cudaGetErrorString(status));
360         GMX_THROW(gmx::InternalError(errorMessage));
361     }
362 }
363
364 #endif