Fix random typos
[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.h"
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                          std::vector<gmx::RVec>* pmeCpuForceBuffer,
65                          const DeviceContext&    deviceContext,
66                          const DeviceStream&     deviceStream) :
67     deviceContext_(deviceContext),
68     pmePpCommStream_(deviceStream),
69     comm_(comm),
70     pmeRank_(pmeRank),
71     pmeCpuForceBuffer_(pmeCpuForceBuffer),
72     d_pmeForces_(nullptr)
73 {
74 }
75
76 PmePpCommGpu::Impl::~Impl() = default;
77
78 void PmePpCommGpu::Impl::reinit(int size)
79 {
80     // Reallocate device buffer used for staging PME force
81     reallocateDeviceBuffer(&d_pmeForces_, size, &d_pmeForcesSize_, &d_pmeForcesSizeAlloc_, deviceContext_);
82
83     // This rank will access PME rank memory directly, so needs to receive the remote PME buffer addresses.
84 #if GMX_MPI
85
86     if (GMX_THREAD_MPI)
87     {
88         // receive device coordinate buffer address from PME rank
89         MPI_Recv(&remotePmeXBuffer_, sizeof(float3*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
90         // send host and device force buffer addresses to PME rank
91         MPI_Send(&d_pmeForces_, sizeof(float3*), MPI_BYTE, pmeRank_, 0, comm_);
92         RVec* pmeCpuForceBufferData = pmeCpuForceBuffer_->data();
93         MPI_Send(&pmeCpuForceBufferData, sizeof(RVec*), MPI_BYTE, pmeRank_, 0, comm_);
94         // Receive address of event and associated flag from PME rank, to allow sync to local stream after force transfer
95         // NOLINTNEXTLINE(bugprone-sizeof-expression)
96         MPI_Recv(&remotePmeForceSendEvent_, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
97         MPI_Recv(&remotePmeForceSendEventRecorded_, sizeof(std::atomic<bool>*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
98     }
99
100 #endif
101 }
102
103 void PmePpCommGpu::Impl::receiveForceFromPmeCudaDirect(bool receivePmeForceToGpu)
104 {
105 #if GMX_MPI
106     // Wait until remote PME task has pushed data, and then enqueue remote event to local stream.
107
108     // Spin until PME rank sets flag
109     while (!(remotePmeForceSendEventRecorded_->load(std::memory_order_acquire))) {};
110
111     // Enqueue remote event
112     remotePmeForceSendEvent_->enqueueWaitEvent(pmePpCommStream_);
113
114     // Reset the flag
115     remotePmeForceSendEventRecorded_->store(false, std::memory_order_release);
116
117     if (receivePmeForceToGpu)
118     {
119         // Record event to be enqueued in the GPU local buffer operations, to
120         // satisfy dependency on receiving the PME force data before
121         // reducing it with the other force contributions.
122         forcesReadySynchronizer_.markEvent(pmePpCommStream_);
123     }
124     else
125     {
126         // Ensure CPU waits for PME forces to be copied before reducing
127         // them with other forces on the CPU
128         pmePpCommStream_.synchronize();
129     }
130 #endif
131 }
132
133 void PmePpCommGpu::Impl::receiveForceFromPmeCudaMpi(float3* pmeForcePtr, int recvSize)
134 {
135 #if GMX_MPI
136     MPI_Recv(pmeForcePtr, recvSize * DIM, MPI_FLOAT, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
137 #else
138     GMX_UNUSED_VALUE(pmeForcePtr);
139     GMX_UNUSED_VALUE(recvSize);
140 #endif
141 }
142
143 void PmePpCommGpu::Impl::receiveForceFromPme(float3* recvPtr, int recvSize, bool receivePmeForceToGpu)
144 {
145     float3* pmeForcePtr = receivePmeForceToGpu ? asFloat3(d_pmeForces_) : recvPtr;
146     if (GMX_THREAD_MPI)
147     {
148         receiveForceFromPmeCudaDirect(receivePmeForceToGpu);
149     }
150     else
151     {
152         receiveForceFromPmeCudaMpi(pmeForcePtr, recvSize);
153     }
154 }
155
156 void PmePpCommGpu::Impl::sendCoordinatesToPmeCudaDirect(float3*               sendPtr,
157                                                         int                   sendSize,
158                                                         GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
159 {
160     // ensure stream waits until coordinate data is available on device
161     coordinatesReadyOnDeviceEvent->enqueueWaitEvent(pmePpCommStream_);
162
163     cudaError_t stat = cudaMemcpyAsync(remotePmeXBuffer_,
164                                        sendPtr,
165                                        sendSize * DIM * sizeof(float),
166                                        cudaMemcpyDefault,
167                                        pmePpCommStream_.stream());
168     CU_RET_ERR(stat, "cudaMemcpyAsync on Send to PME CUDA direct data transfer failed");
169
170 #if GMX_MPI
171     // Record and send event to allow PME task to sync to above transfer before commencing force calculations
172     pmeCoordinatesSynchronizer_.markEvent(pmePpCommStream_);
173     GpuEventSynchronizer* pmeSync = &pmeCoordinatesSynchronizer_;
174     // NOLINTNEXTLINE(bugprone-sizeof-expression)
175     MPI_Send(&pmeSync, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_);
176 #endif
177 }
178
179 void PmePpCommGpu::Impl::sendCoordinatesToPmeCudaMpi(float3*               sendPtr,
180                                                      int                   sendSize,
181                                                      GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
182 {
183     // ensure coordinate data is available on device before we start transfer
184     coordinatesReadyOnDeviceEvent->waitForEvent();
185
186 #if GMX_MPI
187     float3* sendptr_x = sendPtr;
188
189     MPI_Send(sendptr_x, sendSize * DIM, MPI_FLOAT, pmeRank_, eCommType_COORD_GPU, comm_);
190 #else
191     GMX_UNUSED_VALUE(sendPtr);
192     GMX_UNUSED_VALUE(sendSize);
193 #endif
194 }
195
196 void PmePpCommGpu::Impl::sendCoordinatesToPme(float3*               sendPtr,
197                                               int                   sendSize,
198                                               GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
199 {
200     if (GMX_THREAD_MPI)
201     {
202         sendCoordinatesToPmeCudaDirect(sendPtr, sendSize, coordinatesReadyOnDeviceEvent);
203     }
204     else
205     {
206         sendCoordinatesToPmeCudaMpi(sendPtr, sendSize, coordinatesReadyOnDeviceEvent);
207     }
208 }
209 DeviceBuffer<Float3> PmePpCommGpu::Impl::getGpuForceStagingPtr()
210 {
211     return d_pmeForces_;
212 }
213
214 GpuEventSynchronizer* PmePpCommGpu::Impl::getForcesReadySynchronizer()
215 {
216     if (GMX_THREAD_MPI)
217     {
218         return &forcesReadySynchronizer_;
219     }
220     else
221     {
222         return nullptr;
223     }
224 }
225
226 PmePpCommGpu::PmePpCommGpu(MPI_Comm                comm,
227                            int                     pmeRank,
228                            std::vector<gmx::RVec>* pmeCpuForceBuffer,
229                            const DeviceContext&    deviceContext,
230                            const DeviceStream&     deviceStream) :
231     impl_(new Impl(comm, pmeRank, pmeCpuForceBuffer, deviceContext, deviceStream))
232 {
233 }
234
235 PmePpCommGpu::~PmePpCommGpu() = default;
236
237 void PmePpCommGpu::reinit(int size)
238 {
239     impl_->reinit(size);
240 }
241
242 void PmePpCommGpu::receiveForceFromPme(RVec* recvPtr, int recvSize, bool receivePmeForceToGpu)
243 {
244     impl_->receiveForceFromPme(asFloat3(recvPtr), recvSize, receivePmeForceToGpu);
245 }
246
247 void PmePpCommGpu::sendCoordinatesToPmeFromGpu(DeviceBuffer<RVec>    sendPtr,
248                                                int                   sendSize,
249                                                GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
250 {
251     impl_->sendCoordinatesToPme(asFloat3(sendPtr), sendSize, coordinatesReadyOnDeviceEvent);
252 }
253
254 void PmePpCommGpu::sendCoordinatesToPmeFromCpu(RVec*                 sendPtr,
255                                                int                   sendSize,
256                                                GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
257 {
258     impl_->sendCoordinatesToPme(asFloat3(sendPtr), sendSize, coordinatesReadyOnDeviceEvent);
259 }
260
261 DeviceBuffer<Float3> PmePpCommGpu::getGpuForceStagingPtr()
262 {
263     return impl_->getGpuForceStagingPtr();
264 }
265
266 GpuEventSynchronizer* PmePpCommGpu::getForcesReadySynchronizer()
267 {
268     return impl_->getForcesReadySynchronizer();
269 }
270
271 } // namespace gmx