384242f82743230c7aa2bb17ca29d30f9b349ecd
[alexxy/gromacs.git] / src / gromacs / gpu_utils / oclutils.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017, 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 /*! \internal \file
36  *  \brief Define utility routines for OpenCL
37  *
38  *  \author Anca Hamuraru <anca@streamcomputing.eu>
39  */
40 #include "gmxpre.h"
41
42 #include "oclutils.h"
43
44 #include <stdlib.h>
45
46 #include <cassert>
47 #include <cstdio>
48
49 #include <string>
50
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
53
54 /*! \brief Launches synchronous or asynchronous host to device memory copy.
55  *
56  *  If copy_event is not NULL, on return it will contain an event object
57  *  identifying this particular host to device operation. The event can further
58  *  be used to queue a wait for this operation or to query profiling information.
59  */
60 static int ocl_copy_H2D_generic(cl_mem d_dest, void* h_src,
61                                 size_t offset, size_t bytes,
62                                 bool bAsync /* = false*/,
63                                 cl_command_queue command_queue,
64                                 cl_event *copy_event)
65 {
66     cl_int gmx_unused cl_error;
67
68     if (d_dest == NULL || h_src == NULL || bytes == 0)
69     {
70         return -1;
71     }
72
73     if (bAsync)
74     {
75         cl_error = clEnqueueWriteBuffer(command_queue, d_dest, CL_FALSE, offset, bytes, h_src, 0, NULL, copy_event);
76         assert(cl_error == CL_SUCCESS);
77         // TODO: handle errors
78     }
79     else
80     {
81         cl_error = clEnqueueWriteBuffer(command_queue, d_dest, CL_TRUE, offset, bytes, h_src, 0, NULL, copy_event);
82         assert(cl_error == CL_SUCCESS);
83         // TODO: handle errors
84     }
85
86     return 0;
87 }
88
89 /*! \brief Launches asynchronous host to device memory copy.
90  *
91  *  If copy_event is not NULL, on return it will contain an event object
92  *  identifying this particular host to device operation. The event can further
93  *  be used to queue a wait for this operation or to query profiling information.
94  */
95 int ocl_copy_H2D_async(cl_mem d_dest, void * h_src,
96                        size_t offset, size_t bytes,
97                        cl_command_queue command_queue,
98                        cl_event *copy_event)
99 {
100     return ocl_copy_H2D_generic(d_dest, h_src, offset, bytes, true, command_queue, copy_event);
101 }
102
103 /*! \brief Launches synchronous host to device memory copy.
104  */
105 int ocl_copy_H2D_sync(cl_mem d_dest, void * h_src,
106                       size_t offset, size_t bytes,
107                       cl_command_queue command_queue)
108 {
109     return ocl_copy_H2D_generic(d_dest, h_src, offset, bytes, false, command_queue, NULL);
110 }
111
112 /*! \brief Launches synchronous or asynchronous device to host memory copy.
113  *
114  *  If copy_event is not NULL, on return it will contain an event object
115  *  identifying this particular device to host operation. The event can further
116  *  be used to queue a wait for this operation or to query profiling information.
117  */
118 static int ocl_copy_D2H_generic(void * h_dest, cl_mem d_src,
119                                 size_t offset, size_t bytes,
120                                 bool bAsync,
121                                 cl_command_queue command_queue,
122                                 cl_event *copy_event)
123 {
124     cl_int gmx_unused cl_error;
125
126     if (h_dest == NULL || d_src == NULL || bytes == 0)
127     {
128         return -1;
129     }
130
131     if (bAsync)
132     {
133         cl_error = clEnqueueReadBuffer(command_queue, d_src, CL_FALSE, offset, bytes, h_dest, 0, NULL, copy_event);
134         assert(cl_error == CL_SUCCESS);
135         // TODO: handle errors
136     }
137     else
138     {
139         cl_error = clEnqueueReadBuffer(command_queue, d_src, CL_TRUE, offset, bytes, h_dest, 0, NULL, copy_event);
140         assert(cl_error == CL_SUCCESS);
141         // TODO: handle errors
142     }
143
144     return 0;
145 }
146
147 /*! \brief Launches asynchronous device to host memory copy.
148  *
149  *  If copy_event is not NULL, on return it will contain an event object
150  *  identifying this particular host to device operation. The event can further
151  *  be used to queue a wait for this operation or to query profiling information.
152  */
153 int ocl_copy_D2H_async(void * h_dest, cl_mem d_src,
154                        size_t offset, size_t bytes,
155                        cl_command_queue command_queue,
156                        cl_event *copy_event)
157 {
158     return ocl_copy_D2H_generic(h_dest, d_src, offset, bytes, true, command_queue, copy_event);
159 }
160
161 /*! \brief \brief Allocates nbytes of host memory. Use ocl_free to free memory allocated with this function.
162  *
163  *  \todo
164  *  This function should allocate page-locked memory to help reduce D2H and H2D
165  *  transfer times, similar with pmalloc from pmalloc_cuda.cu.
166  *
167  * \param[in,out]    h_ptr   Pointer where to store the address of the newly allocated buffer.
168  * \param[in]        nbytes  Size in bytes of the buffer to be allocated.
169  */
170 void ocl_pmalloc(void **h_ptr, size_t nbytes)
171 {
172     /* Need a temporary type whose size is 1 byte, so that the
173      * implementation of snew_aligned can cope without issuing
174      * warnings. */
175     char **temporary = reinterpret_cast<char **>(h_ptr);
176
177     /* 16-byte alignment is required by the neighbour-searching code,
178      * because it uses four-wide SIMD for bounding-box calculation.
179      * However, when we organize using page-locked memory for
180      * device-host transfers, it will probably need to be aligned to a
181      * 4kb page, like CUDA does. */
182     snew_aligned(*temporary, nbytes, 16);
183 }
184
185 /*! \brief Frees memory allocated with ocl_pmalloc.
186  *
187  * \param[in]    h_ptr   Buffer allocated with ocl_pmalloc that needs to be freed.
188  */
189 void ocl_pfree(void *h_ptr)
190 {
191
192     if (h_ptr)
193     {
194         sfree_aligned(h_ptr);
195     }
196     return;
197 }
198
199 /*! \brief Convert error code to diagnostic string */
200 std::string ocl_get_error_string(cl_int error)
201 {
202     switch (error)
203     {
204         // run-time and JIT compiler errors
205         case 0: return "CL_SUCCESS";
206         case -1: return "CL_DEVICE_NOT_FOUND";
207         case -2: return "CL_DEVICE_NOT_AVAILABLE";
208         case -3: return "CL_COMPILER_NOT_AVAILABLE";
209         case -4: return "CL_MEM_OBJECT_ALLOCATION_FAILURE";
210         case -5: return "CL_OUT_OF_RESOURCES";
211         case -6: return "CL_OUT_OF_HOST_MEMORY";
212         case -7: return "CL_PROFILING_INFO_NOT_AVAILABLE";
213         case -8: return "CL_MEM_COPY_OVERLAP";
214         case -9: return "CL_IMAGE_FORMAT_MISMATCH";
215         case -10: return "CL_IMAGE_FORMAT_NOT_SUPPORTED";
216         case -11: return "CL_BUILD_PROGRAM_FAILURE";
217         case -12: return "CL_MAP_FAILURE";
218         case -13: return "CL_MISALIGNED_SUB_BUFFER_OFFSET";
219         case -14: return "CL_EXEC_STATUS_ERROR_FOR_EVENTS_IN_WAIT_LIST";
220         case -15: return "CL_COMPILE_PROGRAM_FAILURE";
221         case -16: return "CL_LINKER_NOT_AVAILABLE";
222         case -17: return "CL_LINK_PROGRAM_FAILURE";
223         case -18: return "CL_DEVICE_PARTITION_FAILED";
224         case -19: return "CL_KERNEL_ARG_INFO_NOT_AVAILABLE";
225
226         // compile-time errors
227         case -30: return "CL_INVALID_VALUE";
228         case -31: return "CL_INVALID_DEVICE_TYPE";
229         case -32: return "CL_INVALID_PLATFORM";
230         case -33: return "CL_INVALID_DEVICE";
231         case -34: return "CL_INVALID_CONTEXT";
232         case -35: return "CL_INVALID_QUEUE_PROPERTIES";
233         case -36: return "CL_INVALID_COMMAND_QUEUE";
234         case -37: return "CL_INVALID_HOST_PTR";
235         case -38: return "CL_INVALID_MEM_OBJECT";
236         case -39: return "CL_INVALID_IMAGE_FORMAT_DESCRIPTOR";
237         case -40: return "CL_INVALID_IMAGE_SIZE";
238         case -41: return "CL_INVALID_SAMPLER";
239         case -42: return "CL_INVALID_BINARY";
240         case -43: return "CL_INVALID_BUILD_OPTIONS";
241         case -44: return "CL_INVALID_PROGRAM";
242         case -45: return "CL_INVALID_PROGRAM_EXECUTABLE";
243         case -46: return "CL_INVALID_KERNEL_NAME";
244         case -47: return "CL_INVALID_KERNEL_DEFINITION";
245         case -48: return "CL_INVALID_KERNEL";
246         case -49: return "CL_INVALID_ARG_INDEX";
247         case -50: return "CL_INVALID_ARG_VALUE";
248         case -51: return "CL_INVALID_ARG_SIZE";
249         case -52: return "CL_INVALID_KERNEL_ARGS";
250         case -53: return "CL_INVALID_WORK_DIMENSION";
251         case -54: return "CL_INVALID_WORK_GROUP_SIZE";
252         case -55: return "CL_INVALID_WORK_ITEM_SIZE";
253         case -56: return "CL_INVALID_GLOBAL_OFFSET";
254         case -57: return "CL_INVALID_EVENT_WAIT_LIST";
255         case -58: return "CL_INVALID_EVENT";
256         case -59: return "CL_INVALID_OPERATION";
257         case -60: return "CL_INVALID_GL_OBJECT";
258         case -61: return "CL_INVALID_BUFFER_SIZE";
259         case -62: return "CL_INVALID_MIP_LEVEL";
260         case -63: return "CL_INVALID_GLOBAL_WORK_SIZE";
261         case -64: return "CL_INVALID_PROPERTY";
262         case -65: return "CL_INVALID_IMAGE_DESCRIPTOR";
263         case -66: return "CL_INVALID_COMPILER_OPTIONS";
264         case -67: return "CL_INVALID_LINKER_OPTIONS";
265         case -68: return "CL_INVALID_DEVICE_PARTITION_COUNT";
266
267         // extension errors
268         case -1000: return "CL_INVALID_GL_SHAREGROUP_REFERENCE_KHR";
269         case -1001: return "CL_PLATFORM_NOT_FOUND_KHR";
270         case -1002: return "CL_INVALID_D3D10_DEVICE_KHR";
271         case -1003: return "CL_INVALID_D3D10_RESOURCE_KHR";
272         case -1004: return "CL_D3D10_RESOURCE_ALREADY_ACQUIRED_KHR";
273         case -1005: return "CL_D3D10_RESOURCE_NOT_ACQUIRED_KHR";
274         default:    return "Unknown OpenCL error: " +
275                    std::to_string(static_cast<int32_t>(error));
276     }
277 }