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.
35 /*! \libinternal \file
36 * \brief Declaration of class which receives coordinates to GPU memory on PME task
38 * \author Alan Gray <alang@nvidia.com>
40 * \ingroup module_ewald
42 #ifndef GMX_PMECOORDINATERECEIVERGPU_H
43 #define GMX_PMECOORDINATERECEIVERGPU_H
47 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
48 #include "gromacs/math/vectypes.h"
49 #include "gromacs/utility/gmxmpi.h"
62 class PmeCoordinateReceiverGpu
66 /*! \brief Creates PME GPU coordinate receiver object
68 * For multi-GPU runs, the PME GPU can receive coordinates from
69 * multiple PP GPUs. Data from these distinct communications can
70 * be handled separately in the PME spline/spread kernel, allowing
71 * pipelining which overlaps computation and communication. The
72 * class methods are designed to called seperately for each remote
73 * PP rank, and internally a different stream is used for each
74 * remote PP rank to allow overlapping.
76 * \param[in] comm Communicator used for simulation
77 * \param[in] deviceContext GPU context
78 * \param[in] ppRanks List of PP ranks
80 PmeCoordinateReceiverGpu(MPI_Comm comm, const DeviceContext& deviceContext, gmx::ArrayRef<PpRanks> ppRanks);
81 ~PmeCoordinateReceiverGpu();
84 * Re-initialize: set atom ranges and, for thread-MPI case,
85 * send coordinates buffer address to PP rank
86 * This is required after repartitioning since atom ranges and
87 * buffer allocations may have changed.
88 * \param[in] d_x coordinates buffer in GPU memory
90 void reinitCoordinateReceiver(DeviceBuffer<RVec> d_x);
94 * Receive coordinate synchronizer pointer from the PP ranks.
95 * \param[in] ppRank PP rank to receive the synchronizer from.
97 void receiveCoordinatesSynchronizerFromPpCudaDirect(int ppRank);
100 * Used for lib MPI, receives co-ordinates from PP ranks
101 * \param[in] recvbuf coordinates buffer in GPU memory
102 * \param[in] numAtoms starting element in buffer
103 * \param[in] numBytes number of bytes to transfer
104 * \param[in] ppRank PP rank to send data
106 void launchReceiveCoordinatesFromPpCudaMpi(DeviceBuffer<RVec> recvbuf, int numAtoms, int numBytes, int ppRank);
109 * For lib MPI, wait for coordinates from any PP rank
110 * For thread MPI, enqueue PP co-ordinate transfer event received from PP
111 * rank determined from pipeline stage into given stream
112 * \param[in] pipelineStage stage of pipeline corresponding to this transfer
113 * \param[in] deviceStream stream in which to enqueue the wait event.
114 * \returns rank of sending PP task
116 int synchronizeOnCoordinatesFromPpRank(int pipelineStage, const DeviceStream& deviceStream);
118 /*! \brief Perform above synchronizeOnCoordinatesFromPpRanks for all PP ranks,
119 * enqueueing all events to a single stream
120 * \param[in] deviceStream stream in which to enqueue the wait events.
122 void synchronizeOnCoordinatesFromAllPpRanks(const DeviceStream& deviceStream);
125 * Return pointer to stream associated with specific PP rank sender index
126 * \param[in] senderIndex Index of sender PP rank.
128 DeviceStream* ppCommStream(int senderIndex);
131 * Returns range of atoms involved in communication associated with specific PP rank sender
132 * index \param[in] senderIndex Index of sender PP rank.
134 std::tuple<int, int> ppCommAtomRange(int senderIndex);
137 * Return number of PP ranks involved in PME-PP communication
139 int ppCommNumSenderRanks();
143 std::unique_ptr<Impl> impl_;