2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * \brief Implements GPU halo exchange using CUDA.
40 * \author Alan Gray <alang@nvidia.com>
42 * \ingroup module_domdec
46 #include "gpuhaloexchange_impl.cuh"
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/math/vectypes.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/utility/gmxmpi.h"
64 #include "domdec_internal.h"
69 //! Number of CUDA threads in a block
70 // TODO Optimize this through experimentation
71 constexpr static int c_threadsPerBlock = 256;
74 __global__ void packSendBufKernel(float3* __restrict__ dataPacked,
75 const float3* __restrict__ data,
76 const int* __restrict__ map,
78 const float3 coordinateShift)
80 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
81 float3* gm_dataDest = &dataPacked[threadIndex];
82 const float3* gm_dataSrc = &data[map[threadIndex]];
84 if (threadIndex < mapSize)
88 *gm_dataDest = *gm_dataSrc + coordinateShift;
92 *gm_dataDest = *gm_dataSrc;
99 /*! \brief unpack non-local force data buffer on the GPU using pre-populated "map" containing index
100 * information \param[out] data full array of force values \param[in] dataPacked packed
101 * array of force values to be transferred \param[in] map array of indices defining mapping
102 * from full to packed array \param[in] mapSize number of elements in map array
104 template<bool accumulate>
105 __global__ void unpackRecvBufKernel(float3* __restrict__ data,
106 const float3* __restrict__ dataPacked,
107 const int* __restrict__ map,
111 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
112 const float3* gm_dataSrc = &dataPacked[threadIndex];
113 float3* gm_dataDest = &data[map[threadIndex]];
115 if (threadIndex < mapSize)
119 *gm_dataDest += *gm_dataSrc;
123 *gm_dataDest = *gm_dataSrc;
130 void GpuHaloExchange::Impl::reinitHalo(float3* d_coordinatesBuffer, float3* d_forcesBuffer)
133 d_x_ = d_coordinatesBuffer;
134 d_f_ = d_forcesBuffer;
136 cudaStream_t stream = nonLocalStream_;
138 const gmx_domdec_comm_t& comm = *dd_->comm;
139 const gmx_domdec_comm_dim_t& cd = comm.cd[0];
140 const gmx_domdec_ind_t& ind = cd.ind[0];
141 int newSize = ind.nsend[nzone + 1];
143 GMX_RELEASE_ASSERT(cd.numPulses() == 1,
144 "Multiple pulses are not yet supported in GPU halo exchange");
145 GMX_ASSERT(cd.receiveInPlace, "Out-of-place receive is not yet supported in GPU halo exchange");
147 // reallocates only if needed
148 h_indexMap_.resize(newSize);
149 // reallocate on device only if needed
150 if (newSize > maxPackedBufferSize_)
152 reallocateDeviceBuffer(&d_indexMap_, newSize, &indexMapSize_, &indexMapSizeAlloc_, nullptr);
153 reallocateDeviceBuffer(&d_sendBuf_, newSize, &sendBufSize_, &sendBufSizeAlloc_, nullptr);
154 reallocateDeviceBuffer(&d_recvBuf_, newSize, &recvBufSize_, &recvBufSizeAlloc_, nullptr);
155 maxPackedBufferSize_ = newSize;
158 xSendSize_ = newSize;
160 MPI_Sendrecv(&xSendSize_, sizeof(int), MPI_BYTE, sendRankX_, 0, &xRecvSize_, sizeof(int),
161 MPI_BYTE, recvRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
163 fSendSize_ = xRecvSize_;
164 fRecvSize_ = xSendSize_;
166 numHomeAtoms_ = comm.atomRanges.numHomeAtoms(); // offset for data recieved by this rank
168 GMX_ASSERT(ind.index.size() == h_indexMap_.size(), "Size mismatch");
169 std::copy(ind.index.begin(), ind.index.end(), h_indexMap_.begin());
171 copyToDeviceBuffer(&d_indexMap_, h_indexMap_.data(), 0, newSize, stream,
172 GpuApiCallBehavior::Async, nullptr);
174 // This rank will push data to its neighbor, so needs to know
175 // the remote receive address and similarly send its receive
176 // address to other neighbour. We can do this here in reinit fn
177 // since the pointers will not change until the next NS step.
179 // Coordinates buffer:
181 void* recvPtr = static_cast<void*>(&d_coordinatesBuffer[numHomeAtoms_]);
182 MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankX_, 0, &remoteXPtr_, sizeof(void*),
183 MPI_BYTE, sendRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
186 recvPtr = static_cast<void*>(d_recvBuf_);
187 MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankF_, 0, &remoteFPtr_, sizeof(void*),
188 MPI_BYTE, sendRankF_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
195 void GpuHaloExchange::Impl::communicateHaloCoordinates(const matrix box,
196 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
199 // ensure stream waits until coordinate data is available on device
200 coordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
202 // launch kernel to pack send buffer
203 KernelLaunchConfig config;
204 config.blockSize[0] = c_threadsPerBlock;
205 config.blockSize[1] = 1;
206 config.blockSize[2] = 1;
207 config.gridSize[0] = (xSendSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
208 config.gridSize[1] = 1;
209 config.gridSize[2] = 1;
210 config.sharedMemorySize = 0;
211 config.stream = nonLocalStream_;
213 const float3* sendBuf = d_sendBuf_;
214 const float3* d_x = d_x_;
215 const int* indexMap = d_indexMap_;
216 const int size = xSendSize_;
217 // The coordinateShift changes between steps when we have
218 // performed a DD partition, or have updated the box e.g. when
219 // performing pressure coupling. So, for simplicity, the box
220 // is used every step to pass the shift vector as an argument of
221 // the packing kernel.
223 // Because only one-dimensional DD is supported, the coordinate
224 // shift only needs to handle that dimension.
225 const int dimensionIndex = dd_->dim[0];
226 const float3 coordinateShift{ box[dimensionIndex][XX], box[dimensionIndex][YY],
227 box[dimensionIndex][ZZ] };
229 // Avoid launching kernel when there is no work to do
232 auto kernelFn = usePBC_ ? packSendBufKernel<true> : packSendBufKernel<false>;
234 const auto kernelArgs = prepareGpuKernelArguments(kernelFn, config, &sendBuf, &d_x,
235 &indexMap, &size, &coordinateShift);
237 launchGpuKernel(kernelFn, config, nullptr, "Domdec GPU Apply X Halo Exchange", kernelArgs);
240 communicateHaloData(d_x_, HaloQuantity::HaloCoordinates, coordinatesReadyOnDeviceEvent);
245 // The following method should be called after non-local buffer operations,
246 // and before the local buffer operations. It operates in the non-local stream.
247 void GpuHaloExchange::Impl::communicateHaloForces(bool accumulateForces)
250 // Communicate halo data (in non-local stream)
251 communicateHaloData(d_f_, HaloQuantity::HaloForces, nullptr);
255 if (!accumulateForces)
257 // Clear local portion of force array (in local stream)
258 cudaMemsetAsync(d_f, 0, numHomeAtoms_ * sizeof(rvec), localStream_);
261 // ensure non-local stream waits for local stream, due to dependence on
262 // the previous H2D copy of CPU forces (if accumulateForces is true)
263 // or the above clearing.
264 // TODO remove this dependency on localStream - edmine issue #3093
265 GpuEventSynchronizer eventLocal;
266 eventLocal.markEvent(localStream_);
267 eventLocal.enqueueWaitEvent(nonLocalStream_);
269 // Unpack halo buffer into force array
271 KernelLaunchConfig config;
272 config.blockSize[0] = c_threadsPerBlock;
273 config.blockSize[1] = 1;
274 config.blockSize[2] = 1;
275 config.gridSize[0] = (fRecvSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
276 config.gridSize[1] = 1;
277 config.gridSize[2] = 1;
278 config.sharedMemorySize = 0;
279 config.stream = nonLocalStream_;
281 const float3* recvBuf = d_recvBuf_;
282 const int* indexMap = d_indexMap_;
283 const int size = fRecvSize_;
287 auto kernelFn = accumulateForces ? unpackRecvBufKernel<true> : unpackRecvBufKernel<false>;
289 const auto kernelArgs =
290 prepareGpuKernelArguments(kernelFn, config, &d_f, &recvBuf, &indexMap, &size);
292 launchGpuKernel(kernelFn, config, nullptr, "Domdec GPU Apply F Halo Exchange", kernelArgs);
294 fReadyOnDevice_.markEvent(nonLocalStream_);
298 void GpuHaloExchange::Impl::communicateHaloData(float3* d_ptr,
299 HaloQuantity haloQuantity,
300 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
309 if (haloQuantity == HaloQuantity::HaloCoordinates)
311 sendPtr = static_cast<void*>(d_sendBuf_);
312 sendSize = xSendSize_;
313 remotePtr = remoteXPtr_;
314 sendRank = sendRankX_;
315 recvRank = recvRankX_;
318 // Wait for event from receiving task that remote coordinates are ready, and enqueue that event to stream used
319 // for subsequent data push. This avoids a race condition with the remote data being written in the previous timestep.
320 // Similarly send event to task that will push data to this task.
321 GpuEventSynchronizer* remoteCoordinatesReadyOnDeviceEvent;
322 MPI_Sendrecv(&coordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*), MPI_BYTE,
323 recvRank, 0, &remoteCoordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*),
324 MPI_BYTE, sendRank, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
325 remoteCoordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
327 GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
332 sendPtr = static_cast<void*>(&(d_ptr[numHomeAtoms_]));
333 sendSize = fSendSize_;
334 remotePtr = remoteFPtr_;
335 sendRank = sendRankF_;
336 recvRank = recvRankF_;
339 communicateHaloDataWithCudaDirect(sendPtr, sendSize, sendRank, remotePtr, recvRank);
342 void GpuHaloExchange::Impl::communicateHaloDataWithCudaDirect(void* sendPtr,
350 cudaStream_t stream = nonLocalStream_;
352 // We asynchronously push data to remote rank. The remote
353 // destination pointer has already been set in the init fn. We
354 // don't need to worry about overwriting data the remote ranks
355 // still needs since the halo exchange is just done once per
356 // timestep, for each of X and F.
358 // send data to neighbor, if any data exists to send
361 stat = cudaMemcpyAsync(remotePtr, sendPtr, sendSize * DIM * sizeof(float),
362 cudaMemcpyDeviceToDevice, stream);
363 CU_RET_ERR(stat, "cudaMemcpyAsync on GPU Domdec CUDA direct data transfer failed");
367 // ensure pushed data has arrived before remote rank progresses
368 // This rank records an event and sends it to the remote rank which has just been pushed data.
369 // This rank recieves event from remote rank which has pushed data here, and enqueues that event
371 GpuEventSynchronizer* haloDataTransferRemote;
373 haloDataTransferLaunched_->markEvent(stream);
375 MPI_Sendrecv(&haloDataTransferLaunched_, sizeof(GpuEventSynchronizer*), MPI_BYTE, sendRank, 0,
376 &haloDataTransferRemote, sizeof(GpuEventSynchronizer*), MPI_BYTE, recvRank, 0,
377 mpi_comm_mysim_, MPI_STATUS_IGNORE);
379 haloDataTransferRemote->enqueueWaitEvent(stream);
381 GMX_UNUSED_VALUE(sendRank);
382 GMX_UNUSED_VALUE(recvRank);
386 GpuEventSynchronizer* GpuHaloExchange::Impl::getForcesReadyOnDeviceEvent()
388 return &fReadyOnDevice_;
391 /*! \brief Create Domdec GPU object */
392 GpuHaloExchange::Impl::Impl(gmx_domdec_t* dd, MPI_Comm mpi_comm_mysim, void* localStream, void* nonLocalStream) :
394 sendRankX_(dd->neighbor[0][1]),
395 recvRankX_(dd->neighbor[0][0]),
396 sendRankF_(dd->neighbor[0][0]),
397 recvRankF_(dd->neighbor[0][1]),
398 usePBC_(dd->ci[dd->dim[0]] == 0),
399 haloDataTransferLaunched_(new GpuEventSynchronizer()),
400 mpi_comm_mysim_(mpi_comm_mysim),
401 localStream_(*static_cast<cudaStream_t*>(localStream)),
402 nonLocalStream_(*static_cast<cudaStream_t*>(nonLocalStream))
405 GMX_RELEASE_ASSERT(GMX_THREAD_MPI,
406 "GPU Halo exchange is currently only supported with thread-MPI enabled");
410 gmx_fatal(FARGS, "Error: dd->ndim > 1 is not yet supported in GPU halo exchange");
413 if (usePBC_ && dd->unitCellInfo.haveScrewPBC)
415 gmx_fatal(FARGS, "Error: screw is not yet supported in GPU halo exchange\n");
418 changePinningPolicy(&h_indexMap_, gmx::PinningPolicy::PinnedIfSupported);
420 allocateDeviceBuffer(&d_fShift_, 1, nullptr);
423 GpuHaloExchange::Impl::~Impl()
425 freeDeviceBuffer(&d_indexMap_);
426 freeDeviceBuffer(&d_sendBuf_);
427 freeDeviceBuffer(&d_recvBuf_);
428 freeDeviceBuffer(&d_fShift_);
429 delete haloDataTransferLaunched_;
432 GpuHaloExchange::GpuHaloExchange(gmx_domdec_t* dd, MPI_Comm mpi_comm_mysim, void* localStream, void* nonLocalStream) :
433 impl_(new Impl(dd, mpi_comm_mysim, localStream, nonLocalStream))
437 GpuHaloExchange::~GpuHaloExchange() = default;
439 void GpuHaloExchange::reinitHalo(DeviceBuffer<RVec> d_coordinatesBuffer, DeviceBuffer<RVec> d_forcesBuffer)
441 impl_->reinitHalo(asFloat3(d_coordinatesBuffer), asFloat3(d_forcesBuffer));
444 void GpuHaloExchange::communicateHaloCoordinates(const matrix box,
445 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
447 impl_->communicateHaloCoordinates(box, coordinatesReadyOnDeviceEvent);
450 void GpuHaloExchange::communicateHaloForces(bool accumulateForces)
452 impl_->communicateHaloForces(accumulateForces);
455 GpuEventSynchronizer* GpuHaloExchange::getForcesReadyOnDeviceEvent()
457 return impl_->getForcesReadyOnDeviceEvent();