e7d05b01f5fb8f5edc9a84045c71278c0bbb007d
[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                          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     }
95
96 #endif
97 }
98
99 // TODO make this asynchronous by splitting into this into
100 // launchRecvForceFromPmeCudaDirect() and sycnRecvForceFromPmeCudaDirect()
101 void PmePpCommGpu::Impl::receiveForceFromPmeCudaDirect(bool receivePmeForceToGpu)
102 {
103 #if GMX_MPI
104     // Remote PME task pushes GPU data directly data to this PP task.
105
106     // Recieve event from PME task after PME->PP force data push has
107     // been scheduled and enqueue this to PP stream.
108     GpuEventSynchronizer* eventptr;
109     // NOLINTNEXTLINE(bugprone-sizeof-expression)
110     MPI_Recv(&eventptr, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
111     eventptr->enqueueWaitEvent(pmePpCommStream_);
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 #endif
127 }
128
129 void PmePpCommGpu::Impl::receiveForceFromPmeCudaMpi(float3* pmeForcePtr, int recvSize)
130 {
131 #if GMX_MPI
132     MPI_Recv(pmeForcePtr, recvSize * DIM, MPI_FLOAT, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
133 #else
134     GMX_UNUSED_VALUE(pmeForcePtr);
135     GMX_UNUSED_VALUE(recvSize);
136 #endif
137 }
138
139 void PmePpCommGpu::Impl::receiveForceFromPme(float3* recvPtr, int recvSize, bool receivePmeForceToGpu)
140 {
141     float3* pmeForcePtr = receivePmeForceToGpu ? asFloat3(d_pmeForces_) : recvPtr;
142     if (GMX_THREAD_MPI)
143     {
144         receiveForceFromPmeCudaDirect(receivePmeForceToGpu);
145     }
146     else
147     {
148         receiveForceFromPmeCudaMpi(pmeForcePtr, recvSize);
149     }
150 }
151
152 void PmePpCommGpu::Impl::sendCoordinatesToPmeCudaDirect(float3*               sendPtr,
153                                                         int                   sendSize,
154                                                         GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
155 {
156     // ensure stream waits until coordinate data is available on device
157     coordinatesReadyOnDeviceEvent->enqueueWaitEvent(pmePpCommStream_);
158
159     cudaError_t stat = cudaMemcpyAsync(remotePmeXBuffer_,
160                                        sendPtr,
161                                        sendSize * DIM * sizeof(float),
162                                        cudaMemcpyDefault,
163                                        pmePpCommStream_.stream());
164     CU_RET_ERR(stat, "cudaMemcpyAsync on Send to PME CUDA direct data transfer failed");
165
166 #if GMX_MPI
167     // Record and send event to allow PME task to sync to above transfer before commencing force calculations
168     pmeCoordinatesSynchronizer_.markEvent(pmePpCommStream_);
169     GpuEventSynchronizer* pmeSync = &pmeCoordinatesSynchronizer_;
170     // NOLINTNEXTLINE(bugprone-sizeof-expression)
171     MPI_Send(&pmeSync, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_);
172 #endif
173 }
174
175 void PmePpCommGpu::Impl::sendCoordinatesToPmeCudaMpi(float3*               sendPtr,
176                                                      int                   sendSize,
177                                                      GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
178 {
179     // ensure coordinate data is available on device before we start transfer
180     coordinatesReadyOnDeviceEvent->waitForEvent();
181
182 #if GMX_MPI
183     float3* sendptr_x = sendPtr;
184
185     MPI_Send(sendptr_x, sendSize * DIM, MPI_FLOAT, pmeRank_, eCommType_COORD_GPU, comm_);
186 #else
187     GMX_UNUSED_VALUE(sendPtr);
188     GMX_UNUSED_VALUE(sendSize);
189 #endif
190 }
191
192 void PmePpCommGpu::Impl::sendCoordinatesToPme(float3*               sendPtr,
193                                               int                   sendSize,
194                                               GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
195 {
196     if (GMX_THREAD_MPI)
197     {
198         sendCoordinatesToPmeCudaDirect(sendPtr, sendSize, coordinatesReadyOnDeviceEvent);
199     }
200     else
201     {
202         sendCoordinatesToPmeCudaMpi(sendPtr, sendSize, coordinatesReadyOnDeviceEvent);
203     }
204 }
205 DeviceBuffer<Float3> PmePpCommGpu::Impl::getGpuForceStagingPtr()
206 {
207     return d_pmeForces_;
208 }
209
210 GpuEventSynchronizer* PmePpCommGpu::Impl::getForcesReadySynchronizer()
211 {
212     if (GMX_THREAD_MPI)
213     {
214         return &forcesReadySynchronizer_;
215     }
216     else
217     {
218         return nullptr;
219     }
220 }
221
222 PmePpCommGpu::PmePpCommGpu(MPI_Comm                comm,
223                            int                     pmeRank,
224                            std::vector<gmx::RVec>* pmeCpuForceBuffer,
225                            const DeviceContext&    deviceContext,
226                            const DeviceStream&     deviceStream) :
227     impl_(new Impl(comm, pmeRank, pmeCpuForceBuffer, deviceContext, deviceStream))
228 {
229 }
230
231 PmePpCommGpu::~PmePpCommGpu() = default;
232
233 void PmePpCommGpu::reinit(int size)
234 {
235     impl_->reinit(size);
236 }
237
238 void PmePpCommGpu::receiveForceFromPme(RVec* recvPtr, int recvSize, bool receivePmeForceToGpu)
239 {
240     impl_->receiveForceFromPme(asFloat3(recvPtr), recvSize, receivePmeForceToGpu);
241 }
242
243 void PmePpCommGpu::sendCoordinatesToPmeFromGpu(DeviceBuffer<RVec>    sendPtr,
244                                                int                   sendSize,
245                                                GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
246 {
247     impl_->sendCoordinatesToPme(asFloat3(sendPtr), sendSize, coordinatesReadyOnDeviceEvent);
248 }
249
250 void PmePpCommGpu::sendCoordinatesToPmeFromCpu(RVec*                 sendPtr,
251                                                int                   sendSize,
252                                                GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
253 {
254     impl_->sendCoordinatesToPme(asFloat3(sendPtr), sendSize, coordinatesReadyOnDeviceEvent);
255 }
256
257 DeviceBuffer<Float3> PmePpCommGpu::getGpuForceStagingPtr()
258 {
259     return impl_->getGpuForceStagingPtr();
260 }
261
262 GpuEventSynchronizer* PmePpCommGpu::getForcesReadySynchronizer()
263 {
264     return impl_->getForcesReadySynchronizer();
265 }
266
267 } // namespace gmx