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