2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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 Declaration of interfaces for GPU state data propagator object.
39 * This object stores and manages positions, velocities and forces for
40 * all particles in the system on the GPU.
42 * \todo Add cycle counters.
43 * \todo Add synchronization points.
45 * \author Artem Zhmurov <zhmurov@gmail.com>
48 * \ingroup module_mdtypes
50 #ifndef GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_H
51 #define GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_H
53 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
54 #include "gromacs/gpu_utils/gpu_utils.h"
55 #include "gromacs/math/vectypes.h"
56 #include "gromacs/mdtypes/simulation_workload.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/classhelpers.h"
64 class GpuEventSynchronizer;
69 class DeviceStreamManager;
71 class StatePropagatorDataGpu
74 /*! \brief Constructor
76 * The buffers are reallocated only at the reinit call, the padding is
77 * used there for the coordinates buffer. It is needed for PME and added at
78 * the end of the buffer. It is assumed that if the rank has PME duties on the
79 * GPU, all coordinates are copied to the GPU and hence, for this rank, the
80 * coordinates buffer is not split into local and non-local ranges. For other
81 * ranks, the padding size is zero. This works because only one rank ever does
82 * PME work on the GPU, and if that rank also does PP work that is the only
83 * rank. So all coordinates are always transferred.
85 * In OpenCL, only pmeStream is used since it is the only stream created in
86 * PME context. The local and non-local streams are only needed when buffer
87 * ops are offloaded. This feature is currently not available in OpenCL and
88 * hence these streams are not set in these builds.
90 * \param[in] deviceStreamManager Object that owns the DeviceContext and DeviceStreams.
91 * \param[in] transferKind H2D/D2H transfer call behavior (synchronous or not).
92 * \param[in] allocationBlockSizeDivisor Deterines padding size for coordinates buffer.
93 * \param[in] wcycle Wall cycle counter data.
95 StatePropagatorDataGpu(const DeviceStreamManager& deviceStreamManager,
96 GpuApiCallBehavior transferKind,
97 int allocationBlockSizeDivisor,
98 gmx_wallcycle* wcycle);
100 /*! \brief Constructor to use in PME-only rank and in tests.
102 * This constructor should be used if only a coordinate device buffer should be managed
103 * using a single stream. Any operation on force or velocity buffer as well as copy of
104 * non-local coordinates will exit with assertion failure. Note, that the pmeStream can
105 * not be a nullptr and the constructor will exit with an assertion failure.
107 * \todo Currently, unsupported copy operations are blocked by assertion that the stream
108 * not nullptr. This should be improved.
110 * \param[in] pmeStream Device PME stream, nullptr is not allowed.
111 * \param[in] deviceContext Device context, nullptr allowed for non-OpenCL builds.
112 * \param[in] transferKind H2D/D2H transfer call behavior (synchronous or not).
113 * \param[in] allocationBlockSizeDivisor Determines padding size for coordinates buffer.
114 * \param[in] wcycle Wall cycle counter data.
116 StatePropagatorDataGpu(const DeviceStream* pmeStream,
117 const DeviceContext& deviceContext,
118 GpuApiCallBehavior transferKind,
119 int allocationBlockSizeDivisor,
120 gmx_wallcycle* wcycle);
123 StatePropagatorDataGpu(StatePropagatorDataGpu&& other) noexcept;
125 StatePropagatorDataGpu& operator=(StatePropagatorDataGpu&& other) noexcept;
127 ~StatePropagatorDataGpu();
129 /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
131 * Reallocates coordinate, velocities and force buffers on the device.
134 * The coordinates buffer is (re)allocated, when required by PME, with a padding,
135 * the size of which is set by the constructor. The padding region clearing kernel
136 * is scheduled in the \p pmeStream_ (unlike the coordinates H2D) as only the PME
137 * task uses this padding area.
140 * The force buffer is cleared if its size increases, so that previously unused
141 * memory is cleared before forces are accumulated.
143 * \param[in] numAtomsLocal Number of atoms in local domain.
144 * \param[in] numAtomsAll Total number of atoms to handle.
146 void reinit(int numAtomsLocal, int numAtomsAll);
148 /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
150 * \todo There are at least three versions of the function with this functionality in the code:
151 * this one and two more in NBNXM. These should be unified in a shape of a general function
154 * \param[in] atomLocality If all, local or non-local ranges are needed.
156 * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
158 std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality);
161 /*! \brief Get the positions buffer on the GPU.
163 * \returns GPU positions buffer.
165 DeviceBuffer<RVec> getCoordinates();
167 /*! \brief Copy positions to the GPU memory.
169 * \param[in] h_x Positions in the host memory.
170 * \param[in] atomLocality Locality of the particles to copy.
172 void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
174 /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
176 * Returns the event synchronizer which indicates that the coordinates are ready for the
177 * consumption on the device. Takes into account that the producer may be different.
179 * If the update is offloaded, and the current step is not a DD/search step, the returned
180 * synchronizer indicates the completion of GPU update-constraint kernels. Otherwise, on search
181 * steps and if update is not offloaded, the coordinates are provided by the H2D copy and the
182 * returned synchronizer indicates that the copy is complete.
184 * \param[in] atomLocality Locality of the particles to wait for.
185 * \param[in] simulationWork The simulation lifetime flags.
186 * \param[in] stepWork The step lifetime flags.
188 * \returns The event to synchronize the stream that consumes coordinates on device.
190 GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality atomLocality,
191 const SimulationWorkload& simulationWork,
192 const StepWorkload& stepWork);
194 /*! \brief Blocking wait until coordinates are copied to the device.
196 * Synchronizes the stream in which the copy was executed.
198 * \param[in] atomLocality Locality of the particles to wait for.
200 void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
202 /*! \brief Getter for the event synchronizer for the update is done on th GPU
204 * \returns The event to synchronize the stream coordinates wre updated on device.
206 GpuEventSynchronizer* xUpdatedOnDevice();
208 /*! \brief Copy positions from the GPU memory.
210 * \param[in] h_x Positions buffer in the host memory.
211 * \param[in] atomLocality Locality of the particles to copy.
213 void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality);
215 /*! \brief Wait until coordinates are available on the host.
217 * \param[in] atomLocality Locality of the particles to wait for.
219 void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
222 /*! \brief Get the velocities buffer on the GPU.
224 * \returns GPU velocities buffer.
226 DeviceBuffer<RVec> getVelocities();
228 /*! \brief Copy velocities to the GPU memory.
230 * \param[in] h_v Velocities in the host memory.
231 * \param[in] atomLocality Locality of the particles to copy.
233 void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
235 /*! \brief Get the event synchronizer for the H2D velocities copy.
237 * \param[in] atomLocality Locality of the particles to wait for.
239 * \returns The event to synchronize the stream that consumes velocities on device.
241 GpuEventSynchronizer* getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality);
243 /*! \brief Copy velocities from the GPU memory.
245 * \param[in] h_v Velocities buffer in the host memory.
246 * \param[in] atomLocality Locality of the particles to copy.
248 void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
250 /*! \brief Wait until velocities are available on the host.
252 * \param[in] atomLocality Locality of the particles to wait for.
254 void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
257 /*! \brief Get the force buffer on the GPU.
259 * \returns GPU force buffer.
261 DeviceBuffer<RVec> getForces();
263 /*! \brief Copy forces to the GPU memory.
265 * \param[in] h_f Forces in the host memory.
266 * \param[in] atomLocality Locality of the particles to copy.
268 void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
270 /*! \brief Get the event synchronizer for the forces ready on device.
272 * Returns either of the event synchronizers, depending on the offload scenario
273 * for the current simulation timestep:
274 * 1. The forces are copied to the device (when GPU buffer ops are off)
275 * 2. The forces are reduced on the device (GPU buffer ops are on)
277 * \todo Pass step workload instead of the useGpuFBufferOps boolean.
279 * \param[in] atomLocality Locality of the particles to wait for.
280 * \param[in] useGpuFBufferOps If the force buffer ops are offloaded to the GPU.
282 * \returns The event to synchronize the stream that consumes forces on device.
284 GpuEventSynchronizer* getForcesReadyOnDeviceEvent(AtomLocality atomLocality, bool useGpuFBufferOps);
286 /*! \brief Getter for the event synchronizer for the forces are reduced on the GPU.
288 * \returns The event to mark when forces are reduced on the GPU.
290 GpuEventSynchronizer* fReducedOnDevice();
292 /*! \brief Copy forces from the GPU memory.
294 * \param[in] h_f Forces buffer in the host memory.
295 * \param[in] atomLocality Locality of the particles to copy.
297 void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
299 /*! \brief Wait until forces are available on the host.
301 * \param[in] atomLocality Locality of the particles to wait for.
303 void waitForcesReadyOnHost(AtomLocality atomLocality);
305 /*! \brief Getter for the update stream.
307 * \todo This is temporary here, until the management of this stream is taken over.
309 * \returns The device command stream to use in update-constraints.
311 const DeviceStream* getUpdateStream();
313 /*! \brief Getter for the number of local atoms.
315 * \returns The number of local atoms.
319 /*! \brief Getter for the total number of atoms.
321 * \returns The total number of atoms.
327 gmx::PrivateImplPointer<Impl> impl_;
328 GMX_DISALLOW_COPY_AND_ASSIGN(StatePropagatorDataGpu);
333 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_H