cb9e787c446836d27549d3e1c6775017520bc07f
[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 "pme_pp_comm_gpu_impl.h"
47
48 #include "config.h"
49
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/gpu_utils/device_context.h"
52 #include "gromacs/gpu_utils/device_stream.h"
53 #include "gromacs/gpu_utils/devicebuffer.h"
54 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
55 #include "gromacs/utility/gmxmpi.h"
56
57 namespace gmx
58 {
59
60 PmePpCommGpu::Impl::Impl(MPI_Comm             comm,
61                          int                  pmeRank,
62                          const DeviceContext& deviceContext,
63                          const DeviceStream&  deviceStream) :
64     deviceContext_(deviceContext),
65     pmePpCommStream_(deviceStream),
66     comm_(comm),
67     pmeRank_(pmeRank),
68     d_pmeForces_(nullptr)
69 {
70     GMX_RELEASE_ASSERT(
71             GMX_THREAD_MPI,
72             "PME-PP GPU Communication is currently only supported with thread-MPI enabled");
73 }
74
75 PmePpCommGpu::Impl::~Impl() = default;
76
77 void PmePpCommGpu::Impl::reinit(int size)
78 {
79     // This rank will access PME rank memory directly, so needs to receive the remote PME buffer addresses.
80 #if GMX_MPI
81     MPI_Recv(&remotePmeXBuffer_, sizeof(void**), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
82     MPI_Recv(&remotePmeFBuffer_, sizeof(void**), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
83
84     // Reallocate buffer used for staging PME force on GPU
85     reallocateDeviceBuffer(&d_pmeForces_, size, &d_pmeForcesSize_, &d_pmeForcesSizeAlloc_, deviceContext_);
86 #else
87     GMX_UNUSED_VALUE(size);
88 #endif
89     return;
90 }
91
92 // TODO make this asynchronous by splitting into this into
93 // launchRecvForceFromPmeCudaDirect() and sycnRecvForceFromPmeCudaDirect()
94 void PmePpCommGpu::Impl::receiveForceFromPmeCudaDirect(void* recvPtr, 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     MPI_Recv(&pmeSync, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_, MPI_STATUS_IGNORE);
101     pmeSync->enqueueWaitEvent(pmePpCommStream_);
102
103     // Pull force data from remote GPU
104     void*       pmeForcePtr = receivePmeForceToGpu ? static_cast<void*>(d_pmeForces_) : recvPtr;
105     cudaError_t stat        = cudaMemcpyAsync(pmeForcePtr,
106                                        remotePmeFBuffer_,
107                                        recvSize * DIM * sizeof(float),
108                                        cudaMemcpyDefault,
109                                        pmePpCommStream_.stream());
110     CU_RET_ERR(stat, "cudaMemcpyAsync on Recv from PME CUDA direct data transfer failed");
111
112     if (receivePmeForceToGpu)
113     {
114         // Record event to be enqueued in the GPU local buffer operations, to
115         // satisfy dependency on receiving the PME force data before
116         // reducing it with the other force contributions.
117         forcesReadySynchronizer_.markEvent(pmePpCommStream_);
118     }
119     else
120     {
121         // Ensure CPU waits for PME forces to be copied before reducing
122         // them with other forces on the CPU
123         cudaStreamSynchronize(pmePpCommStream_.stream());
124     }
125 #else
126     GMX_UNUSED_VALUE(recvPtr);
127     GMX_UNUSED_VALUE(recvSize);
128     GMX_UNUSED_VALUE(receivePmeForceToGpu);
129 #endif
130 }
131
132 void PmePpCommGpu::Impl::sendCoordinatesToPmeCudaDirect(void* sendPtr,
133                                                         int   sendSize,
134                                                         bool gmx_unused sendPmeCoordinatesFromGpu,
135                                                         GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
136 {
137 #if GMX_MPI
138     // ensure stream waits until coordinate data is available on device
139     coordinatesReadyOnDeviceEvent->enqueueWaitEvent(pmePpCommStream_);
140
141     cudaError_t stat = cudaMemcpyAsync(remotePmeXBuffer_,
142                                        sendPtr,
143                                        sendSize * DIM * sizeof(float),
144                                        cudaMemcpyDefault,
145                                        pmePpCommStream_.stream());
146     CU_RET_ERR(stat, "cudaMemcpyAsync on Send to PME CUDA direct data transfer failed");
147
148     // Record and send event to allow PME task to sync to above transfer before commencing force calculations
149     pmeCoordinatesSynchronizer_.markEvent(pmePpCommStream_);
150     GpuEventSynchronizer* pmeSync = &pmeCoordinatesSynchronizer_;
151     MPI_Send(&pmeSync, sizeof(GpuEventSynchronizer*), MPI_BYTE, pmeRank_, 0, comm_);
152 #else
153     GMX_UNUSED_VALUE(sendPtr);
154     GMX_UNUSED_VALUE(sendSize);
155     GMX_UNUSED_VALUE(sendPmeCoordinatesFromGpu);
156     GMX_UNUSED_VALUE(coordinatesReadyOnDeviceEvent);
157 #endif
158 }
159
160 DeviceBuffer<Float3> PmePpCommGpu::Impl::getGpuForceStagingPtr()
161 {
162     return d_pmeForces_;
163 }
164
165 GpuEventSynchronizer* PmePpCommGpu::Impl::getForcesReadySynchronizer()
166 {
167     return &forcesReadySynchronizer_;
168 }
169
170 PmePpCommGpu::PmePpCommGpu(MPI_Comm             comm,
171                            int                  pmeRank,
172                            const DeviceContext& deviceContext,
173                            const DeviceStream&  deviceStream) :
174     impl_(new Impl(comm, pmeRank, deviceContext, deviceStream))
175 {
176 }
177
178 PmePpCommGpu::~PmePpCommGpu() = default;
179
180 void PmePpCommGpu::reinit(int size)
181 {
182     impl_->reinit(size);
183 }
184
185 void PmePpCommGpu::receiveForceFromPmeCudaDirect(void* recvPtr, int recvSize, bool receivePmeForceToGpu)
186 {
187     impl_->receiveForceFromPmeCudaDirect(recvPtr, recvSize, receivePmeForceToGpu);
188 }
189
190 void PmePpCommGpu::sendCoordinatesToPmeCudaDirect(void*                 sendPtr,
191                                                   int                   sendSize,
192                                                   bool                  sendPmeCoordinatesFromGpu,
193                                                   GpuEventSynchronizer* coordinatesReadyOnDeviceEvent)
194 {
195     impl_->sendCoordinatesToPmeCudaDirect(
196             sendPtr, sendSize, sendPmeCoordinatesFromGpu, coordinatesReadyOnDeviceEvent);
197 }
198
199 DeviceBuffer<gmx::RVec> PmePpCommGpu::getGpuForceStagingPtr()
200 {
201     return impl_->getGpuForceStagingPtr();
202 }
203
204 GpuEventSynchronizer* PmePpCommGpu::getForcesReadySynchronizer()
205 {
206     return impl_->getForcesReadySynchronizer();
207 }
208
209 } // namespace gmx