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 Declaration of low-level functions and fields of GPU state propagator object.
39 * \author Artem Zhmurov <zhmurov@gmail.com>
41 * \ingroup module_mdtypes
43 #ifndef GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H
44 #define GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H
52 #include "gromacs/gpu_utils/devicebuffer.h"
54 # include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
56 # include "gromacs/gpu_utils/gpueventsynchronizer_ocl.h"
58 # include "gromacs/gpu_utils/gpueventsynchronizer_sycl.h"
61 #include "gromacs/math/vectypes.h"
62 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
63 #include "gromacs/utility/enumerationhelpers.h"
70 class StatePropagatorDataGpu::Impl
76 /*! \brief Constructor
78 * The buffers are reallocated only at the reinit call, the padding is
79 * used there for the coordinates buffer. It is needed for PME and added at
80 * the end of the buffer. It is assumed that if the rank has PME duties on the
81 * GPU, all coordinates are copied to the GPU and hence, for this rank, the
82 * coordinates buffer is not split into local and non-local ranges. For other
83 * ranks, the padding size is zero. This works because only one rank ever does
84 * PME work on the GPU, and if that rank also does PP work that is the only
85 * rank. So all coordinates are always transferred.
87 * In OpenCL, only pmeStream is used since it is the only stream created in
88 * PME context. The local and non-local streams are only needed when buffer
89 * ops are offloaded. This feature is currently not available in OpenCL and
90 * hence these streams are not set in these builds.
92 * \param[in] deviceStreamManager Object that owns the DeviceContext and DeviceStreams.
93 * \param[in] transferKind H2D/D2H transfer call behavior (synchronous or not).
94 * \param[in] allocationBlockSizeDivisor Determines the padding size for coordinates buffer.
95 * \param[in] wcycle Wall cycle counter data.
97 Impl(const DeviceStreamManager& deviceStreamManager,
98 GpuApiCallBehavior transferKind,
99 int allocationBlockSizeDivisor,
100 gmx_wallcycle* wcycle);
102 /*! \brief Constructor to use in PME-only rank and in tests.
104 * This constructor should be used if only a coordinate device buffer should be managed
105 * using a single stream. Any operation on force or velocity buffer as well as copy of
106 * non-local coordinates will exit with assertion failure. Note, that the pmeStream can
107 * not be a nullptr and the constructor will exit with an assertion failure.
109 * \todo Currently, unsupported copy operations are blocked by assertion that the stream
110 * not nullptr. This should be improved.
112 * \param[in] pmeStream Device PME stream, nullptr is not allowed.
113 * \param[in] deviceContext Device context, nullptr allowed for non-OpenCL builds.
114 * \param[in] transferKind H2D/D2H transfer call behavior (synchronous or not).
115 * \param[in] allocationBlockSizeDivisor Determines the padding size for coordinates buffer.
116 * \param[in] wcycle Wall cycle counter data.
118 Impl(const DeviceStream* pmeStream,
119 const DeviceContext& deviceContext,
120 GpuApiCallBehavior transferKind,
121 int allocationBlockSizeDivisor,
122 gmx_wallcycle* wcycle);
127 /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
129 * Reallocates coordinate, velocities and force buffers on the device.
132 * The coordinates buffer is (re)allocated, when required by PME, with a padding,
133 * the size of which is set by the constructor. The padding region clearing kernel
134 * is scheduled in the \p pmeStream_ (unlike the coordinates H2D) as only the PME
135 * task uses this padding area.
138 * The force buffer is cleared if its size increases, so that previously unused
139 * memory is cleared before forces are accumulated.
141 * \param[in] numAtomsLocal Number of atoms in local domain.
142 * \param[in] numAtomsAll Total number of atoms to handle.
144 void reinit(int numAtomsLocal, int numAtomsAll);
146 /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
148 * \todo There are at least three versions of the function with this functionality in the code:
149 * this one and two more in NBNXM. These should be unified in a shape of a general function
152 * \param[in] atomLocality If all, local or non-local ranges are needed.
154 * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
156 std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality);
159 /*! \brief Get the positions buffer on the GPU.
161 * \returns GPU positions buffer.
163 DeviceBuffer<RVec> getCoordinates();
165 /*! \brief Copy positions to the GPU memory.
167 * \param[in] h_x Positions in the host memory.
168 * \param[in] atomLocality Locality of the particles to copy.
170 void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
172 /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
174 * Returns the event synchronizer which indicates that the coordinates are ready for the
175 * consumption on the device. Takes into account that the producer may be different.
177 * If the update is offloaded, and the current step is not a DD/search step, the returned
178 * synchronizer indicates the completion of GPU update-constraint kernels. Otherwise, on search
179 * steps and if update is not offloaded, the coordinates are provided by the H2D copy and the
180 * returned synchronizer indicates that the copy is complete.
182 * \param[in] atomLocality Locality of the particles to wait for.
183 * \param[in] simulationWork The simulation lifetime flags.
184 * \param[in] stepWork The step lifetime flags.
186 * \returns The event to synchronize the stream that consumes coordinates on device.
188 GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality atomLocality,
189 const SimulationWorkload& simulationWork,
190 const StepWorkload& stepWork);
192 /*! \brief Blocking wait until coordinates are copied to the device.
194 * Synchronizes the stream in which the copy was executed.
196 * \param[in] atomLocality Locality of the particles to wait for.
198 void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
200 /*! \brief Getter for the event synchronizer for the update is done on th GPU
202 * \returns The event to synchronize the stream coordinates wre updated on device.
204 GpuEventSynchronizer* xUpdatedOnDevice();
206 /*! \brief Copy positions from the GPU memory.
208 * \param[in] h_x Positions buffer in the host memory.
209 * \param[in] atomLocality Locality of the particles to copy.
211 void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality);
213 /*! \brief Wait until coordinates are available on the host.
215 * \param[in] atomLocality Locality of the particles to wait for.
217 void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
220 /*! \brief Get the velocities buffer on the GPU.
222 * \returns GPU velocities buffer.
224 DeviceBuffer<RVec> getVelocities();
226 /*! \brief Copy velocities to the GPU memory.
228 * \param[in] h_v Velocities in the host memory.
229 * \param[in] atomLocality Locality of the particles to copy.
231 void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
233 /*! \brief Get the event synchronizer on the H2D velocities copy.
235 * \param[in] atomLocality Locality of the particles to wait for.
237 * \returns The event to synchronize the stream that consumes velocities on device.
239 GpuEventSynchronizer* getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality);
241 /*! \brief Copy velocities from the GPU memory.
243 * \param[in] h_v Velocities buffer in the host memory.
244 * \param[in] atomLocality Locality of the particles to copy.
246 void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
248 /*! \brief Wait until velocities are available on the host.
250 * \param[in] atomLocality Locality of the particles to wait for.
252 void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
255 /*! \brief Get the force buffer on the GPU.
257 * \returns GPU force buffer.
259 DeviceBuffer<RVec> getForces();
261 /*! \brief Copy forces to the GPU memory.
263 * \param[in] h_f Forces in the host memory.
264 * \param[in] atomLocality Locality of the particles to copy.
266 void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
268 /*! \brief Clear forces in the GPU memory.
270 * \param[in] atomLocality Locality of the particles to clear.
272 void clearForcesOnGpu(AtomLocality atomLocality);
274 /*! \brief Get the event synchronizer for the forces ready on device.
276 * Returns either of the event synchronizers, depending on the offload scenario
277 * for the current simulation timestep:
278 * 1. The forces are copied to the device (when GPU buffer ops are off)
279 * 2. The forces are reduced on the device (GPU buffer ops are on)
281 * \todo Pass step workload instead of the useGpuFBufferOps boolean.
283 * \param[in] atomLocality Locality of the particles to wait for.
284 * \param[in] useGpuFBufferOps If the force buffer ops are offloaded to the GPU.
286 * \returns The event to synchronize the stream that consumes forces on device.
288 GpuEventSynchronizer* getForcesReadyOnDeviceEvent(AtomLocality atomLocality, bool useGpuFBufferOps);
290 /*! \brief Getter for the event synchronizer for the forces are reduced on the GPU.
292 * \returns The event to mark when forces are reduced on the GPU.
294 GpuEventSynchronizer* fReducedOnDevice();
296 /*! \brief Copy forces from the GPU memory.
298 * \param[in] h_f Forces buffer in the host memory.
299 * \param[in] atomLocality Locality of the particles to copy.
301 void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
303 /*! \brief Wait until forces are available on the host.
305 * \param[in] atomLocality Locality of the particles to wait for.
307 void waitForcesReadyOnHost(AtomLocality atomLocality);
309 /*! \brief Getter for the update stream.
311 * \todo This is temporary here, until the management of this stream is taken over.
313 * \returns The device command stream to use in update-constraints.
315 const DeviceStream* getUpdateStream();
317 /*! \brief Getter for the number of local atoms.
319 * \returns The number of local atoms.
323 /*! \brief Getter for the total number of atoms.
325 * \returns The total number of atoms.
331 const DeviceStream* pmeStream_;
332 //! GPU NBNXM local stream.
333 const DeviceStream* localStream_;
334 //! GPU NBNXM non-local stream.
335 const DeviceStream* nonLocalStream_;
336 //! GPU Update-constreaints stream.
337 const DeviceStream* updateStream_;
339 // Streams to use for coordinates H2D and D2H copies (one event for each atom locality)
340 EnumerationArray<AtomLocality, const DeviceStream*> xCopyStreams_ = { { nullptr } };
341 // Streams to use for velocities H2D and D2H copies (one event for each atom locality)
342 EnumerationArray<AtomLocality, const DeviceStream*> vCopyStreams_ = { { nullptr } };
343 // Streams to use for forces H2D and D2H copies (one event for each atom locality)
344 EnumerationArray<AtomLocality, const DeviceStream*> fCopyStreams_ = { { nullptr } };
346 /*! \brief An array of events that indicate H2D copy is complete (one event for each atom locality)
348 * \todo Reconsider naming. It should be xCopiedToDevice or xH2DCopyComplete, etc.
350 EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnDevice_;
351 //! An event that the coordinates are ready after update-constraints execution
352 GpuEventSynchronizer xUpdatedOnDevice_;
353 //! An array of events that indicate D2H copy of coordinates is complete (one event for each atom locality)
354 EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnHost_;
356 //! An array of events that indicate H2D copy of velocities is complete (one event for each atom locality)
357 EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnDevice_;
358 //! An array of events that indicate D2H copy of velocities is complete (one event for each atom locality)
359 EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnHost_;
361 //! An array of events that indicate H2D copy of forces is complete (one event for each atom locality)
362 EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnDevice_;
363 //! An event that the forces were reduced on the GPU
364 GpuEventSynchronizer fReducedOnDevice_;
365 //! An array of events that indicate D2H copy of forces is complete (one event for each atom locality)
366 EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnHost_;
368 //! GPU context (for OpenCL builds)
369 const DeviceContext& deviceContext_;
370 //! Default GPU calls behavior
371 GpuApiCallBehavior transferKind_ = GpuApiCallBehavior::Async;
372 //! Required minimum divisor of the allocation size of the coordinates buffer
373 int allocationBlockSizeDivisor_ = 0;
375 //! Number of local atoms
376 int numAtomsLocal_ = -1;
377 //! Total number of atoms
378 int numAtomsAll_ = -1;
380 //! Device positions buffer
381 DeviceBuffer<RVec> d_x_;
382 //! Number of particles saved in the positions buffer
384 //! Allocation size for the positions buffer
385 int d_xCapacity_ = -1;
387 //! Device velocities buffer
388 DeviceBuffer<RVec> d_v_;
389 //! Number of particles saved in the velocities buffer
391 //! Allocation size for the velocities buffer
392 int d_vCapacity_ = -1;
394 //! Device force buffer
395 DeviceBuffer<RVec> d_f_;
396 //! Number of particles saved in the force buffer
398 //! Allocation size for the force buffer
399 int d_fCapacity_ = -1;
401 //! \brief Pointer to wallcycle structure.
402 gmx_wallcycle* wcycle_;
404 /*! \brief Performs the copy of data from host to device buffer.
406 * \todo Template on locality.
408 * \param[out] d_data Device-side buffer.
409 * \param[in] h_data Host-side buffer.
410 * \param[in] dataSize Device-side data allocation size.
411 * \param[in] atomLocality If all, local or non-local ranges should be copied.
412 * \param[in] deviceStream GPU stream to execute copy in.
414 void copyToDevice(DeviceBuffer<RVec> d_data,
415 gmx::ArrayRef<const gmx::RVec> h_data,
417 AtomLocality atomLocality,
418 const DeviceStream& deviceStream);
420 /*! \brief Performs the copy of data from device to host buffer.
422 * \param[out] h_data Host-side buffer.
423 * \param[in] d_data Device-side buffer.
424 * \param[in] dataSize Device-side data allocation size.
425 * \param[in] atomLocality If all, local or non-local ranges should be copied.
426 * \param[in] deviceStream GPU stream to execute copy in.
428 void copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
429 DeviceBuffer<RVec> d_data,
431 AtomLocality atomLocality,
432 const DeviceStream& deviceStream);
434 /*! \brief Performs the clearing of data in device buffer.
436 * \todo Template on locality.
438 * \param[out] d_data Device-side buffer.
439 * \param[in] dataSize Device-side data allocation size.
440 * \param[in] atomLocality If all, local or non-local ranges should be cleared.
441 * \param[in] deviceStream GPU stream to execute copy in.
443 void clearOnDevice(DeviceBuffer<RVec> d_data,
445 AtomLocality atomLocality,
446 const DeviceStream& deviceStream);
451 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H