Remove MPI comm from GPU PME-PP force transfer initiation
[alexxy/gromacs.git] / src / gromacs / ewald / pme_force_sender_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_force_sender_gpu_impl.h"
47
48 #include "config.h"
49
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
52 #include "gromacs/utility/gmxmpi.h"
53
54 namespace gmx
55 {
56
57 /*! \brief Create PME-PP GPU communication object */
58 PmeForceSenderGpu::Impl::Impl(GpuEventSynchronizer*  pmeForcesReady,
59                               MPI_Comm               comm,
60                               const DeviceContext&   deviceContext,
61                               gmx::ArrayRef<PpRanks> ppRanks) :
62     pmeForcesReady_(pmeForcesReady),
63     comm_(comm),
64     ppRanks_(ppRanks),
65     deviceContext_(deviceContext),
66     ppCommStream_(ppRanks.size()),
67     ppCommEvent_(ppRanks.size()),
68     pmeRemoteGpuForcePtr_(ppRanks.size()),
69     pmeRemoteCpuForcePtr_(ppRanks.size())
70 {
71     // Create streams and events to manage pushing of force buffers to remote PP ranks
72     std::unique_ptr<DeviceStream>         stream;
73     std::unique_ptr<GpuEventSynchronizer> event;
74     size_t                                i = 0;
75     for (i = 0; i < ppRanks_.size(); i++)
76     {
77         stream = std::make_unique<DeviceStream>(deviceContext_, DeviceStreamPriority::High, false);
78         ppCommStream_[i] = std::move(stream);
79         event            = std::make_unique<GpuEventSynchronizer>();
80         ppCommEvent_[i]  = std::move(event);
81     }
82 }
83
84 PmeForceSenderGpu::Impl::~Impl() = default;
85
86 /*! \brief Sets location of force to be sent to each PP rank  */
87 void PmeForceSenderGpu::Impl::setForceSendBuffer(DeviceBuffer<Float3> d_f)
88 {
89
90     // Need to send address to PP rank only for thread-MPI as PP rank pulls
91     // data using cudamemcpy
92     if (!GMX_THREAD_MPI)
93     {
94         return;
95     }
96
97 #if GMX_MPI
98
99     if (localForcePtr_.empty())
100     {
101         localForcePtr_.resize(ppRanks_.size());
102     }
103     int ind_start = 0;
104     int ind_end   = 0;
105     int i         = 0;
106     for (const auto& receiver : ppRanks_)
107     {
108         ind_start = ind_end;
109         ind_end   = ind_start + receiver.numAtoms;
110
111         localForcePtr_[i] = &d_f[ind_start];
112         // NOLINTNEXTLINE(bugprone-sizeof-expression)
113         MPI_Recv(&pmeRemoteGpuForcePtr_[i], sizeof(float3*), MPI_BYTE, receiver.rankId, 0, comm_, MPI_STATUS_IGNORE);
114         // NOLINTNEXTLINE(bugprone-sizeof-expression)
115         MPI_Recv(&pmeRemoteCpuForcePtr_[i], sizeof(float3*), MPI_BYTE, receiver.rankId, 0, comm_, MPI_STATUS_IGNORE);
116         i++;
117     }
118
119 #else
120     GMX_UNUSED_VALUE(d_f);
121 #endif
122 }
123
124
125 /*! \brief Send PME synchronizer directly using CUDA memory copy */
126 void PmeForceSenderGpu::Impl::sendFToPpCudaDirect(int ppRank, int numAtoms, bool sendForcesDirectToPpGpu)
127 {
128
129     GMX_ASSERT(GMX_THREAD_MPI, "sendFToPpCudaDirect is expected to be called only for Thread-MPI");
130
131
132 #if GMX_MPI
133     float3* pmeRemoteForcePtr =
134             sendForcesDirectToPpGpu ? pmeRemoteGpuForcePtr_[ppRank] : pmeRemoteCpuForcePtr_[ppRank];
135
136     pmeForcesReady_->enqueueWaitEvent(*ppCommStream_[ppRank]);
137
138     cudaError_t stat = cudaMemcpyAsync(pmeRemoteForcePtr,
139                                        localForcePtr_[ppRank],
140                                        numAtoms * sizeof(rvec),
141                                        cudaMemcpyDefault,
142                                        ppCommStream_[ppRank]->stream());
143     CU_RET_ERR(stat, "cudaMemcpyAsync on Recv from PME CUDA direct data transfer failed");
144     ppCommEvent_[ppRank]->markEvent(*ppCommStream_[ppRank]);
145     // NOLINTNEXTLINE(bugprone-sizeof-expression)
146     MPI_Send(&ppCommEvent_[ppRank], sizeof(GpuEventSynchronizer*), MPI_BYTE, ppRank, 0, comm_);
147 #else
148     GMX_UNUSED_VALUE(ppRank);
149     GMX_UNUSED_VALUE(numAtoms);
150 #endif
151 }
152
153 /*! \brief Send PME data directly using CUDA-aware MPI */
154 void PmeForceSenderGpu::Impl::sendFToPpCudaMpi(DeviceBuffer<RVec> sendbuf,
155                                                int                offset,
156                                                int                numBytes,
157                                                int                ppRank,
158                                                MPI_Request*       request)
159 {
160     GMX_ASSERT(GMX_LIB_MPI, "sendFToPpCudaMpi is expected to be called only for Lib-MPI");
161
162 #if GMX_MPI
163     // if using GPU direct comm with CUDA-aware MPI, make sure forces are ready on device
164     // before sending it to PP ranks
165     pmeForcesReady_->waitForEvent();
166
167     MPI_Isend(sendbuf[offset], numBytes, MPI_BYTE, ppRank, 0, comm_, request);
168
169 #else
170     GMX_UNUSED_VALUE(sendbuf);
171     GMX_UNUSED_VALUE(offset);
172     GMX_UNUSED_VALUE(numBytes);
173     GMX_UNUSED_VALUE(ppRank);
174     GMX_UNUSED_VALUE(request);
175 #endif
176 }
177
178 PmeForceSenderGpu::PmeForceSenderGpu(GpuEventSynchronizer*  pmeForcesReady,
179                                      MPI_Comm               comm,
180                                      const DeviceContext&   deviceContext,
181                                      gmx::ArrayRef<PpRanks> ppRanks) :
182     impl_(new Impl(pmeForcesReady, comm, deviceContext, ppRanks))
183 {
184 }
185
186 PmeForceSenderGpu::~PmeForceSenderGpu() = default;
187
188
189 void PmeForceSenderGpu::setForceSendBuffer(DeviceBuffer<RVec> d_f)
190 {
191     impl_->setForceSendBuffer(d_f);
192 }
193
194 void PmeForceSenderGpu::sendFToPpCudaMpi(DeviceBuffer<RVec> sendbuf,
195                                          int                offset,
196                                          int                numBytes,
197                                          int                ppRank,
198                                          MPI_Request*       request)
199 {
200     impl_->sendFToPpCudaMpi(sendbuf, offset, numBytes, ppRank, request);
201 }
202
203 void PmeForceSenderGpu::sendFToPpCudaDirect(int ppRank, int numAtoms, bool sendForcesDirectToPpGpu)
204 {
205     impl_->sendFToPpCudaDirect(ppRank, numAtoms, sendForcesDirectToPpGpu);
206 }
207
208
209 } // namespace gmx