cbcff3defbc7274da10cbbbd2119fdb4cf129476
[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/timing/wallcycle.h"
65 #include "gromacs/utility/gmxmpi.h"
66
67 #include "domdec_internal.h"
68
69 namespace gmx
70 {
71
72 //! Number of CUDA threads in a block
73 // TODO Optimize this through experimentation
74 constexpr static int c_threadsPerBlock = 256;
75
76 template<bool usePBC>
77 __global__ void packSendBufKernel(float3* __restrict__ dataPacked,
78                                   const float3* __restrict__ data,
79                                   const int* __restrict__ map,
80                                   const int    mapSize,
81                                   const float3 coordinateShift)
82 {
83     int           threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
84     float3*       gm_dataDest = &dataPacked[threadIndex];
85     const float3* gm_dataSrc  = &data[map[threadIndex]];
86
87     if (threadIndex < mapSize)
88     {
89         if (usePBC)
90         {
91             *gm_dataDest = *gm_dataSrc + coordinateShift;
92         }
93         else
94         {
95             *gm_dataDest = *gm_dataSrc;
96         }
97     }
98
99     return;
100 }
101
102 /*! \brief unpack non-local force data buffer on the GPU using pre-populated "map" containing index
103  * information \param[out] data        full array of force values \param[in]  dataPacked  packed
104  * array of force values to be transferred \param[in]  map         array of indices defining mapping
105  * from full to packed array \param[in]  mapSize     number of elements in map array
106  */
107 template<bool accumulate>
108 __global__ void unpackRecvBufKernel(float3* __restrict__ data,
109                                     const float3* __restrict__ dataPacked,
110                                     const int* __restrict__ map,
111                                     const int mapSize)
112 {
113
114     int           threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
115     const float3* gm_dataSrc  = &dataPacked[threadIndex];
116     float3*       gm_dataDest = &data[map[threadIndex]];
117
118     if (threadIndex < mapSize)
119     {
120         if (accumulate)
121         {
122             *gm_dataDest += *gm_dataSrc;
123         }
124         else
125         {
126             *gm_dataDest = *gm_dataSrc;
127         }
128     }
129
130     return;
131 }
132
133 void GpuHaloExchange::Impl::reinitHalo(float3* d_coordinatesBuffer, float3* d_forcesBuffer)
134 {
135
136     d_x_ = d_coordinatesBuffer;
137     d_f_ = d_forcesBuffer;
138
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, nonLocalStream_,
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     wallcycle_start(wcycle_, ewcLAUNCH_GPU);
210     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_MOVEX);
211
212     // launch kernel to pack send buffer
213     KernelLaunchConfig config;
214     config.blockSize[0]     = c_threadsPerBlock;
215     config.blockSize[1]     = 1;
216     config.blockSize[2]     = 1;
217     config.gridSize[0]      = (xSendSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
218     config.gridSize[1]      = 1;
219     config.gridSize[2]      = 1;
220     config.sharedMemorySize = 0;
221
222     const float3* sendBuf  = d_sendBuf_;
223     const float3* d_x      = d_x_;
224     const int*    indexMap = d_indexMap_;
225     const int     size     = xSendSize_;
226     // The coordinateShift changes between steps when we have
227     // performed a DD partition, or have updated the box e.g. when
228     // performing pressure coupling. So, for simplicity, the box
229     // is used every step to pass the shift vector as an argument of
230     // the packing kernel.
231     //
232     // Because only one-dimensional DD is supported, the coordinate
233     // shift only needs to handle that dimension.
234     const int    dimensionIndex = dd_->dim[0];
235     const float3 coordinateShift{ box[dimensionIndex][XX], box[dimensionIndex][YY],
236                                   box[dimensionIndex][ZZ] };
237
238     // Avoid launching kernel when there is no work to do
239     if (size > 0)
240     {
241         auto kernelFn = usePBC_ ? packSendBufKernel<true> : packSendBufKernel<false>;
242
243         const auto kernelArgs = prepareGpuKernelArguments(kernelFn, config, &sendBuf, &d_x,
244                                                           &indexMap, &size, &coordinateShift);
245
246         launchGpuKernel(kernelFn, config, nonLocalStream_, nullptr,
247                         "Domdec GPU Apply X Halo Exchange", kernelArgs);
248     }
249
250     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_MOVEX);
251     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
252
253     // Consider time spent in communicateHaloData as Comm.X counter
254     // ToDo: We need further refinement here as communicateHaloData includes launch time for cudamemcpyasync
255     wallcycle_start(wcycle_, ewcMOVEX);
256
257     communicateHaloData(d_x_, HaloQuantity::HaloCoordinates, coordinatesReadyOnDeviceEvent);
258
259     wallcycle_stop(wcycle_, ewcMOVEX);
260
261     return;
262 }
263
264 // The following method should be called after non-local buffer operations,
265 // and before the local buffer operations. It operates in the non-local stream.
266 void GpuHaloExchange::Impl::communicateHaloForces(bool accumulateForces)
267 {
268     // Consider time spent in communicateHaloData as Comm.F counter
269     // ToDo: We need further refinement here as communicateHaloData includes launch time for cudamemcpyasync
270     wallcycle_start(wcycle_, ewcMOVEF);
271
272     // Communicate halo data (in non-local stream)
273     communicateHaloData(d_f_, HaloQuantity::HaloForces, nullptr);
274
275     wallcycle_stop(wcycle_, ewcMOVEF);
276
277     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
278     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_MOVEF);
279
280     float3* d_f = d_f_;
281
282     if (pulse_ == (dd_->comm->cd[0].numPulses() - 1))
283     {
284         if (!accumulateForces)
285         {
286             // Clear local portion of force array (in local stream)
287             cudaMemsetAsync(d_f, 0, numHomeAtoms_ * sizeof(rvec), localStream_.stream());
288         }
289
290         // ensure non-local stream waits for local stream, due to dependence on
291         // the previous H2D copy of CPU forces (if accumulateForces is true)
292         // or the above clearing.
293         // TODO remove this dependency on localStream - edmine Issue #3093
294         GpuEventSynchronizer eventLocal;
295         eventLocal.markEvent(localStream_);
296         eventLocal.enqueueWaitEvent(nonLocalStream_);
297     }
298
299     // Unpack halo buffer into force array
300
301     KernelLaunchConfig config;
302     config.blockSize[0]     = c_threadsPerBlock;
303     config.blockSize[1]     = 1;
304     config.blockSize[2]     = 1;
305     config.gridSize[0]      = (fRecvSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
306     config.gridSize[1]      = 1;
307     config.gridSize[2]      = 1;
308     config.sharedMemorySize = 0;
309
310     const float3* recvBuf  = d_recvBuf_;
311     const int*    indexMap = d_indexMap_;
312     const int     size     = fRecvSize_;
313
314     if (pulse_ > 0)
315     {
316         // We need to accumulate rather than set, since it is possible
317         // that, in this pulse, a value could be written to a location
318         // corresponding to the halo region of a following pulse.
319         accumulateForces = true;
320     }
321
322     if (size > 0)
323     {
324         auto kernelFn = accumulateForces ? unpackRecvBufKernel<true> : unpackRecvBufKernel<false>;
325
326         const auto kernelArgs =
327                 prepareGpuKernelArguments(kernelFn, config, &d_f, &recvBuf, &indexMap, &size);
328
329         launchGpuKernel(kernelFn, config, nonLocalStream_, nullptr,
330                         "Domdec GPU Apply F Halo Exchange", kernelArgs);
331     }
332
333     if (pulse_ == 0)
334     {
335         fReadyOnDevice_.markEvent(nonLocalStream_);
336     }
337
338     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_MOVEF);
339     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
340 }
341
342
343 void GpuHaloExchange::Impl::communicateHaloData(float3*               d_ptr,
344                                                 HaloQuantity          haloQuantity,
345                                                 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
346 {
347
348     void* sendPtr;
349     int   sendSize;
350     void* remotePtr;
351     int   sendRank;
352     int   recvRank;
353
354     if (haloQuantity == HaloQuantity::HaloCoordinates)
355     {
356         sendPtr   = static_cast<void*>(d_sendBuf_);
357         sendSize  = xSendSize_;
358         remotePtr = remoteXPtr_;
359         sendRank  = sendRankX_;
360         recvRank  = recvRankX_;
361
362 #if GMX_MPI
363         // Wait for event from receiving task that remote coordinates are ready, and enqueue that event to stream used
364         // for subsequent data push. This avoids a race condition with the remote data being written in the previous timestep.
365         // Similarly send event to task that will push data to this task.
366         GpuEventSynchronizer* remoteCoordinatesReadyOnDeviceEvent;
367         MPI_Sendrecv(&coordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*), MPI_BYTE,
368                      recvRank, 0, &remoteCoordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*),
369                      MPI_BYTE, sendRank, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
370         remoteCoordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
371 #else
372         GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
373 #endif
374     }
375     else
376     {
377         int recvOffset = dd_->comm->atomRanges.end(DDAtomRanges::Type::Zones);
378         for (int p = pulse_; p < dd_->comm->cd[0].numPulses(); p++)
379         {
380             recvOffset -= dd_->comm->cd[0].ind[p].nrecv[nzone_ + 1];
381         }
382         sendPtr   = static_cast<void*>(&(d_ptr[recvOffset]));
383         sendSize  = fSendSize_;
384         remotePtr = remoteFPtr_;
385         sendRank  = sendRankF_;
386         recvRank  = recvRankF_;
387     }
388
389     communicateHaloDataWithCudaDirect(sendPtr, sendSize, sendRank, remotePtr, recvRank);
390 }
391
392 void GpuHaloExchange::Impl::communicateHaloDataWithCudaDirect(void* sendPtr,
393                                                               int   sendSize,
394                                                               int   sendRank,
395                                                               void* remotePtr,
396                                                               int   recvRank)
397 {
398
399     cudaError_t stat;
400
401     // We asynchronously push data to remote rank. The remote
402     // destination pointer has already been set in the init fn.  We
403     // don't need to worry about overwriting data the remote ranks
404     // still needs since the halo exchange is just done once per
405     // timestep, for each of X and F.
406
407     // send data to neighbor, if any data exists to send
408     if (sendSize > 0)
409     {
410         stat = cudaMemcpyAsync(remotePtr, sendPtr, sendSize * DIM * sizeof(float),
411                                cudaMemcpyDeviceToDevice, nonLocalStream_.stream());
412
413         CU_RET_ERR(stat, "cudaMemcpyAsync on GPU Domdec CUDA direct data transfer failed");
414     }
415
416 #if GMX_MPI
417     // ensure pushed data has arrived before remote rank progresses
418     // This rank records an event and sends it to the remote rank which has just been pushed data.
419     // This rank recieves event from remote rank which has pushed data here, and enqueues that event
420     // to its stream.
421     GpuEventSynchronizer* haloDataTransferRemote;
422
423     haloDataTransferLaunched_->markEvent(nonLocalStream_);
424
425     MPI_Sendrecv(&haloDataTransferLaunched_, sizeof(GpuEventSynchronizer*), MPI_BYTE, sendRank, 0,
426                  &haloDataTransferRemote, sizeof(GpuEventSynchronizer*), MPI_BYTE, recvRank, 0,
427                  mpi_comm_mysim_, MPI_STATUS_IGNORE);
428
429     haloDataTransferRemote->enqueueWaitEvent(nonLocalStream_);
430 #else
431     GMX_UNUSED_VALUE(sendRank);
432     GMX_UNUSED_VALUE(recvRank);
433 #endif
434 }
435
436 GpuEventSynchronizer* GpuHaloExchange::Impl::getForcesReadyOnDeviceEvent()
437 {
438     return &fReadyOnDevice_;
439 }
440
441 /*! \brief Create Domdec GPU object */
442 GpuHaloExchange::Impl::Impl(gmx_domdec_t*        dd,
443                             MPI_Comm             mpi_comm_mysim,
444                             const DeviceContext& deviceContext,
445                             const DeviceStream&  localStream,
446                             const DeviceStream&  nonLocalStream,
447                             int                  pulse,
448                             gmx_wallcycle*       wcycle) :
449     dd_(dd),
450     sendRankX_(dd->neighbor[0][1]),
451     recvRankX_(dd->neighbor[0][0]),
452     sendRankF_(dd->neighbor[0][0]),
453     recvRankF_(dd->neighbor[0][1]),
454     usePBC_(dd->ci[dd->dim[0]] == 0),
455     haloDataTransferLaunched_(new GpuEventSynchronizer()),
456     mpi_comm_mysim_(mpi_comm_mysim),
457     deviceContext_(deviceContext),
458     localStream_(localStream),
459     nonLocalStream_(nonLocalStream),
460     pulse_(pulse),
461     wcycle_(wcycle)
462 {
463
464     GMX_RELEASE_ASSERT(GMX_THREAD_MPI,
465                        "GPU Halo exchange is currently only supported with thread-MPI enabled");
466
467     if (dd->ndim > 1)
468     {
469         gmx_fatal(FARGS, "Error: dd->ndim > 1 is not yet supported in GPU halo exchange");
470     }
471
472     if (usePBC_ && dd->unitCellInfo.haveScrewPBC)
473     {
474         gmx_fatal(FARGS, "Error: screw is not yet supported in GPU halo exchange\n");
475     }
476
477     changePinningPolicy(&h_indexMap_, gmx::PinningPolicy::PinnedIfSupported);
478
479     allocateDeviceBuffer(&d_fShift_, 1, deviceContext_);
480 }
481
482 GpuHaloExchange::Impl::~Impl()
483 {
484     freeDeviceBuffer(&d_indexMap_);
485     freeDeviceBuffer(&d_sendBuf_);
486     freeDeviceBuffer(&d_recvBuf_);
487     freeDeviceBuffer(&d_fShift_);
488     delete haloDataTransferLaunched_;
489 }
490
491 GpuHaloExchange::GpuHaloExchange(gmx_domdec_t*        dd,
492                                  MPI_Comm             mpi_comm_mysim,
493                                  const DeviceContext& deviceContext,
494                                  const DeviceStream&  localStream,
495                                  const DeviceStream&  nonLocalStream,
496                                  int                  pulse,
497                                  gmx_wallcycle*       wcycle) :
498     impl_(new Impl(dd, mpi_comm_mysim, deviceContext, localStream, nonLocalStream, pulse, wcycle))
499 {
500 }
501
502 GpuHaloExchange::~GpuHaloExchange() = default;
503
504 void GpuHaloExchange::reinitHalo(DeviceBuffer<RVec> d_coordinatesBuffer, DeviceBuffer<RVec> d_forcesBuffer)
505 {
506     impl_->reinitHalo(asFloat3(d_coordinatesBuffer), asFloat3(d_forcesBuffer));
507 }
508
509 void GpuHaloExchange::communicateHaloCoordinates(const matrix          box,
510                                                  GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
511 {
512     impl_->communicateHaloCoordinates(box, coordinatesReadyOnDeviceEvent);
513 }
514
515 void GpuHaloExchange::communicateHaloForces(bool accumulateForces)
516 {
517     impl_->communicateHaloForces(accumulateForces);
518 }
519
520 GpuEventSynchronizer* GpuHaloExchange::getForcesReadyOnDeviceEvent()
521 {
522     return impl_->getForcesReadyOnDeviceEvent();
523 }
524 } // namespace gmx