4a62a895324c1a98ea742cf109137e22a22f192c
[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,2021, 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 <utility>
54
55 #include "gromacs/domdec/domdec.h"
56 #include "gromacs/domdec/domdec_struct.h"
57 #include "gromacs/domdec/gpuhaloexchange.h"
58 #include "gromacs/gpu_utils/cudautils.cuh"
59 #include "gromacs/gpu_utils/device_context.h"
60 #include "gromacs/gpu_utils/devicebuffer.h"
61 #include "gromacs/gpu_utils/gpueventsynchronizer.h"
62 #include "gromacs/gpu_utils/typecasts.cuh"
63 #include "gromacs/gpu_utils/vectype_ops.cuh"
64 #include "gromacs/math/vectypes.h"
65 #include "gromacs/pbcutil/ishift.h"
66 #include "gromacs/timing/wallcycle.h"
67 #include "gromacs/utility/gmxmpi.h"
68
69 #include "domdec_internal.h"
70
71 namespace gmx
72 {
73
74 //! Number of CUDA threads in a block
75 // TODO Optimize this through experimentation
76 constexpr static int c_threadsPerBlock = 256;
77
78 template<bool usePBC>
79 __global__ void packSendBufKernel(float3* __restrict__ dataPacked,
80                                   const float3* __restrict__ data,
81                                   const int* __restrict__ map,
82                                   const int    mapSize,
83                                   const float3 coordinateShift)
84 {
85     int           threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
86     float3*       gm_dataDest = &dataPacked[threadIndex];
87     const float3* gm_dataSrc  = &data[map[threadIndex]];
88
89     if (threadIndex < mapSize)
90     {
91         if (usePBC)
92         {
93             *gm_dataDest = *gm_dataSrc + coordinateShift;
94         }
95         else
96         {
97             *gm_dataDest = *gm_dataSrc;
98         }
99     }
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
131 void GpuHaloExchange::Impl::reinitHalo(float3* d_coordinatesBuffer, float3* d_forcesBuffer)
132 {
133     wallcycle_start(wcycle_, WallCycleCounter::Domdec);
134     wallcycle_sub_start(wcycle_, WallCycleSubCounter::DDGpu);
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[dimIndex_];
141     const gmx_domdec_ind_t&      ind  = cd.ind[pulse_];
142
143     numHomeAtoms_ = comm.atomRanges.numHomeAtoms(); // offset for data recieved by this rank
144
145     // Determine receive offset for the dimension index and pulse of this halo exchange object
146     int numZoneTemp   = 1;
147     int numZone       = 0;
148     int numAtomsTotal = numHomeAtoms_;
149     for (int i = 0; i <= dimIndex_; i++)
150     {
151         int pulseMax = (i == dimIndex_) ? pulse_ : (comm.cd[i].numPulses() - 1);
152         for (int p = 0; p <= pulseMax; p++)
153         {
154             atomOffset_                     = numAtomsTotal;
155             const gmx_domdec_ind_t& indTemp = comm.cd[i].ind[p];
156             numAtomsTotal += indTemp.nrecv[numZoneTemp + 1];
157         }
158         numZone = numZoneTemp;
159         numZoneTemp += numZoneTemp;
160     }
161
162     int newSize = ind.nsend[numZone + 1];
163
164     GMX_ASSERT(cd.receiveInPlace, "Out-of-place receive is not yet supported in GPU halo exchange");
165
166     // reallocates only if needed
167     h_indexMap_.resize(newSize);
168     // reallocate on device only if needed
169     if (newSize > maxPackedBufferSize_)
170     {
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;
175     }
176
177     xSendSize_ = newSize;
178 #if GMX_MPI
179     MPI_Sendrecv(&xSendSize_,
180                  sizeof(int),
181                  MPI_BYTE,
182                  sendRankX_,
183                  0,
184                  &xRecvSize_,
185                  sizeof(int),
186                  MPI_BYTE,
187                  recvRankX_,
188                  0,
189                  mpi_comm_mysim_,
190                  MPI_STATUS_IGNORE);
191 #endif
192     fSendSize_ = xRecvSize_;
193     fRecvSize_ = xSendSize_;
194
195     if (newSize > 0)
196     {
197         GMX_ASSERT(ind.index.size() == h_indexMap_.size(),
198                    "Size mismatch between domain decomposition communication index array and GPU "
199                    "halo exchange index mapping array");
200         std::copy(ind.index.begin(), ind.index.end(), h_indexMap_.begin());
201
202         copyToDeviceBuffer(
203                 &d_indexMap_, h_indexMap_.data(), 0, newSize, *haloStream_, GpuApiCallBehavior::Async, nullptr);
204     }
205
206 #if GMX_MPI
207     // Exchange of remote addresses from neighboring ranks is needed only with CUDA-direct as cudamemcpy needs both src/dst pointer
208     // MPI calls such as MPI_send doesn't worry about receiving address, that is taken care by MPI_recv call in neighboring rank
209     if (GMX_THREAD_MPI)
210     {
211         // This rank will push data to its neighbor, so needs to know
212         // the remote receive address and similarly send its receive
213         // address to other neighbour. We can do this here in reinit fn
214         // since the pointers will not change until the next NS step.
215
216         // Coordinates buffer:
217         float3* recvPtr = &d_x_[atomOffset_];
218         MPI_Sendrecv(&recvPtr,
219                      sizeof(void*),
220                      MPI_BYTE,
221                      recvRankX_,
222                      0,
223                      &remoteXPtr_,
224                      sizeof(void*),
225                      MPI_BYTE,
226                      sendRankX_,
227                      0,
228                      mpi_comm_mysim_,
229                      MPI_STATUS_IGNORE);
230
231         // Force buffer:
232         recvPtr = d_recvBuf_;
233         MPI_Sendrecv(&recvPtr,
234                      sizeof(void*),
235                      MPI_BYTE,
236                      recvRankF_,
237                      0,
238                      &remoteFPtr_,
239                      sizeof(void*),
240                      MPI_BYTE,
241                      sendRankF_,
242                      0,
243                      mpi_comm_mysim_,
244                      MPI_STATUS_IGNORE);
245     }
246 #endif
247
248     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::DDGpu);
249     wallcycle_stop(wcycle_, WallCycleCounter::Domdec);
250 }
251
252 void GpuHaloExchange::Impl::enqueueWaitRemoteCoordinatesReadyEvent(GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
253 {
254     GMX_ASSERT(coordinatesReadyOnDeviceEvent != nullptr,
255                "Co-ordinate Halo exchange requires valid co-ordinate ready event");
256
257     // Wait for event from receiving task that remote coordinates are ready, and enqueue that event to stream used
258     // for subsequent data push. This avoids a race condition with the remote data being written in the previous timestep.
259     // Similarly send event to task that will push data to this task.
260     GpuEventSynchronizer* remoteCoordinatesReadyOnDeviceEvent;
261     MPI_Sendrecv(&coordinatesReadyOnDeviceEvent,
262                  sizeof(GpuEventSynchronizer*), //NOLINT(bugprone-sizeof-expression)
263                  MPI_BYTE,
264                  recvRankX_,
265                  0,
266                  &remoteCoordinatesReadyOnDeviceEvent,
267                  sizeof(GpuEventSynchronizer*), //NOLINT(bugprone-sizeof-expression)
268                  MPI_BYTE,
269                  sendRankX_,
270                  0,
271                  mpi_comm_mysim_,
272                  MPI_STATUS_IGNORE);
273     remoteCoordinatesReadyOnDeviceEvent->enqueueWaitEvent(*haloStream_);
274 }
275
276 GpuEventSynchronizer* GpuHaloExchange::Impl::communicateHaloCoordinates(const matrix box,
277                                                                         GpuEventSynchronizer* dependencyEvent)
278 {
279     wallcycle_start(wcycle_, WallCycleCounter::LaunchGpu);
280
281     // ensure stream waits until dependency has been satisfied
282     dependencyEvent->enqueueWaitEvent(*haloStream_);
283
284     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuMoveX);
285
286     // launch kernel to pack send buffer
287     KernelLaunchConfig config;
288     config.blockSize[0]     = c_threadsPerBlock;
289     config.blockSize[1]     = 1;
290     config.blockSize[2]     = 1;
291     config.gridSize[0]      = (xSendSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
292     config.gridSize[1]      = 1;
293     config.gridSize[2]      = 1;
294     config.sharedMemorySize = 0;
295
296     const float3* sendBuf  = d_sendBuf_;
297     const float3* d_x      = d_x_;
298     const int*    indexMap = d_indexMap_;
299     const int     size     = xSendSize_;
300     // The coordinateShift changes between steps when we have
301     // performed a DD partition, or have updated the box e.g. when
302     // performing pressure coupling. So, for simplicity, the box
303     // is used every step to pass the shift vector as an argument of
304     // the packing kernel.
305     const int    boxDimensionIndex = dd_->dim[dimIndex_];
306     const float3 coordinateShift{ box[boxDimensionIndex][XX],
307                                   box[boxDimensionIndex][YY],
308                                   box[boxDimensionIndex][ZZ] };
309
310     // Avoid launching kernel when there is no work to do
311     if (size > 0)
312     {
313         auto kernelFn = usePBC_ ? packSendBufKernel<true> : packSendBufKernel<false>;
314
315         const auto kernelArgs = prepareGpuKernelArguments(
316                 kernelFn, config, &sendBuf, &d_x, &indexMap, &size, &coordinateShift);
317
318         launchGpuKernel(kernelFn, config, *haloStream_, nullptr, "Domdec GPU Apply X Halo Exchange", kernelArgs);
319     }
320
321     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuMoveX);
322     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
323
324     // Consider time spent in communicateHaloData as Comm.X counter
325     // ToDo: We need further refinement here as communicateHaloData includes launch time for cudamemcpyasync
326     wallcycle_start(wcycle_, WallCycleCounter::MoveX);
327
328     // wait for remote co-ordinates is implicit with process-MPI as non-local stream is synchronized before MPI calls
329     // and MPI_Waitall call makes sure both neighboring ranks' non-local stream is synchronized before data transfer is initiated
330     if (GMX_THREAD_MPI && dimIndex_ == 0 && pulse_ == 0)
331     {
332         enqueueWaitRemoteCoordinatesReadyEvent(dependencyEvent);
333     }
334
335     float3* recvPtr = GMX_THREAD_MPI ? remoteXPtr_ : &d_x_[atomOffset_];
336     communicateHaloData(d_sendBuf_, xSendSize_, sendRankX_, recvPtr, xRecvSize_, recvRankX_);
337
338     coordinateHaloLaunched_.markEvent(*haloStream_);
339
340     wallcycle_stop(wcycle_, WallCycleCounter::MoveX);
341
342     return &coordinateHaloLaunched_;
343 }
344
345 // The following method should be called after non-local buffer operations,
346 // and before the local buffer operations.
347 void GpuHaloExchange::Impl::communicateHaloForces(bool accumulateForces,
348                                                   FixedCapacityVector<GpuEventSynchronizer*, 2>* dependencyEvents)
349 {
350
351     // Consider time spent in communicateHaloData as Comm.F counter
352     // ToDo: We need further refinement here as communicateHaloData includes launch time for cudamemcpyasync
353     wallcycle_start(wcycle_, WallCycleCounter::MoveF);
354
355     while (!dependencyEvents->empty())
356     {
357         auto* dependency = dependencyEvents->back();
358         dependency->enqueueWaitEvent(*haloStream_);
359         dependencyEvents->pop_back();
360     }
361
362     float3* recvPtr = GMX_THREAD_MPI ? remoteFPtr_ : d_recvBuf_;
363
364     // Communicate halo data
365     communicateHaloData(&(d_f_[atomOffset_]), fSendSize_, sendRankF_, recvPtr, fRecvSize_, recvRankF_);
366
367     wallcycle_stop(wcycle_, WallCycleCounter::MoveF);
368
369     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
370     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuMoveF);
371
372     float3* d_f = d_f_;
373
374     // Unpack halo buffer into force array
375
376     KernelLaunchConfig config;
377     config.blockSize[0]     = c_threadsPerBlock;
378     config.blockSize[1]     = 1;
379     config.blockSize[2]     = 1;
380     config.gridSize[0]      = (fRecvSize_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
381     config.gridSize[1]      = 1;
382     config.gridSize[2]      = 1;
383     config.sharedMemorySize = 0;
384
385     const float3* recvBuf  = d_recvBuf_;
386     const int*    indexMap = d_indexMap_;
387     const int     size     = fRecvSize_;
388
389     if (pulse_ > 0 || dd_->ndim > 1)
390     {
391         // We need to accumulate rather than set, since it is possible
392         // that, in this pulse/dim, a value could be written to a location
393         // corresponding to the halo region of a following pulse/dim.
394         accumulateForces = true;
395     }
396
397     if (size > 0)
398     {
399         auto kernelFn = accumulateForces ? unpackRecvBufKernel<true> : unpackRecvBufKernel<false>;
400
401         const auto kernelArgs =
402                 prepareGpuKernelArguments(kernelFn, config, &d_f, &recvBuf, &indexMap, &size);
403
404         launchGpuKernel(kernelFn, config, *haloStream_, nullptr, "Domdec GPU Apply F Halo Exchange", kernelArgs);
405     }
406
407     fReadyOnDevice_.markEvent(*haloStream_);
408
409     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuMoveF);
410     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
411 }
412
413 void GpuHaloExchange::Impl::communicateHaloData(float3* sendPtr,
414                                                 int     sendSize,
415                                                 int     sendRank,
416                                                 float3* recvPtr,
417                                                 int     recvSize,
418                                                 int     recvRank)
419 {
420     if (GMX_THREAD_MPI)
421     {
422         // no need to explicitly sync with GMX_THREAD_MPI as all operations are
423         // anyway launched in correct stream
424         communicateHaloDataWithCudaDirect(sendPtr, sendSize, sendRank, recvPtr, recvRank);
425     }
426     else
427     {
428         communicateHaloDataWithCudaMPI(sendPtr, sendSize, sendRank, recvPtr, recvSize, recvRank);
429     }
430 }
431
432 void GpuHaloExchange::Impl::communicateHaloDataWithCudaMPI(float3* sendPtr,
433                                                            int     sendSize,
434                                                            int     sendRank,
435                                                            float3* recvPtr,
436                                                            int     recvSize,
437                                                            int     recvRank)
438 {
439     // no need to wait for haloDataReadyOnDevice event if this rank is not sending any data
440     if (sendSize > 0)
441     {
442         // wait for halo stream to complete all outstanding
443         // activities, to ensure that buffer is up-to-date in GPU memory
444         // before transferring to remote rank
445
446         // ToDo: Replace stream synchronize with event synchronize
447         haloStream_->synchronize();
448     }
449
450     // perform halo exchange directly in device buffers
451 #if GMX_MPI
452     MPI_Request request;
453
454     // recv remote data into halo region
455     MPI_Irecv(recvPtr, recvSize * DIM, MPI_FLOAT, recvRank, 0, mpi_comm_mysim_, &request);
456
457     // send data to remote halo region
458     MPI_Send(sendPtr, sendSize * DIM, MPI_FLOAT, sendRank, 0, mpi_comm_mysim_);
459
460     MPI_Wait(&request, MPI_STATUS_IGNORE);
461 #endif
462 }
463
464 void GpuHaloExchange::Impl::communicateHaloDataWithCudaDirect(float3* sendPtr,
465                                                               int     sendSize,
466                                                               int     sendRank,
467                                                               float3* remotePtr,
468                                                               int     recvRank)
469 {
470
471     cudaError_t stat;
472
473     // We asynchronously push data to remote rank. The remote
474     // destination pointer has already been set in the init fn.  We
475     // don't need to worry about overwriting data the remote ranks
476     // still needs since the halo exchange is just done once per
477     // timestep, for each of X and F.
478
479     // send data to neighbor, if any data exists to send
480     if (sendSize > 0)
481     {
482         stat = cudaMemcpyAsync(remotePtr,
483                                sendPtr,
484                                sendSize * DIM * sizeof(float),
485                                cudaMemcpyDeviceToDevice,
486                                haloStream_->stream());
487
488         CU_RET_ERR(stat, "cudaMemcpyAsync on GPU Domdec CUDA direct data transfer failed");
489     }
490
491 #if GMX_MPI
492     // ensure pushed data has arrived before remote rank progresses
493     // This rank records an event and sends it to the remote rank which has just been pushed data.
494     // This rank recieves event from remote rank which has pushed data here, and enqueues that event
495     // to its stream.
496     GpuEventSynchronizer* haloDataTransferRemote;
497
498     GMX_ASSERT(haloDataTransferLaunched_ != nullptr,
499                "Halo exchange requires valid event to synchronize data transfer initiated in "
500                "remote rank");
501     haloDataTransferLaunched_->markEvent(*haloStream_);
502
503     MPI_Sendrecv(&haloDataTransferLaunched_,
504                  sizeof(GpuEventSynchronizer*), //NOLINT(bugprone-sizeof-expression)
505                  MPI_BYTE,
506                  sendRank,
507                  0,
508                  &haloDataTransferRemote,
509                  sizeof(GpuEventSynchronizer*), //NOLINT(bugprone-sizeof-expression)
510                  MPI_BYTE,
511                  recvRank,
512                  0,
513                  mpi_comm_mysim_,
514                  MPI_STATUS_IGNORE);
515
516     haloDataTransferRemote->enqueueWaitEvent(*haloStream_);
517 #else
518     GMX_UNUSED_VALUE(sendRank);
519     GMX_UNUSED_VALUE(recvRank);
520 #endif
521 }
522
523 GpuEventSynchronizer* GpuHaloExchange::Impl::getForcesReadyOnDeviceEvent()
524 {
525     return &fReadyOnDevice_;
526 }
527
528 /*! \brief Create Domdec GPU object */
529 GpuHaloExchange::Impl::Impl(gmx_domdec_t*        dd,
530                             int                  dimIndex,
531                             MPI_Comm             mpi_comm_mysim,
532                             const DeviceContext& deviceContext,
533                             int                  pulse,
534                             gmx_wallcycle*       wcycle) :
535     dd_(dd),
536     sendRankX_(dd->neighbor[dimIndex][1]),
537     recvRankX_(dd->neighbor[dimIndex][0]),
538     sendRankF_(dd->neighbor[dimIndex][0]),
539     recvRankF_(dd->neighbor[dimIndex][1]),
540     usePBC_(dd->ci[dd->dim[dimIndex]] == 0),
541     haloDataTransferLaunched_(GMX_THREAD_MPI ? new GpuEventSynchronizer() : nullptr),
542     mpi_comm_mysim_(mpi_comm_mysim),
543     deviceContext_(deviceContext),
544     haloStream_(new DeviceStream(deviceContext, DeviceStreamPriority::High, false)),
545     dimIndex_(dimIndex),
546     pulse_(pulse),
547     wcycle_(wcycle)
548 {
549     if (usePBC_ && dd->unitCellInfo.haveScrewPBC)
550     {
551         gmx_fatal(FARGS, "Error: screw is not yet supported in GPU halo exchange\n");
552     }
553
554     changePinningPolicy(&h_indexMap_, gmx::PinningPolicy::PinnedIfSupported);
555
556     allocateDeviceBuffer(&d_fShift_, 1, deviceContext_);
557 }
558
559 GpuHaloExchange::Impl::~Impl()
560 {
561     freeDeviceBuffer(&d_indexMap_);
562     freeDeviceBuffer(&d_sendBuf_);
563     freeDeviceBuffer(&d_recvBuf_);
564     freeDeviceBuffer(&d_fShift_);
565     delete haloDataTransferLaunched_;
566 }
567
568 GpuHaloExchange::GpuHaloExchange(gmx_domdec_t*        dd,
569                                  int                  dimIndex,
570                                  MPI_Comm             mpi_comm_mysim,
571                                  const DeviceContext& deviceContext,
572                                  int                  pulse,
573                                  gmx_wallcycle*       wcycle) :
574     impl_(new Impl(dd, dimIndex, mpi_comm_mysim, deviceContext, pulse, wcycle))
575 {
576 }
577
578 GpuHaloExchange::GpuHaloExchange(GpuHaloExchange&&) noexcept = default;
579
580 GpuHaloExchange& GpuHaloExchange::operator=(GpuHaloExchange&& other) noexcept
581 {
582     std::swap(impl_, other.impl_);
583     return *this;
584 }
585
586 GpuHaloExchange::~GpuHaloExchange() = default;
587
588 void GpuHaloExchange::reinitHalo(DeviceBuffer<RVec> d_coordinatesBuffer, DeviceBuffer<RVec> d_forcesBuffer)
589 {
590     impl_->reinitHalo(asFloat3(d_coordinatesBuffer), asFloat3(d_forcesBuffer));
591 }
592
593 GpuEventSynchronizer* GpuHaloExchange::communicateHaloCoordinates(const matrix          box,
594                                                                   GpuEventSynchronizer* dependencyEvent)
595 {
596     return impl_->communicateHaloCoordinates(box, dependencyEvent);
597 }
598
599 void GpuHaloExchange::communicateHaloForces(bool accumulateForces,
600                                             FixedCapacityVector<GpuEventSynchronizer*, 2>* dependencyEvents)
601 {
602     impl_->communicateHaloForces(accumulateForces, dependencyEvents);
603 }
604
605 GpuEventSynchronizer* GpuHaloExchange::getForcesReadyOnDeviceEvent()
606 {
607     return impl_->getForcesReadyOnDeviceEvent();
608 }
609 } // namespace gmx