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/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"
67 #include "domdec_internal.h"
72 //! Number of CUDA threads in a block
73 // TODO Optimize this through experimentation
74 constexpr static int c_threadsPerBlock = 256;
77 __global__ void packSendBufKernel(float3* __restrict__ dataPacked,
78 const float3* __restrict__ data,
79 const int* __restrict__ map,
81 const float3 coordinateShift)
83 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
84 float3* gm_dataDest = &dataPacked[threadIndex];
85 const float3* gm_dataSrc = &data[map[threadIndex]];
87 if (threadIndex < mapSize)
91 *gm_dataDest = *gm_dataSrc + coordinateShift;
95 *gm_dataDest = *gm_dataSrc;
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
107 template<bool accumulate>
108 __global__ void unpackRecvBufKernel(float3* __restrict__ data,
109 const float3* __restrict__ dataPacked,
110 const int* __restrict__ map,
114 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
115 const float3* gm_dataSrc = &dataPacked[threadIndex];
116 float3* gm_dataDest = &data[map[threadIndex]];
118 if (threadIndex < mapSize)
122 *gm_dataDest += *gm_dataSrc;
126 *gm_dataDest = *gm_dataSrc;
133 void GpuHaloExchange::Impl::reinitHalo(float3* d_coordinatesBuffer, float3* d_forcesBuffer)
136 d_x_ = d_coordinatesBuffer;
137 d_f_ = d_forcesBuffer;
139 const gmx_domdec_comm_t& comm = *dd_->comm;
140 const gmx_domdec_comm_dim_t& cd = comm.cd[dimIndex_];
141 const gmx_domdec_ind_t& ind = cd.ind[pulse_];
143 numHomeAtoms_ = comm.atomRanges.numHomeAtoms(); // offset for data recieved by this rank
145 // Determine receive offset for the dimension index and pulse of this halo exchange object
148 int numAtomsTotal = numHomeAtoms_;
149 for (int i = 0; i <= dimIndex_; i++)
151 int pulseMax = (i == dimIndex_) ? pulse_ : (comm.cd[i].numPulses() - 1);
152 for (int p = 0; p <= pulseMax; p++)
154 atomOffset_ = numAtomsTotal;
155 const gmx_domdec_ind_t& indTemp = comm.cd[i].ind[p];
156 numAtomsTotal += indTemp.nrecv[numZoneTemp + 1];
158 numZone = numZoneTemp;
159 numZoneTemp += numZoneTemp;
162 int newSize = ind.nsend[numZone + 1];
164 GMX_ASSERT(cd.receiveInPlace, "Out-of-place receive is not yet supported in GPU halo exchange");
166 // reallocates only if needed
167 h_indexMap_.resize(newSize);
168 // reallocate on device only if needed
169 if (newSize > maxPackedBufferSize_)
171 reallocateDeviceBuffer(&d_indexMap_, newSize, &indexMapSize_, &indexMapSizeAlloc_, deviceContext_);
172 reallocateDeviceBuffer(&d_sendBuf_, newSize, &sendBufSize_, &sendBufSizeAlloc_, deviceContext_);
173 reallocateDeviceBuffer(&d_recvBuf_, newSize, &recvBufSize_, &recvBufSizeAlloc_, deviceContext_);
174 maxPackedBufferSize_ = newSize;
177 xSendSize_ = newSize;
179 MPI_Sendrecv(&xSendSize_, sizeof(int), MPI_BYTE, sendRankX_, 0, &xRecvSize_, sizeof(int),
180 MPI_BYTE, recvRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
182 fSendSize_ = xRecvSize_;
183 fRecvSize_ = xSendSize_;
187 GMX_ASSERT(ind.index.size() == h_indexMap_.size(),
188 "Size mismatch between domain decomposition communication index array and GPU "
189 "halo exchange index mapping array");
190 std::copy(ind.index.begin(), ind.index.end(), h_indexMap_.begin());
192 copyToDeviceBuffer(&d_indexMap_, h_indexMap_.data(), 0, newSize, nonLocalStream_,
193 GpuApiCallBehavior::Async, nullptr);
195 // This rank will push data to its neighbor, so needs to know
196 // the remote receive address and similarly send its receive
197 // address to other neighbour. We can do this here in reinit fn
198 // since the pointers will not change until the next NS step.
200 // Coordinates buffer:
201 void* recvPtr = static_cast<void*>(&d_x_[atomOffset_]);
203 MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankX_, 0, &remoteXPtr_, sizeof(void*),
204 MPI_BYTE, sendRankX_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
207 recvPtr = static_cast<void*>(d_recvBuf_);
208 MPI_Sendrecv(&recvPtr, sizeof(void*), MPI_BYTE, recvRankF_, 0, &remoteFPtr_, sizeof(void*),
209 MPI_BYTE, sendRankF_, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
215 void GpuHaloExchange::Impl::communicateHaloCoordinates(const matrix box,
216 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
221 // ensure stream waits until coordinate data is available on device
222 coordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
225 wallcycle_start(wcycle_, ewcLAUNCH_GPU);
226 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_MOVEX);
228 // launch kernel to pack send buffer
229 KernelLaunchConfig config;
230 config.blockSize[0] = c_threadsPerBlock;
231 config.blockSize[1] = 1;
232 config.blockSize[2] = 1;
233 config.gridSize[0] = (xSendSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
234 config.gridSize[1] = 1;
235 config.gridSize[2] = 1;
236 config.sharedMemorySize = 0;
238 const float3* sendBuf = d_sendBuf_;
239 const float3* d_x = d_x_;
240 const int* indexMap = d_indexMap_;
241 const int size = xSendSize_;
242 // The coordinateShift changes between steps when we have
243 // performed a DD partition, or have updated the box e.g. when
244 // performing pressure coupling. So, for simplicity, the box
245 // is used every step to pass the shift vector as an argument of
246 // the packing kernel.
247 const int boxDimensionIndex = dd_->dim[dimIndex_];
248 const float3 coordinateShift{ box[boxDimensionIndex][XX], box[boxDimensionIndex][YY],
249 box[boxDimensionIndex][ZZ] };
251 // Avoid launching kernel when there is no work to do
254 auto kernelFn = usePBC_ ? packSendBufKernel<true> : packSendBufKernel<false>;
256 const auto kernelArgs = prepareGpuKernelArguments(kernelFn, config, &sendBuf, &d_x,
257 &indexMap, &size, &coordinateShift);
259 launchGpuKernel(kernelFn, config, nonLocalStream_, nullptr,
260 "Domdec GPU Apply X Halo Exchange", kernelArgs);
263 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_MOVEX);
264 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
266 // Consider time spent in communicateHaloData as Comm.X counter
267 // ToDo: We need further refinement here as communicateHaloData includes launch time for cudamemcpyasync
268 wallcycle_start(wcycle_, ewcMOVEX);
270 communicateHaloData(d_x_, HaloQuantity::HaloCoordinates, coordinatesReadyOnDeviceEvent);
272 wallcycle_stop(wcycle_, ewcMOVEX);
277 // The following method should be called after non-local buffer operations,
278 // and before the local buffer operations. It operates in the non-local stream.
279 void GpuHaloExchange::Impl::communicateHaloForces(bool accumulateForces)
281 // Consider time spent in communicateHaloData as Comm.F counter
282 // ToDo: We need further refinement here as communicateHaloData includes launch time for cudamemcpyasync
283 wallcycle_start(wcycle_, ewcMOVEF);
285 // Communicate halo data (in non-local stream)
286 communicateHaloData(d_f_, HaloQuantity::HaloForces, nullptr);
288 wallcycle_stop(wcycle_, ewcMOVEF);
290 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
291 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_MOVEF);
294 // If this is the last pulse and index (noting the force halo
295 // exchanges across multiple pulses and indices are called in
296 // reverse order) then perform the following preparation
298 if ((pulse_ == (dd_->comm->cd[dimIndex_].numPulses() - 1)) && (dimIndex_ == (dd_->ndim - 1)))
300 if (!accumulateForces)
302 // Clear local portion of force array (in local stream)
303 cudaMemsetAsync(d_f, 0, numHomeAtoms_ * sizeof(rvec), localStream_.stream());
306 // ensure non-local stream waits for local stream, due to dependence on
307 // the previous H2D copy of CPU forces (if accumulateForces is true)
308 // or the above clearing.
309 // TODO remove this dependency on localStream - edmine Issue #3093
310 GpuEventSynchronizer eventLocal;
311 eventLocal.markEvent(localStream_);
312 eventLocal.enqueueWaitEvent(nonLocalStream_);
315 // Unpack halo buffer into force array
317 KernelLaunchConfig config;
318 config.blockSize[0] = c_threadsPerBlock;
319 config.blockSize[1] = 1;
320 config.blockSize[2] = 1;
321 config.gridSize[0] = (fRecvSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
322 config.gridSize[1] = 1;
323 config.gridSize[2] = 1;
324 config.sharedMemorySize = 0;
326 const float3* recvBuf = d_recvBuf_;
327 const int* indexMap = d_indexMap_;
328 const int size = fRecvSize_;
330 if (pulse_ > 0 || dd_->ndim > 1)
332 // We need to accumulate rather than set, since it is possible
333 // that, in this pulse/dim, a value could be written to a location
334 // corresponding to the halo region of a following pulse/dim.
335 accumulateForces = true;
340 auto kernelFn = accumulateForces ? unpackRecvBufKernel<true> : unpackRecvBufKernel<false>;
342 const auto kernelArgs =
343 prepareGpuKernelArguments(kernelFn, config, &d_f, &recvBuf, &indexMap, &size);
345 launchGpuKernel(kernelFn, config, nonLocalStream_, nullptr,
346 "Domdec GPU Apply F Halo Exchange", kernelArgs);
351 fReadyOnDevice_.markEvent(nonLocalStream_);
354 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_MOVEF);
355 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
359 void GpuHaloExchange::Impl::communicateHaloData(float3* d_ptr,
360 HaloQuantity haloQuantity,
361 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
370 if (haloQuantity == HaloQuantity::HaloCoordinates)
372 sendPtr = static_cast<void*>(d_sendBuf_);
373 sendSize = xSendSize_;
374 remotePtr = remoteXPtr_;
375 sendRank = sendRankX_;
376 recvRank = recvRankX_;
379 // Wait for event from receiving task that remote coordinates are ready, and enqueue that event to stream used
380 // for subsequent data push. This avoids a race condition with the remote data being written in the previous timestep.
381 // Similarly send event to task that will push data to this task.
382 GpuEventSynchronizer* remoteCoordinatesReadyOnDeviceEvent;
383 MPI_Sendrecv(&coordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*), MPI_BYTE,
384 recvRank, 0, &remoteCoordinatesReadyOnDeviceEvent, sizeof(GpuEventSynchronizer*),
385 MPI_BYTE, sendRank, 0, mpi_comm_mysim_, MPI_STATUS_IGNORE);
386 remoteCoordinatesReadyOnDeviceEvent->enqueueWaitEvent(nonLocalStream_);
388 GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
393 sendPtr = static_cast<void*>(&(d_ptr[atomOffset_]));
394 sendSize = fSendSize_;
395 remotePtr = remoteFPtr_;
396 sendRank = sendRankF_;
397 recvRank = recvRankF_;
400 communicateHaloDataWithCudaDirect(sendPtr, sendSize, sendRank, remotePtr, recvRank);
403 void GpuHaloExchange::Impl::communicateHaloDataWithCudaDirect(void* sendPtr,
412 // We asynchronously push data to remote rank. The remote
413 // destination pointer has already been set in the init fn. We
414 // don't need to worry about overwriting data the remote ranks
415 // still needs since the halo exchange is just done once per
416 // timestep, for each of X and F.
418 // send data to neighbor, if any data exists to send
421 stat = cudaMemcpyAsync(remotePtr, sendPtr, sendSize * DIM * sizeof(float),
422 cudaMemcpyDeviceToDevice, nonLocalStream_.stream());
424 CU_RET_ERR(stat, "cudaMemcpyAsync on GPU Domdec CUDA direct data transfer failed");
428 // ensure pushed data has arrived before remote rank progresses
429 // This rank records an event and sends it to the remote rank which has just been pushed data.
430 // This rank recieves event from remote rank which has pushed data here, and enqueues that event
432 GpuEventSynchronizer* haloDataTransferRemote;
434 haloDataTransferLaunched_->markEvent(nonLocalStream_);
436 MPI_Sendrecv(&haloDataTransferLaunched_, sizeof(GpuEventSynchronizer*), MPI_BYTE, sendRank, 0,
437 &haloDataTransferRemote, sizeof(GpuEventSynchronizer*), MPI_BYTE, recvRank, 0,
438 mpi_comm_mysim_, MPI_STATUS_IGNORE);
440 haloDataTransferRemote->enqueueWaitEvent(nonLocalStream_);
442 GMX_UNUSED_VALUE(sendRank);
443 GMX_UNUSED_VALUE(recvRank);
447 GpuEventSynchronizer* GpuHaloExchange::Impl::getForcesReadyOnDeviceEvent()
449 return &fReadyOnDevice_;
452 /*! \brief Create Domdec GPU object */
453 GpuHaloExchange::Impl::Impl(gmx_domdec_t* dd,
455 MPI_Comm mpi_comm_mysim,
456 const DeviceContext& deviceContext,
457 const DeviceStream& localStream,
458 const DeviceStream& nonLocalStream,
460 gmx_wallcycle* wcycle) :
463 sendRankX_(dd->neighbor[dimIndex][1]),
464 recvRankX_(dd->neighbor[dimIndex][0]),
465 sendRankF_(dd->neighbor[dimIndex][0]),
466 recvRankF_(dd->neighbor[dimIndex][1]),
467 usePBC_(dd->ci[dd->dim[dimIndex]] == 0),
468 haloDataTransferLaunched_(new GpuEventSynchronizer()),
469 mpi_comm_mysim_(mpi_comm_mysim),
470 deviceContext_(deviceContext),
471 localStream_(localStream),
472 nonLocalStream_(nonLocalStream),
477 GMX_RELEASE_ASSERT(GMX_THREAD_MPI,
478 "GPU Halo exchange is currently only supported with thread-MPI enabled");
480 if (usePBC_ && dd->unitCellInfo.haveScrewPBC)
482 gmx_fatal(FARGS, "Error: screw is not yet supported in GPU halo exchange\n");
485 changePinningPolicy(&h_indexMap_, gmx::PinningPolicy::PinnedIfSupported);
487 allocateDeviceBuffer(&d_fShift_, 1, deviceContext_);
490 GpuHaloExchange::Impl::~Impl()
492 freeDeviceBuffer(&d_indexMap_);
493 freeDeviceBuffer(&d_sendBuf_);
494 freeDeviceBuffer(&d_recvBuf_);
495 freeDeviceBuffer(&d_fShift_);
496 delete haloDataTransferLaunched_;
499 GpuHaloExchange::GpuHaloExchange(gmx_domdec_t* dd,
501 MPI_Comm mpi_comm_mysim,
502 const DeviceContext& deviceContext,
503 const DeviceStream& localStream,
504 const DeviceStream& nonLocalStream,
506 gmx_wallcycle* wcycle) :
507 impl_(new Impl(dd, dimIndex, mpi_comm_mysim, deviceContext, localStream, nonLocalStream, pulse, wcycle))
511 GpuHaloExchange::~GpuHaloExchange() = default;
513 void GpuHaloExchange::reinitHalo(DeviceBuffer<RVec> d_coordinatesBuffer, DeviceBuffer<RVec> d_forcesBuffer)
515 impl_->reinitHalo(asFloat3(d_coordinatesBuffer), asFloat3(d_forcesBuffer));
518 void GpuHaloExchange::communicateHaloCoordinates(const matrix box,
519 GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
521 impl_->communicateHaloCoordinates(box, coordinatesReadyOnDeviceEvent);
524 void GpuHaloExchange::communicateHaloForces(bool accumulateForces)
526 impl_->communicateHaloForces(accumulateForces);
529 GpuEventSynchronizer* GpuHaloExchange::getForcesReadyOnDeviceEvent()
531 return impl_->getForcesReadyOnDeviceEvent();