2 * This file is part of the GROMACS molecular simulation package.
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.
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/pbcutil/ishift.h"
61 #include "gromacs/utility/gmxmpi.h"
63 #include "domdec_internal.h"
68 //! Number of CUDA threads in a block
69 // TODO Optimize this through experimentation
70 constexpr static int c_threadsPerBlock = 256;
73 __global__ void packSendBufKernel(float3* __restrict__ dataPacked,
74 const float3* __restrict__ data,
75 const int* __restrict__ map,
77 const float3 coordinateShift)
79 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
80 float3* gm_dataDest = &dataPacked[threadIndex];
81 const float3* gm_dataSrc = &data[map[threadIndex]];
83 if (threadIndex < mapSize)
87 *gm_dataDest = *gm_dataSrc + coordinateShift;
91 *gm_dataDest = *gm_dataSrc;
98 /*! \brief unpack non-local force data buffer on the GPU using pre-populated "map" containing index
99 * information \param[out] data full array of force values \param[in] dataPacked packed
100 * array of force values to be transferred \param[in] map array of indices defining mapping
101 * from full to packed array \param[in] mapSize number of elements in map array
103 template<bool accumulate>
104 __global__ void unpackRecvBufKernel(float3* __restrict__ data,
105 const float3* __restrict__ dataPacked,
106 const int* __restrict__ map,
110 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
111 const float3* gm_dataSrc = &dataPacked[threadIndex];
112 float3* gm_dataDest = &data[map[threadIndex]];
114 if (threadIndex < mapSize)
118 *gm_dataDest += *gm_dataSrc;
122 *gm_dataDest = *gm_dataSrc;
129 void GpuHaloExchange::Impl::reinitHalo(float3* d_coordinatesBuffer, float3* d_forcesBuffer)
132 d_x_ = d_coordinatesBuffer;
133 d_f_ = d_forcesBuffer;
135 cudaStream_t stream = nonLocalStream_;
137 const gmx_domdec_comm_t& comm = *dd_->comm;
138 const gmx_domdec_comm_dim_t& cd = comm.cd[0];
139 const gmx_domdec_ind_t& ind = cd.ind[0];
140 int newSize = ind.nsend[nzone + 1];
142 GMX_RELEASE_ASSERT(cd.numPulses() == 1,
143 "Multiple pulses are not yet supported in GPU halo exchange");
144 GMX_ASSERT(cd.receiveInPlace, "Out-of-place receive is not yet supported in GPU halo exchange");
146 // reallocates only if needed
147 h_indexMap_.resize(newSize);
148 // reallocate on device only if needed
149 if (newSize > maxPackedBufferSize_)
151 reallocateDeviceBuffer(&d_indexMap_, newSize, &indexMapSize_, &indexMapSizeAlloc_, nullptr);
152 reallocateDeviceBuffer(&d_sendBuf_, newSize, &sendBufSize_, &sendBufSizeAlloc_, nullptr);
153 reallocateDeviceBuffer(&d_recvBuf_, newSize, &recvBufSize_, &recvBufSizeAlloc_, nullptr);
154 maxPackedBufferSize_ = newSize;
157 xSendSize_ = newSize;
159 MPI_Sendrecv(&xSendSize_, sizeof(int), MPI_BYTE, sendRankX_, 0, &xRecvSize_, sizeof(int),
160 MPI_BYTE, recvRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
162 fSendSize_ = xRecvSize_;
163 fRecvSize_ = xSendSize_;
165 numHomeAtoms_ = comm.atomRanges.numHomeAtoms(); // offset for data recieved by this rank
167 GMX_ASSERT(ind.index.size() == h_indexMap_.size(), "Size mismatch");
168 std::copy(ind.index.begin(), ind.index.end(), h_indexMap_.begin());
170 copyToDeviceBuffer(&d_indexMap_, h_indexMap_.data(), 0, newSize, stream,
171 GpuApiCallBehavior::Async, nullptr);
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.
178 // Coordinates buffer:
180 void* recvPtr = static_cast<void*>(&d_coordinatesBuffer[numHomeAtoms_]);
181 MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankX_, 0, &remoteXPtr_, sizeof(void*),
182 MPI_BYTE, sendRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
185 recvPtr = static_cast<void*>(d_recvBuf_);
186 MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankF_, 0, &remoteFPtr_, sizeof(void*),
187 MPI_BYTE, sendRankF_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
194 void GpuHaloExchange::Impl::communicateHaloCoordinates(const matrix box,
195 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
198 // ensure stream waits until coordinate data is available on device
199 coordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
201 // launch kernel to pack send buffer
202 KernelLaunchConfig config;
203 config.blockSize[0] = c_threadsPerBlock;
204 config.blockSize[1] = 1;
205 config.blockSize[2] = 1;
206 config.gridSize[0] = (xSendSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
207 config.gridSize[1] = 1;
208 config.gridSize[2] = 1;
209 config.sharedMemorySize = 0;
210 config.stream = nonLocalStream_;
212 const float3* sendBuf = d_sendBuf_;
213 const float3* d_x = d_x_;
214 const int* indexMap = d_indexMap_;
215 const int size = xSendSize_;
216 // The coordinateShift changes between steps when we have
217 // performed a DD partition, or have updated the box e.g. when
218 // performing pressure coupling. So, for simplicity, the box
219 // is used every step to pass the shift vector as an argument of
220 // the packing kernel.
222 // Because only one-dimensional DD is supported, the coordinate
223 // shift only needs to handle that dimension.
224 const int dimensionIndex = dd_->dim[0];
225 const float3 coordinateShift{ box[dimensionIndex][XX], box[dimensionIndex][YY],
226 box[dimensionIndex][ZZ] };
228 // Avoid launching kernel when there is no work to do
231 auto kernelFn = usePBC_ ? packSendBufKernel<true> : packSendBufKernel<false>;
233 const auto kernelArgs = prepareGpuKernelArguments(kernelFn, config, &sendBuf, &d_x,
234 &indexMap, &size, &coordinateShift);
236 launchGpuKernel(kernelFn, config, nullptr, "Domdec GPU Apply X Halo Exchange", kernelArgs);
239 communicateHaloData(d_x_, HaloQuantity::HaloCoordinates, coordinatesReadyOnDeviceEvent);
244 // The following method should be called after non-local buffer operations,
245 // and before the local buffer operations. It operates in the non-local stream.
246 void GpuHaloExchange::Impl::communicateHaloForces(bool accumulateForces)
249 // Communicate halo data (in non-local stream)
250 communicateHaloData(d_f_, HaloQuantity::HaloForces, nullptr);
254 if (!accumulateForces)
256 // Clear local portion of force array (in local stream)
257 cudaMemsetAsync(d_f, 0, numHomeAtoms_ * sizeof(rvec), localStream_);
260 // ensure non-local stream waits for local stream, due to dependence on
261 // the previous H2D copy of CPU forces (if accumulateForces is true)
262 // or the above clearing.
263 // TODO remove this dependency on localStream - edmine issue #3093
264 GpuEventSynchronizer eventLocal;
265 eventLocal.markEvent(localStream_);
266 eventLocal.enqueueWaitEvent(nonLocalStream_);
268 // Unpack halo buffer into force array
270 KernelLaunchConfig config;
271 config.blockSize[0] = c_threadsPerBlock;
272 config.blockSize[1] = 1;
273 config.blockSize[2] = 1;
274 config.gridSize[0] = (fRecvSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
275 config.gridSize[1] = 1;
276 config.gridSize[2] = 1;
277 config.sharedMemorySize = 0;
278 config.stream = nonLocalStream_;
280 const float3* recvBuf = d_recvBuf_;
281 const int* indexMap = d_indexMap_;
282 const int size = fRecvSize_;
286 auto kernelFn = accumulateForces ? unpackRecvBufKernel<true> : unpackRecvBufKernel<false>;
288 const auto kernelArgs =
289 prepareGpuKernelArguments(kernelFn, config, &d_f, &recvBuf, &indexMap, &size);
291 launchGpuKernel(kernelFn, config, nullptr, "Domdec GPU Apply F Halo Exchange", kernelArgs);
293 fReadyOnDevice_.markEvent(nonLocalStream_);
297 void GpuHaloExchange::Impl::communicateHaloData(float3* d_ptr,
298 HaloQuantity haloQuantity,
299 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
308 if (haloQuantity == HaloQuantity::HaloCoordinates)
310 sendPtr = static_cast<void*>(d_sendBuf_);
311 sendSize = xSendSize_;
312 remotePtr = remoteXPtr_;
313 sendRank = sendRankX_;
314 recvRank = recvRankX_;
317 // Wait for event from receiving task that remote coordinates are ready, and enqueue that event to stream used
318 // for subsequent data push. This avoids a race condition with the remote data being written in the previous timestep.
319 // Similarly send event to task that will push data to this task.
320 GpuEventSynchronizer* remoteCoordinatesReadyOnDeviceEvent;
321 MPI_Sendrecv(&coordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*), MPI_BYTE,
322 recvRank, 0, &remoteCoordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*),
323 MPI_BYTE, sendRank, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
324 remoteCoordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
326 GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
331 sendPtr = static_cast<void*>(&(d_ptr[numHomeAtoms_]));
332 sendSize = fSendSize_;
333 remotePtr = remoteFPtr_;
334 sendRank = sendRankF_;
335 recvRank = recvRankF_;
338 communicateHaloDataWithCudaDirect(sendPtr, sendSize, sendRank, remotePtr, recvRank);
341 void GpuHaloExchange::Impl::communicateHaloDataWithCudaDirect(void* sendPtr,
349 cudaStream_t stream = nonLocalStream_;
351 // We asynchronously push data to remote rank. The remote
352 // destination pointer has already been set in the init fn. We
353 // don't need to worry about overwriting data the remote ranks
354 // still needs since the halo exchange is just done once per
355 // timestep, for each of X and F.
357 // send data to neighbor, if any data exists to send
360 stat = cudaMemcpyAsync(remotePtr, sendPtr, sendSize * DIM * sizeof(float),
361 cudaMemcpyDeviceToDevice, stream);
362 CU_RET_ERR(stat, "cudaMemcpyAsync on GPU Domdec CUDA direct data transfer failed");
366 // ensure pushed data has arrived before remote rank progresses
367 // This rank records an event and sends it to the remote rank which has just been pushed data.
368 // This rank recieves event from remote rank which has pushed data here, and enqueues that event
370 GpuEventSynchronizer* haloDataTransferRemote;
372 haloDataTransferLaunched_->markEvent(stream);
374 MPI_Sendrecv(&haloDataTransferLaunched_, sizeof(GpuEventSynchronizer*), MPI_BYTE, sendRank, 0,
375 &haloDataTransferRemote, sizeof(GpuEventSynchronizer*), MPI_BYTE, recvRank, 0,
376 mpi_comm_mysim_, MPI_STATUS_IGNORE);
378 haloDataTransferRemote->enqueueWaitEvent(stream);
380 GMX_UNUSED_VALUE(sendRank);
381 GMX_UNUSED_VALUE(recvRank);
385 GpuEventSynchronizer* GpuHaloExchange::Impl::getForcesReadyOnDeviceEvent()
387 return &fReadyOnDevice_;
390 /*! \brief Create Domdec GPU object */
391 GpuHaloExchange::Impl::Impl(gmx_domdec_t* dd, MPI_Comm mpi_comm_mysim, void* localStream, void* nonLocalStream) :
393 sendRankX_(dd->neighbor[0][1]),
394 recvRankX_(dd->neighbor[0][0]),
395 sendRankF_(dd->neighbor[0][0]),
396 recvRankF_(dd->neighbor[0][1]),
397 usePBC_(dd->ci[dd->dim[0]] == 0),
398 haloDataTransferLaunched_(new GpuEventSynchronizer()),
399 mpi_comm_mysim_(mpi_comm_mysim),
400 localStream_(*static_cast<cudaStream_t*>(localStream)),
401 nonLocalStream_(*static_cast<cudaStream_t*>(nonLocalStream))
404 GMX_RELEASE_ASSERT(GMX_THREAD_MPI,
405 "GPU Halo exchange is currently only supported with thread-MPI enabled");
409 gmx_fatal(FARGS, "Error: dd->ndim > 1 is not yet supported in GPU halo exchange");
412 if (usePBC_ && dd->unitCellInfo.haveScrewPBC)
414 gmx_fatal(FARGS, "Error: screw is not yet supported in GPU halo exchange\n");
417 changePinningPolicy(&h_indexMap_, gmx::PinningPolicy::PinnedIfSupported);
419 allocateDeviceBuffer(&d_fShift_, 1, nullptr);
422 GpuHaloExchange::Impl::~Impl()
424 freeDeviceBuffer(&d_indexMap_);
425 freeDeviceBuffer(&d_sendBuf_);
426 freeDeviceBuffer(&d_recvBuf_);
427 freeDeviceBuffer(&d_fShift_);
428 delete haloDataTransferLaunched_;
431 GpuHaloExchange::GpuHaloExchange(gmx_domdec_t* dd, MPI_Comm mpi_comm_mysim, void* localStream, void* nonLocalStream) :
432 impl_(new Impl(dd, mpi_comm_mysim, localStream, nonLocalStream))
436 GpuHaloExchange::~GpuHaloExchange() = default;
438 void GpuHaloExchange::reinitHalo(DeviceBuffer<float> d_coordinatesBuffer, DeviceBuffer<float> d_forcesBuffer)
440 impl_->reinitHalo(reinterpret_cast<float3*>(d_coordinatesBuffer),
441 reinterpret_cast<float3*>(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();