Apply clang-tidy-11 fixes to CUDA files
[alexxy/gromacs.git] / src / gromacs / ewald / pme_pp_comm_gpu_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 PME-PP communication using CUDA
38  *
39  *
40  * \author Alan Gray <alang@nvidia.com>
41  *
42  * \ingroup module_ewald
43  */
44 #include "gmxpre.h"
45
46 #include "gromacs/ewald/pme_pp_communication.h"
47 #include "pme_pp_comm_gpu_impl.h"
48
49 #include "config.h"
50
51 #include "gromacs/gpu_utils/cudautils.cuh"
52 #include "gromacs/gpu_utils/device_context.h"
53 #include "gromacs/gpu_utils/device_stream.h"
54 #include "gromacs/gpu_utils/devicebuffer.h"
55 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
56 #include "gromacs/gpu_utils/typecasts.cuh"
57 #include "gromacs/utility/gmxmpi.h"
58
59 namespace gmx
60 {
61
62 PmePpCommGpu::Impl::Impl(MPI_Comm             comm,
63                          int                  pmeRank,
64                          const DeviceContext& deviceContext,
65                          const DeviceStream&  deviceStream) :
66     deviceContext_(deviceContext),
67     pmePpCommStream_(deviceStream),
68     comm_(comm),
69     pmeRank_(pmeRank),
70     d_pmeForces_(nullptr)
71 {
72 }
73
74 PmePpCommGpu::Impl::~Impl() = default;
75
76 void PmePpCommGpu::Impl::reinit(int size)
77 {
78     // This rank will access PME rank memory directly, so needs to receive the remote PME buffer addresses.
79 #if GMX_MPI
80
81     if (GMX_THREAD_MPI)
82     {
83         // receive device buffer address from PME rank
84         MPI_Recv(&remotePmeXBuffer_, sizeof(float3*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
85         MPI_Recv(&remotePmeFBuffer_, sizeof(float3*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
86     }
87
88 #endif
89
90     // Reallocate buffer used for staging PME force on GPU
91     reallocateDeviceBuffer(&d_pmeForces_, size, &d_pmeForcesSize_, &d_pmeForcesSizeAlloc_, deviceContext_);
92 }
93
94 void PmePpCommGpu::Impl::receiveForceFromPmeCudaDirect(float3* pmeForcePtr, int recvSize, bool receivePmeForceToGpu)
95 {
96 #if GMX_MPI
97     // Receive event from PME task and add to stream, to ensure pull of data doesn't
98     // occur before PME force calc is completed
99     GpuEventSynchronizer* pmeSync;
100     // NOLINTNEXTLINE(bugprone-sizeof-expression)
101     MPI_Recv(&pmeSync, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
102     pmeSync->enqueueWaitEvent(pmePpCommStream_);
103 #endif
104
105     // Pull force data from remote GPU
106     cudaError_t stat = cudaMemcpyAsync(pmeForcePtr,
107                                        remotePmeFBuffer_,
108                                        recvSize * DIM * sizeof(float),
109                                        cudaMemcpyDefault,
110                                        pmePpCommStream_.stream());
111     CU_RET_ERR(stat, "cudaMemcpyAsync on Recv from PME CUDA direct data transfer failed");
112
113     if (receivePmeForceToGpu)
114     {
115         // Record event to be enqueued in the GPU local buffer operations, to
116         // satisfy dependency on receiving the PME force data before
117         // reducing it with the other force contributions.
118         forcesReadySynchronizer_.markEvent(pmePpCommStream_);
119     }
120     else
121     {
122         // Ensure CPU waits for PME forces to be copied before reducing
123         // them with other forces on the CPU
124         pmePpCommStream_.synchronize();
125     }
126 }
127
128 void PmePpCommGpu::Impl::receiveForceFromPmeCudaMpi(float3* pmeForcePtr, int recvSize)
129 {
130 #if GMX_MPI
131     MPI_Recv(pmeForcePtr, recvSize * DIM, MPI_FLOAT, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
132 #else
133     GMX_UNUSED_VALUE(pmeForcePtr);
134     GMX_UNUSED_VALUE(recvSize);
135 #endif
136 }
137
138 void PmePpCommGpu::Impl::receiveForceFromPme(float3* recvPtr, int recvSize, bool receivePmeForceToGpu)
139 {
140     float3* pmeForcePtr = receivePmeForceToGpu ? asFloat3(d_pmeForces_) : recvPtr;
141     if (GMX_THREAD_MPI)
142     {
143         receiveForceFromPmeCudaDirect(pmeForcePtr, recvSize, receivePmeForceToGpu);
144     }
145     else
146     {
147         receiveForceFromPmeCudaMpi(pmeForcePtr, recvSize);
148     }
149 }
150
151 void PmePpCommGpu::Impl::sendCoordinatesToPmeCudaDirect(float3*               sendPtr,
152                                                         int                   sendSize,
153                                                         GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
154 {
155     // ensure stream waits until coordinate data is available on device
156     coordinatesReadyOnDeviceEvent->enqueueWaitEvent(pmePpCommStream_);
157
158     cudaError_t stat = cudaMemcpyAsync(remotePmeXBuffer_,
159                                        sendPtr,
160                                        sendSize * DIM * sizeof(float),
161                                        cudaMemcpyDefault,
162                                        pmePpCommStream_.stream());
163     CU_RET_ERR(stat, "cudaMemcpyAsync on Send to PME CUDA direct data transfer failed");
164
165 #if GMX_MPI
166     // Record and send event to allow PME task to sync to above transfer before commencing force calculations
167     pmeCoordinatesSynchronizer_.markEvent(pmePpCommStream_);
168     GpuEventSynchronizer* pmeSync = &pmeCoordinatesSynchronizer_;
169     // NOLINTNEXTLINE(bugprone-sizeof-expression)
170     MPI_Send(&pmeSync, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_);
171 #endif
172 }
173
174 void PmePpCommGpu::Impl::sendCoordinatesToPmeCudaMpi(float3*               sendPtr,
175                                                      int                   sendSize,
176                                                      GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
177 {
178     // ensure coordinate data is available on device before we start transfer
179     coordinatesReadyOnDeviceEvent->waitForEvent();
180
181 #if GMX_MPI
182     float3* sendptr_x = sendPtr;
183
184     MPI_Send(sendptr_x, sendSize * DIM, MPI_FLOAT, pmeRank_, eCommType_COORD_GPU, comm_);
185 #else
186     GMX_UNUSED_VALUE(sendPtr);
187     GMX_UNUSED_VALUE(sendSize);
188 #endif
189 }
190
191 void PmePpCommGpu::Impl::sendCoordinatesToPme(float3*               sendPtr,
192                                               int                   sendSize,
193                                               GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
194 {
195     if (GMX_THREAD_MPI)
196     {
197         sendCoordinatesToPmeCudaDirect(sendPtr, sendSize, coordinatesReadyOnDeviceEvent);
198     }
199     else
200     {
201         sendCoordinatesToPmeCudaMpi(sendPtr, sendSize, coordinatesReadyOnDeviceEvent);
202     }
203 }
204 DeviceBuffer<Float3> PmePpCommGpu::Impl::getGpuForceStagingPtr()
205 {
206     return d_pmeForces_;
207 }
208
209 GpuEventSynchronizer* PmePpCommGpu::Impl::getForcesReadySynchronizer()
210 {
211     if (GMX_THREAD_MPI)
212     {
213         return &forcesReadySynchronizer_;
214     }
215     else
216     {
217         return nullptr;
218     }
219 }
220
221 PmePpCommGpu::PmePpCommGpu(MPI_Comm             comm,
222                            int                  pmeRank,
223                            const DeviceContext& deviceContext,
224                            const DeviceStream&  deviceStream) :
225     impl_(new Impl(comm, pmeRank, deviceContext, deviceStream))
226 {
227 }
228
229 PmePpCommGpu::~PmePpCommGpu() = default;
230
231 void PmePpCommGpu::reinit(int size)
232 {
233     impl_->reinit(size);
234 }
235
236 void PmePpCommGpu::receiveForceFromPme(RVec* recvPtr, int recvSize, bool receivePmeForceToGpu)
237 {
238     impl_->receiveForceFromPme(asFloat3(recvPtr), recvSize, receivePmeForceToGpu);
239 }
240
241 void PmePpCommGpu::sendCoordinatesToPmeFromGpu(DeviceBuffer<RVec>    sendPtr,
242                                                int                   sendSize,
243                                                GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
244 {
245     impl_->sendCoordinatesToPme(asFloat3(sendPtr), sendSize, coordinatesReadyOnDeviceEvent);
246 }
247
248 void PmePpCommGpu::sendCoordinatesToPmeFromCpu(RVec*                 sendPtr,
249                                                int                   sendSize,
250                                                GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
251 {
252     impl_->sendCoordinatesToPme(asFloat3(sendPtr), sendSize, coordinatesReadyOnDeviceEvent);
253 }
254
255 DeviceBuffer<Float3> PmePpCommGpu::getGpuForceStagingPtr()
256 {
257     return impl_->getGpuForceStagingPtr();
258 }
259
260 GpuEventSynchronizer* PmePpCommGpu::getForcesReadySynchronizer()
261 {
262     return impl_->getForcesReadySynchronizer();
263 }
264
265 } // namespace gmx