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