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