Fix race condition in GPU X halo exchange with GPU update
[alexxy/gromacs.git] / src / gromacs / domdec / gpuhaloexchange_impl.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \internal \file
36  *
37  * \brief Implements GPU halo exchange using CUDA.
38  *
39  *
40  * \author Alan Gray <alang@nvidia.com>
41  *
42  * \ingroup module_domdec
43  */
44 #include "gmxpre.h"
45
46 #include "gpuhaloexchange_impl.cuh"
47
48 #include "config.h"
49
50 #include <assert.h>
51 #include <stdio.h>
52
53 #include "gromacs/domdec/domdec.h"
54 #include "gromacs/domdec/domdec_struct.h"
55 #include "gromacs/domdec/gpuhaloexchange.h"
56 #include "gromacs/gpu_utils/cudautils.cuh"
57 #include "gromacs/gpu_utils/devicebuffer.h"
58 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
59 #include "gromacs/gpu_utils/vectype_ops.cuh"
60 #include "gromacs/pbcutil/ishift.h"
61
62 #include "domdec_internal.h"
63
64 namespace gmx
65 {
66
67 //! Number of CUDA threads in a block
68 // TODO Optimize this through experimentation
69 constexpr static int c_threadsPerBlock = 256;
70
71 template<bool usePBC>
72 __global__ void packSendBufKernel(float3* __restrict__ dataPacked,
73                                   const float3* __restrict__ data,
74                                   const int* __restrict__ map,
75                                   const int    mapSize,
76                                   const float3 coordinateShift)
77 {
78     int           threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
79     float3*       gm_dataDest = &dataPacked[threadIndex];
80     const float3* gm_dataSrc  = &data[map[threadIndex]];
81
82     if (threadIndex < mapSize)
83     {
84         if (usePBC)
85         {
86             *gm_dataDest = *gm_dataSrc + coordinateShift;
87         }
88         else
89         {
90             *gm_dataDest = *gm_dataSrc;
91         }
92     }
93
94     return;
95 }
96
97 /*! \brief unpack non-local force data buffer on the GPU using pre-populated "map" containing index
98  * information \param[out] data        full array of force values \param[in]  dataPacked  packed
99  * array of force values to be transferred \param[in]  map         array of indices defining mapping
100  * from full to packed array \param[in]  mapSize     number of elements in map array
101  */
102 template<bool accumulate>
103 __global__ void unpackRecvBufKernel(float3* __restrict__ data,
104                                     const float3* __restrict__ dataPacked,
105                                     const int* __restrict__ map,
106                                     const int mapSize)
107 {
108
109     int           threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
110     const float3* gm_dataSrc  = &dataPacked[threadIndex];
111     float3*       gm_dataDest = &data[map[threadIndex]];
112
113     if (threadIndex < mapSize)
114     {
115         if (accumulate)
116         {
117             *gm_dataDest += *gm_dataSrc;
118         }
119         else
120         {
121             *gm_dataDest = *gm_dataSrc;
122         }
123     }
124
125     return;
126 }
127
128 void GpuHaloExchange::Impl::reinitHalo(float3* d_coordinatesBuffer, float3* d_forcesBuffer)
129 {
130
131     d_x_ = d_coordinatesBuffer;
132     d_f_ = d_forcesBuffer;
133
134     cudaStream_t                 stream  = nonLocalStream_;
135     int                          nzone   = 1;
136     const gmx_domdec_comm_t&     comm    = *dd_->comm;
137     const gmx_domdec_comm_dim_t& cd      = comm.cd[0];
138     const gmx_domdec_ind_t&      ind     = cd.ind[0];
139     int                          newSize = ind.nsend[nzone + 1];
140
141     GMX_RELEASE_ASSERT(cd.numPulses() == 1,
142                        "Multiple pulses are not yet supported in GPU halo exchange");
143     GMX_ASSERT(cd.receiveInPlace, "Out-of-place receive is not yet supported in GPU halo exchange");
144
145     // reallocates only if needed
146     h_indexMap_.resize(newSize);
147     // reallocate on device only if needed
148     if (newSize > maxPackedBufferSize_)
149     {
150         reallocateDeviceBuffer(&d_indexMap_, newSize, &indexMapSize_, &indexMapSizeAlloc_, nullptr);
151         reallocateDeviceBuffer(&d_sendBuf_, newSize, &sendBufSize_, &sendBufSizeAlloc_, nullptr);
152         reallocateDeviceBuffer(&d_recvBuf_, newSize, &recvBufSize_, &recvBufSizeAlloc_, nullptr);
153         maxPackedBufferSize_ = newSize;
154     }
155
156     xSendSize_ = newSize;
157 #if GMX_MPI
158     MPI_Sendrecv(&xSendSize_, sizeof(int), MPI_BYTE, sendRankX_, 0, &xRecvSize_, sizeof(int),
159                  MPI_BYTE, recvRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
160 #endif
161     fSendSize_ = xRecvSize_;
162     fRecvSize_ = xSendSize_;
163
164     numHomeAtoms_ = comm.atomRanges.numHomeAtoms(); // offset for data recieved by this rank
165
166     GMX_ASSERT(ind.index.size() == h_indexMap_.size(), "Size mismatch");
167     std::copy(ind.index.begin(), ind.index.end(), h_indexMap_.begin());
168
169     copyToDeviceBuffer(&d_indexMap_, h_indexMap_.data(), 0, newSize, stream,
170                        GpuApiCallBehavior::Async, nullptr);
171
172     // This rank will push data to its neighbor, so needs to know
173     // the remote receive address and similarly send its receive
174     // address to other neighbour. We can do this here in reinit fn
175     // since the pointers will not change until the next NS step.
176
177     // Coordinates buffer:
178 #if GMX_MPI
179     void* recvPtr = static_cast<void*>(&d_coordinatesBuffer[numHomeAtoms_]);
180     MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankX_, 0, &remoteXPtr_, sizeof(void*),
181                  MPI_BYTE, sendRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
182
183     // Force buffer:
184     recvPtr = static_cast<void*>(d_recvBuf_);
185     MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankF_, 0, &remoteFPtr_, sizeof(void*),
186                  MPI_BYTE, sendRankF_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
187 #endif
188
189
190     return;
191 }
192
193 void GpuHaloExchange::Impl::communicateHaloCoordinates(const matrix          box,
194                                                        GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
195 {
196
197     // ensure stream waits until coordinate data is available on device
198     coordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
199
200     // launch kernel to pack send buffer
201     KernelLaunchConfig config;
202     config.blockSize[0]     = c_threadsPerBlock;
203     config.blockSize[1]     = 1;
204     config.blockSize[2]     = 1;
205     config.gridSize[0]      = (xSendSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
206     config.gridSize[1]      = 1;
207     config.gridSize[2]      = 1;
208     config.sharedMemorySize = 0;
209     config.stream           = nonLocalStream_;
210
211     const float3* sendBuf  = d_sendBuf_;
212     const float3* d_x      = d_x_;
213     const int*    indexMap = d_indexMap_;
214     const int     size     = xSendSize_;
215     // The coordinateShift changes between steps when we have
216     // performed a DD partition, or have updated the box e.g. when
217     // performing pressure coupling. So, for simplicity, the box
218     // is used every step to pass the shift vector as an argument of
219     // the packing kernel.
220     //
221     // Because only one-dimensional DD is supported, the coordinate
222     // shift only needs to handle that dimension.
223     const int    dimensionIndex = dd_->dim[0];
224     const float3 coordinateShift{ box[dimensionIndex][XX], box[dimensionIndex][YY],
225                                   box[dimensionIndex][ZZ] };
226
227     // Avoid launching kernel when there is no work to do
228     if (size > 0)
229     {
230         auto kernelFn = usePBC_ ? packSendBufKernel<true> : packSendBufKernel<false>;
231
232         const auto kernelArgs = prepareGpuKernelArguments(kernelFn, config, &sendBuf, &d_x,
233                                                           &indexMap, &size, &coordinateShift);
234
235         launchGpuKernel(kernelFn, config, nullptr, "Domdec GPU Apply X Halo Exchange", kernelArgs);
236     }
237
238     communicateHaloData(d_x_, HaloQuantity::HaloCoordinates, coordinatesReadyOnDeviceEvent);
239
240     return;
241 }
242
243 // The following method should be called after non-local buffer operations,
244 // and before the local buffer operations. It operates in the non-local stream.
245 void GpuHaloExchange::Impl::communicateHaloForces(bool accumulateForces)
246 {
247
248     // Communicate halo data (in non-local stream)
249     communicateHaloData(d_f_, HaloQuantity::HaloForces, nullptr);
250
251     float3* d_f = d_f_;
252
253     if (!accumulateForces)
254     {
255         // Clear local portion of force array (in local stream)
256         cudaMemsetAsync(d_f, 0, numHomeAtoms_ * sizeof(rvec), localStream_);
257     }
258
259     // ensure non-local stream waits for local stream, due to dependence on
260     // the previous H2D copy of CPU forces (if accumulateForces is true)
261     // or the above clearing.
262     // TODO remove this dependency on localStream - edmine issue #3093
263     GpuEventSynchronizer eventLocal;
264     eventLocal.markEvent(localStream_);
265     eventLocal.enqueueWaitEvent(nonLocalStream_);
266
267     // Unpack halo buffer into force array
268
269     KernelLaunchConfig config;
270     config.blockSize[0]     = c_threadsPerBlock;
271     config.blockSize[1]     = 1;
272     config.blockSize[2]     = 1;
273     config.gridSize[0]      = (fRecvSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
274     config.gridSize[1]      = 1;
275     config.gridSize[2]      = 1;
276     config.sharedMemorySize = 0;
277     config.stream           = nonLocalStream_;
278
279     const float3* recvBuf  = d_recvBuf_;
280     const int*    indexMap = d_indexMap_;
281     const int     size     = fRecvSize_;
282
283     if (size > 0)
284     {
285         auto kernelFn = accumulateForces ? unpackRecvBufKernel<true> : unpackRecvBufKernel<false>;
286
287         const auto kernelArgs =
288                 prepareGpuKernelArguments(kernelFn, config, &d_f, &recvBuf, &indexMap, &size);
289
290         launchGpuKernel(kernelFn, config, nullptr, "Domdec GPU Apply F Halo Exchange", kernelArgs);
291     }
292 }
293
294
295 void GpuHaloExchange::Impl::communicateHaloData(float3*               d_ptr,
296                                                 HaloQuantity          haloQuantity,
297                                                 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
298 {
299
300     void* sendPtr;
301     int   sendSize;
302     void* remotePtr;
303     int   sendRank;
304     int   recvRank;
305
306     if (haloQuantity == HaloQuantity::HaloCoordinates)
307     {
308         sendPtr   = static_cast<void*>(d_sendBuf_);
309         sendSize  = xSendSize_;
310         remotePtr = remoteXPtr_;
311         sendRank  = sendRankX_;
312         recvRank  = recvRankX_;
313
314 #if GMX_MPI
315         // Wait for event from receiving task that remote coordinates are ready, and enqueue that event to stream used
316         // for subsequent data push. This avoids a race condition with the remote data being written in the previous timestep.
317         // Similarly send event to task that will push data to this task.
318         GpuEventSynchronizer* remoteCoordinatesReadyOnDeviceEvent;
319         MPI_Sendrecv(&coordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*), MPI_BYTE,
320                      recvRank, 0, &remoteCoordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*),
321                      MPI_BYTE, sendRank, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
322         remoteCoordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
323 #endif
324     }
325     else
326     {
327         sendPtr   = static_cast<void*>(&(d_ptr[numHomeAtoms_]));
328         sendSize  = fSendSize_;
329         remotePtr = remoteFPtr_;
330         sendRank  = sendRankF_;
331         recvRank  = recvRankF_;
332     }
333
334     communicateHaloDataWithCudaDirect(sendPtr, sendSize, sendRank, remotePtr, recvRank);
335 }
336
337 void GpuHaloExchange::Impl::communicateHaloDataWithCudaDirect(void* sendPtr,
338                                                               int   sendSize,
339                                                               int   sendRank,
340                                                               void* remotePtr,
341                                                               int   recvRank)
342 {
343
344     cudaError_t  stat;
345     cudaStream_t stream = nonLocalStream_;
346
347     // We asynchronously push data to remote rank. The remote
348     // destination pointer has already been set in the init fn.  We
349     // don't need to worry about overwriting data the remote ranks
350     // still needs since the halo exchange is just done once per
351     // timestep, for each of X and F.
352
353     // send data to neighbor, if any data exists to send
354     if (sendSize > 0)
355     {
356         stat = cudaMemcpyAsync(remotePtr, sendPtr, sendSize * DIM * sizeof(float),
357                                cudaMemcpyDeviceToDevice, stream);
358         CU_RET_ERR(stat, "cudaMemcpyAsync on GPU Domdec CUDA direct data transfer failed");
359     }
360
361 #if GMX_MPI
362     // ensure pushed data has arrived before remote rank progresses
363     // This rank records an event and sends it to the remote rank which has just been pushed data.
364     // This rank recieves event from remote rank which has pushed data here, and enqueues that event
365     // to its stream.
366     GpuEventSynchronizer* haloDataTransferRemote;
367
368     haloDataTransferLaunched_->markEvent(stream);
369
370     MPI_Sendrecv(&haloDataTransferLaunched_, sizeof(GpuEventSynchronizer*), MPI_BYTE, sendRank, 0,
371                  &haloDataTransferRemote, sizeof(GpuEventSynchronizer*), MPI_BYTE, recvRank, 0,
372                  mpi_comm_mysim_, MPI_STATUS_IGNORE);
373
374     haloDataTransferRemote->enqueueWaitEvent(stream);
375 #else
376     GMX_UNUSED_VALUE(sendRank);
377     GMX_UNUSED_VALUE(recvRank);
378 #endif
379 }
380
381 /*! \brief Create Domdec GPU object */
382 GpuHaloExchange::Impl::Impl(gmx_domdec_t* dd, MPI_Comm mpi_comm_mysim, void* localStream, void* nonLocalStream) :
383     dd_(dd),
384     sendRankX_(dd->neighbor[0][1]),
385     recvRankX_(dd->neighbor[0][0]),
386     sendRankF_(dd->neighbor[0][0]),
387     recvRankF_(dd->neighbor[0][1]),
388     usePBC_(dd->ci[dd->dim[0]] == 0),
389     haloDataTransferLaunched_(new GpuEventSynchronizer()),
390     mpi_comm_mysim_(mpi_comm_mysim),
391     localStream_(*static_cast<cudaStream_t*>(localStream)),
392     nonLocalStream_(*static_cast<cudaStream_t*>(nonLocalStream))
393 {
394
395     GMX_RELEASE_ASSERT(GMX_THREAD_MPI,
396                        "GPU Halo exchange is currently only supported with thread-MPI enabled");
397
398     if (dd->ndim > 1)
399     {
400         gmx_fatal(FARGS, "Error: dd->ndim > 1 is not yet supported in GPU halo exchange");
401     }
402
403     if (usePBC_ && dd->unitCellInfo.haveScrewPBC)
404     {
405         gmx_fatal(FARGS, "Error: screw is not yet supported in GPU halo exchange\n");
406     }
407
408     changePinningPolicy(&h_indexMap_, gmx::PinningPolicy::PinnedIfSupported);
409
410     allocateDeviceBuffer(&d_fShift_, 1, nullptr);
411 }
412
413 GpuHaloExchange::Impl::~Impl()
414 {
415     freeDeviceBuffer(&d_indexMap_);
416     freeDeviceBuffer(&d_sendBuf_);
417     freeDeviceBuffer(&d_recvBuf_);
418     freeDeviceBuffer(&d_fShift_);
419     delete haloDataTransferLaunched_;
420 }
421
422 GpuHaloExchange::GpuHaloExchange(gmx_domdec_t* dd, MPI_Comm mpi_comm_mysim, void* localStream, void* nonLocalStream) :
423     impl_(new Impl(dd, mpi_comm_mysim, localStream, nonLocalStream))
424 {
425 }
426
427 GpuHaloExchange::~GpuHaloExchange() = default;
428
429 void GpuHaloExchange::reinitHalo(DeviceBuffer<float> d_coordinatesBuffer, DeviceBuffer<float> d_forcesBuffer)
430 {
431     impl_->reinitHalo(reinterpret_cast<float3*>(d_coordinatesBuffer),
432                       reinterpret_cast<float3*>(d_forcesBuffer));
433 }
434
435 void GpuHaloExchange::communicateHaloCoordinates(const matrix          box,
436                                                  GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
437 {
438     impl_->communicateHaloCoordinates(box, coordinatesReadyOnDeviceEvent);
439 }
440
441 void GpuHaloExchange::communicateHaloForces(bool accumulateForces)
442 {
443     impl_->communicateHaloForces(accumulateForces);
444 }
445
446 } // namespace gmx