2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * \brief Implements PME-PP communication using CUDA
40 * \author Alan Gray <alang@nvidia.com>
42 * \ingroup module_ewald
46 #include "pme_force_sender_gpu_impl.h"
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/gpu_utils/gpueventsynchronizer.h"
52 #include "gromacs/utility/gmxmpi.h"
57 /*! \brief Create PME-PP GPU communication object */
58 PmeForceSenderGpu::Impl::Impl(GpuEventSynchronizer* pmeForcesReady,
60 const DeviceContext& deviceContext,
61 gmx::ArrayRef<PpRanks> ppRanks) :
62 pmeForcesReady_(pmeForcesReady),
65 ppCommStream_(ppRanks.size()),
66 ppCommEvent_(ppRanks.size()),
67 ppCommEventRecorded_(ppRanks.size()),
68 deviceContext_(deviceContext),
69 pmeRemoteCpuForcePtr_(ppRanks.size()),
70 pmeRemoteGpuForcePtr_(ppRanks.size())
72 // Create streams and events to manage pushing of force buffers to remote PP ranks
73 std::unique_ptr<DeviceStream> stream;
74 std::unique_ptr<GpuEventSynchronizer> event;
76 for (i = 0; i < ppRanks_.size(); i++)
78 stream = std::make_unique<DeviceStream>(deviceContext_, DeviceStreamPriority::High, false);
79 ppCommStream_[i] = std::move(stream);
80 event = std::make_unique<GpuEventSynchronizer>();
81 ppCommEvent_[i] = std::move(event);
85 PmeForceSenderGpu::Impl::~Impl() = default;
87 /*! \brief Sets location of force to be sent to each PP rank */
88 void PmeForceSenderGpu::Impl::setForceSendBuffer(DeviceBuffer<Float3> d_f)
91 // Need to send address to PP rank only for thread-MPI as PP rank pulls
92 // data using cudamemcpy
100 if (localForcePtr_.empty())
102 localForcePtr_.resize(ppRanks_.size());
107 for (const auto& receiver : ppRanks_)
110 ind_end = ind_start + receiver.numAtoms;
112 localForcePtr_[i] = &d_f[ind_start];
113 // NOLINTNEXTLINE(bugprone-sizeof-expression)
114 MPI_Recv(&pmeRemoteGpuForcePtr_[i], sizeof(float3*), MPI_BYTE, receiver.rankId, 0, comm_, MPI_STATUS_IGNORE);
115 // NOLINTNEXTLINE(bugprone-sizeof-expression)
116 MPI_Recv(&pmeRemoteCpuForcePtr_[i], sizeof(float3*), MPI_BYTE, receiver.rankId, 0, comm_, MPI_STATUS_IGNORE);
117 // Send address of event and associated flag to PP rank, to allow remote enqueueing
118 // NOLINTNEXTLINE(bugprone-sizeof-expression)
119 MPI_Send(&ppCommEvent_[i], sizeof(GpuEventSynchronizer*), MPI_BYTE, receiver.rankId, 0, comm_);
121 std::atomic<bool>* tmpPpCommEventRecordedPtr =
122 reinterpret_cast<std::atomic<bool>*>(&(ppCommEventRecorded_[i]));
123 tmpPpCommEventRecordedPtr->store(false, std::memory_order_release);
124 // NOLINTNEXTLINE(bugprone-sizeof-expression)
125 MPI_Send(&tmpPpCommEventRecordedPtr, sizeof(std::atomic<bool>*), MPI_BYTE, receiver.rankId, 0, comm_);
130 GMX_UNUSED_VALUE(d_f);
135 /*! \brief Send PME synchronizer directly using CUDA memory copy */
136 void PmeForceSenderGpu::Impl::sendFToPpCudaDirect(int ppRank, int numAtoms, bool sendForcesDirectToPpGpu)
139 GMX_ASSERT(GMX_THREAD_MPI, "sendFToPpCudaDirect is expected to be called only for Thread-MPI");
143 float3* pmeRemoteForcePtr =
144 sendForcesDirectToPpGpu ? pmeRemoteGpuForcePtr_[ppRank] : pmeRemoteCpuForcePtr_[ppRank];
146 pmeForcesReady_->enqueueWaitEvent(*ppCommStream_[ppRank]);
148 cudaError_t stat = cudaMemcpyAsync(pmeRemoteForcePtr,
149 localForcePtr_[ppRank],
150 numAtoms * sizeof(rvec),
152 ppCommStream_[ppRank]->stream());
153 CU_RET_ERR(stat, "cudaMemcpyAsync on Recv from PME CUDA direct data transfer failed");
154 ppCommEvent_[ppRank]->markEvent(*ppCommStream_[ppRank]);
155 std::atomic<bool>* tmpPpCommEventRecordedPtr =
156 reinterpret_cast<std::atomic<bool>*>(&(ppCommEventRecorded_[ppRank]));
157 tmpPpCommEventRecordedPtr->store(true, std::memory_order_release);
159 GMX_UNUSED_VALUE(ppRank);
160 GMX_UNUSED_VALUE(numAtoms);
164 /*! \brief Send PME data directly using CUDA-aware MPI */
165 void PmeForceSenderGpu::Impl::sendFToPpCudaMpi(DeviceBuffer<RVec> sendbuf,
169 MPI_Request* request)
171 GMX_ASSERT(GMX_LIB_MPI, "sendFToPpCudaMpi is expected to be called only for Lib-MPI");
174 // if using GPU direct comm with CUDA-aware MPI, make sure forces are ready on device
175 // before sending it to PP ranks
176 pmeForcesReady_->waitForEvent();
178 MPI_Isend(sendbuf[offset], numBytes, MPI_BYTE, ppRank, 0, comm_, request);
181 GMX_UNUSED_VALUE(sendbuf);
182 GMX_UNUSED_VALUE(offset);
183 GMX_UNUSED_VALUE(numBytes);
184 GMX_UNUSED_VALUE(ppRank);
185 GMX_UNUSED_VALUE(request);
189 PmeForceSenderGpu::PmeForceSenderGpu(GpuEventSynchronizer* pmeForcesReady,
191 const DeviceContext& deviceContext,
192 gmx::ArrayRef<PpRanks> ppRanks) :
193 impl_(new Impl(pmeForcesReady, comm, deviceContext, ppRanks))
197 PmeForceSenderGpu::~PmeForceSenderGpu() = default;
200 void PmeForceSenderGpu::setForceSendBuffer(DeviceBuffer<RVec> d_f)
202 impl_->setForceSendBuffer(d_f);
205 void PmeForceSenderGpu::sendFToPpCudaMpi(DeviceBuffer<RVec> sendbuf,
209 MPI_Request* request)
211 impl_->sendFToPpCudaMpi(sendbuf, offset, numBytes, ppRank, request);
214 void PmeForceSenderGpu::sendFToPpCudaDirect(int ppRank, int numAtoms, bool sendForcesDirectToPpGpu)
216 impl_->sendFToPpCudaDirect(ppRank, numAtoms, sendForcesDirectToPpGpu);