Fix compilation issues with ARM SIMD
[alexxy/gromacs.git] / src / gromacs / gpu_utils / cudautils.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2014,2015,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
36 #include "gmxpre.h"
37
38 #include "cudautils.cuh"
39
40 #include <stdlib.h>
41
42 #include "gromacs/utility/smalloc.h"
43
44 /*** Generic CUDA data operation wrappers ***/
45
46 /*! Launches synchronous or asynchronous host to device memory copy.
47  *
48  *  The copy is launched in stream s or if not specified, in stream 0.
49  */
50 static int cu_copy_D2H_generic(void * h_dest, void * d_src, size_t bytes,
51                                bool bAsync = false, cudaStream_t s = 0)
52 {
53     cudaError_t stat;
54
55     if (h_dest == NULL || d_src == NULL || bytes == 0)
56     {
57         return -1;
58     }
59
60     if (bAsync)
61     {
62         stat = cudaMemcpyAsync(h_dest, d_src, bytes, cudaMemcpyDeviceToHost, s);
63         CU_RET_ERR(stat, "DtoH cudaMemcpyAsync failed");
64
65     }
66     else
67     {
68         stat = cudaMemcpy(h_dest, d_src, bytes, cudaMemcpyDeviceToHost);
69         CU_RET_ERR(stat, "DtoH cudaMemcpy failed");
70     }
71
72     return 0;
73 }
74
75 int cu_copy_D2H(void * h_dest, void * d_src, size_t bytes)
76 {
77     return cu_copy_D2H_generic(h_dest, d_src, bytes, false);
78 }
79
80 /*!
81  *  The copy is launched in stream s or if not specified, in stream 0.
82  */
83 int cu_copy_D2H_async(void * h_dest, void * d_src, size_t bytes, cudaStream_t s = 0)
84 {
85     return cu_copy_D2H_generic(h_dest, d_src, bytes, true, s);
86 }
87
88 /*! Launches synchronous or asynchronous device to host memory copy.
89  *
90  *  The copy is launched in stream s or if not specified, in stream 0.
91  */
92 static int cu_copy_H2D_generic(void * d_dest, void * h_src, size_t bytes,
93                                bool bAsync = false, cudaStream_t s = 0)
94 {
95     cudaError_t stat;
96
97     if (d_dest == NULL || h_src == NULL || bytes == 0)
98     {
99         return -1;
100     }
101
102     if (bAsync)
103     {
104         stat = cudaMemcpyAsync(d_dest, h_src, bytes, cudaMemcpyHostToDevice, s);
105         CU_RET_ERR(stat, "HtoD cudaMemcpyAsync failed");
106     }
107     else
108     {
109         stat = cudaMemcpy(d_dest, h_src, bytes, cudaMemcpyHostToDevice);
110         CU_RET_ERR(stat, "HtoD cudaMemcpy failed");
111     }
112
113     return 0;
114 }
115
116 int cu_copy_H2D(void * d_dest, void * h_src, size_t bytes)
117 {
118     return cu_copy_H2D_generic(d_dest, h_src, bytes, false);
119 }
120
121 /*!
122  *  The copy is launched in stream s or if not specified, in stream 0.
123  */
124 int cu_copy_H2D_async(void * d_dest, void * h_src, size_t bytes, cudaStream_t s = 0)
125 {
126     return cu_copy_H2D_generic(d_dest, h_src, bytes, true, s);
127 }
128
129 int cu_copy_H2D_alloc(void ** d_dest, void * h_src, size_t bytes)
130 {
131     cudaError_t stat;
132
133     if (d_dest == NULL || h_src == NULL || bytes == 0)
134     {
135         return -1;
136     }
137
138     stat = cudaMalloc(d_dest, bytes);
139     CU_RET_ERR(stat, "cudaMalloc failed in cu_copy_H2D_alloc");
140
141     return cu_copy_H2D(*d_dest, h_src, bytes);
142 }
143
144 float cu_event_elapsed(cudaEvent_t start, cudaEvent_t end)
145 {
146     float       t = 0.0;
147     cudaError_t stat;
148
149     stat = cudaEventElapsedTime(&t, start, end);
150     CU_RET_ERR(stat, "cudaEventElapsedTime failed in cu_event_elapsed");
151
152     return t;
153 }
154
155 int cu_wait_event(cudaEvent_t e)
156 {
157     cudaError_t s;
158
159     s = cudaEventSynchronize(e);
160     CU_RET_ERR(s, "cudaEventSynchronize failed in cu_wait_event");
161
162     return 0;
163 }
164
165 /*!
166  *  If time != NULL it also calculates the time elapsed between start and end and
167  *  return this is milliseconds.
168  */
169 int cu_wait_event_time(cudaEvent_t end, cudaEvent_t start, float *time)
170 {
171     cudaError_t s;
172
173     s = cudaEventSynchronize(end);
174     CU_RET_ERR(s, "cudaEventSynchronize failed in cu_wait_event");
175
176     if (time)
177     {
178         *time = cu_event_elapsed(start, end);
179     }
180
181     return 0;
182 }
183
184 /**** Operation on buffered arrays (arrays with "over-allocation" in gmx wording) *****/
185
186 /*!
187  * If the pointers to the size variables are NULL no resetting happens.
188  */
189 void cu_free_buffered(void *d_ptr, int *n, int *nalloc)
190 {
191     cudaError_t stat;
192
193     if (d_ptr)
194     {
195         stat = cudaFree(d_ptr);
196         CU_RET_ERR(stat, "cudaFree failed");
197     }
198
199     if (n)
200     {
201         *n = -1;
202     }
203
204     if (nalloc)
205     {
206         *nalloc = -1;
207     }
208 }
209
210 /*!
211  *  Reallocation of the memory pointed by d_ptr and copying of the data from
212  *  the location pointed by h_src host-side pointer is done. Allocation is
213  *  buffered and therefore freeing is only needed if the previously allocated
214  *  space is not enough.
215  *  The H2D copy is launched in stream s and can be done synchronously or
216  *  asynchronously (the default is the latter).
217  */
218 void cu_realloc_buffered(void **d_dest, void *h_src,
219                          size_t type_size,
220                          int *curr_size, int *curr_alloc_size,
221                          int req_size,
222                          cudaStream_t s,
223                          bool bAsync = true)
224 {
225     cudaError_t stat;
226
227     if (d_dest == NULL || req_size < 0)
228     {
229         return;
230     }
231
232     /* reallocate only if the data does not fit = allocation size is smaller
233        than the current requested size */
234     if (req_size > *curr_alloc_size)
235     {
236         /* only free if the array has already been initialized */
237         if (*curr_alloc_size >= 0)
238         {
239             cu_free_buffered(*d_dest, curr_size, curr_alloc_size);
240         }
241
242         *curr_alloc_size = over_alloc_large(req_size);
243
244         stat = cudaMalloc(d_dest, *curr_alloc_size * type_size);
245         CU_RET_ERR(stat, "cudaMalloc failed in cu_free_buffered");
246     }
247
248     /* size could have changed without actual reallocation */
249     *curr_size = req_size;
250
251     /* upload to device */
252     if (h_src)
253     {
254         if (bAsync)
255         {
256             cu_copy_H2D_async(*d_dest, h_src, *curr_size * type_size, s);
257         }
258         else
259         {
260             cu_copy_H2D(*d_dest, h_src,  *curr_size * type_size);
261         }
262     }
263 }