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