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"
53 #include "gromacs/gpu_utils/gpueventsynchronizer.h"
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
56 #include "gromacs/utility/enumerationhelpers.h"
63 class StatePropagatorDataGpu::Impl
69 /*! \brief Constructor
71 * The buffers are reallocated only at the reinit call, the padding is
72 * used there for the coordinates buffer. It is needed for PME and added at
73 * the end of the buffer. It is assumed that if the rank has PME duties on the
74 * GPU, all coordinates are copied to the GPU and hence, for this rank, the
75 * coordinates buffer is not split into local and non-local ranges. For other
76 * ranks, the padding size is zero. This works because only one rank ever does
77 * PME work on the GPU, and if that rank also does PP work that is the only
78 * rank. So all coordinates are always transferred.
80 * In OpenCL, only pmeStream is used since it is the only stream created in
81 * PME context. The local and non-local streams are only needed when buffer
82 * ops are offloaded. This feature is currently not available in OpenCL and
83 * hence these streams are not set in these builds.
85 * \param[in] deviceStreamManager Object that owns the DeviceContext and DeviceStreams.
86 * \param[in] transferKind H2D/D2H transfer call behavior (synchronous or not).
87 * \param[in] allocationBlockSizeDivisor Determines the padding size for coordinates buffer.
88 * \param[in] wcycle Wall cycle counter data.
90 Impl(const DeviceStreamManager& deviceStreamManager,
91 GpuApiCallBehavior transferKind,
92 int allocationBlockSizeDivisor,
93 gmx_wallcycle* wcycle);
95 /*! \brief Constructor to use in PME-only rank and in tests.
97 * This constructor should be used if only a coordinate device buffer should be managed
98 * using a single stream. Any operation on force or velocity buffer as well as copy of
99 * non-local coordinates will exit with assertion failure. Note, that the pmeStream can
100 * not be a nullptr and the constructor will exit with an assertion failure.
102 * \todo Currently, unsupported copy operations are blocked by assertion that the stream
103 * not nullptr. This should be improved.
105 * \param[in] pmeStream Device PME stream, nullptr is not allowed.
106 * \param[in] deviceContext Device context, nullptr allowed for non-OpenCL builds.
107 * \param[in] transferKind H2D/D2H transfer call behavior (synchronous or not).
108 * \param[in] allocationBlockSizeDivisor Determines the padding size for coordinates buffer.
109 * \param[in] wcycle Wall cycle counter data.
111 Impl(const DeviceStream* pmeStream,
112 const DeviceContext& deviceContext,
113 GpuApiCallBehavior transferKind,
114 int allocationBlockSizeDivisor,
115 gmx_wallcycle* wcycle);
120 /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
122 * Reallocates coordinate, velocities and force buffers on the device.
125 * The coordinates buffer is (re)allocated, when required by PME, with a padding,
126 * the size of which is set by the constructor. The padding region clearing kernel
127 * is scheduled in the \p pmeStream_ (unlike the coordinates H2D) as only the PME
128 * task uses this padding area.
131 * The force buffer is cleared if its size increases, so that previously unused
132 * memory is cleared before forces are accumulated.
134 * \param[in] numAtomsLocal Number of atoms in local domain.
135 * \param[in] numAtomsAll Total number of atoms to handle.
137 void reinit(int numAtomsLocal, int numAtomsAll);
139 /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
141 * \todo There are at least three versions of the function with this functionality in the code:
142 * this one and two more in NBNXM. These should be unified in a shape of a general function
145 * \param[in] atomLocality If all, local or non-local ranges are needed.
147 * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
149 std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality) const;
152 /*! \brief Get the positions buffer on the GPU.
154 * \returns GPU positions buffer.
156 DeviceBuffer<RVec> getCoordinates();
158 /*! \brief Copy positions to the GPU memory.
160 * \param[in] h_x Positions in the host memory.
161 * \param[in] atomLocality Locality of the particles to copy.
163 void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
165 /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
167 * Returns the event synchronizer which indicates that the coordinates are ready for the
168 * consumption on the device. Takes into account that the producer may be different.
170 * If the update is offloaded, and the current step is not a DD/search step, the returned
171 * synchronizer indicates the completion of GPU update-constraint kernels. Otherwise, on search
172 * steps and if update is not offloaded, the coordinates are provided by the H2D copy and the
173 * returned synchronizer indicates that the copy is complete.
175 * \param[in] atomLocality Locality of the particles to wait for.
176 * \param[in] simulationWork The simulation lifetime flags.
177 * \param[in] stepWork The step lifetime flags.
178 * \param[in] gpuCoordinateHaloLaunched Event recorded when GPU coordinate halo has been launched.
180 * \returns The event to synchronize the stream that consumes coordinates on device.
182 GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality atomLocality,
183 const SimulationWorkload& simulationWork,
184 const StepWorkload& stepWork,
185 GpuEventSynchronizer* gpuCoordinateHaloLaunched = nullptr);
187 /*! \brief Blocking wait until coordinates are copied to the device.
189 * Synchronizes the stream in which the copy was executed.
191 * \param[in] atomLocality Locality of the particles to wait for.
193 void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
195 /*! \brief Consume the event for copying coordinates to the device.
197 * Used for manual event consumption. Does nothing except changing the internal event counter.
199 * \param[in] atomLocality Locality of the particles.
201 void consumeCoordinatesCopiedToDeviceEvent(AtomLocality atomLocality);
203 /*! \brief Reset the event for copying coordinates to the device.
205 * Used for manual event consumption. Does nothing except resetting the event.
207 * \param[in] atomLocality Locality of the particles.
209 void resetCoordinatesCopiedToDeviceEvent(AtomLocality atomLocality);
211 /*! \brief Setter for the event synchronizer for the update is done on th GPU
213 * \param[in] xUpdatedOnDeviceEvent The event to synchronize the stream coordinates wre updated on device.
215 void setXUpdatedOnDeviceEvent(GpuEventSynchronizer* xUpdatedOnDeviceEvent);
217 /*! \brief Copy positions from the GPU memory, with an optional explicit dependency.
219 * \param[in] h_x Positions buffer in the host memory.
220 * \param[in] atomLocality Locality of the particles to copy.
221 * \param[in] dependency Dependency event for this operation.
223 void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x,
224 AtomLocality atomLocality,
225 GpuEventSynchronizer* dependency = nullptr);
227 /*! \brief Wait until coordinates are available on the host.
229 * \param[in] atomLocality Locality of the particles to wait for.
231 void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
234 /*! \brief Get the velocities buffer on the GPU.
236 * \returns GPU velocities buffer.
238 DeviceBuffer<RVec> getVelocities();
240 /*! \brief Copy velocities to the GPU memory.
242 * \param[in] h_v Velocities in the host memory.
243 * \param[in] atomLocality Locality of the particles to copy.
245 void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
247 /*! \brief Copy velocities from the GPU memory.
249 * \param[in] h_v Velocities buffer in the host memory.
250 * \param[in] atomLocality Locality of the particles to copy.
252 void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
254 /*! \brief Wait until velocities are available on the host.
256 * \param[in] atomLocality Locality of the particles to wait for.
258 void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
261 /*! \brief Get the force buffer on the GPU.
263 * \returns GPU force buffer.
265 DeviceBuffer<RVec> getForces();
267 /*! \brief Copy forces to the GPU memory.
269 * \param[in] h_f Forces in the host memory.
270 * \param[in] atomLocality Locality of the particles to copy.
272 void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
274 /*! \brief Clear forces in the GPU memory.
276 * \param[in] atomLocality Locality of the particles to clear.
277 * \param[in] dependency Dependency event for this operation.
279 void clearForcesOnGpu(AtomLocality atomLocality, GpuEventSynchronizer* dependency);
281 /*! \brief Get the event synchronizer for the forces ready on device.
283 * Returns either of the event synchronizers, depending on the offload scenario
284 * for the current simulation timestep:
285 * 1. The forces are copied to the device (when GPU buffer ops are off)
286 * 2. The forces are reduced on the device (GPU buffer ops are on)
288 * \param[in] stepWork Step workload flags
289 * \param[in] simulationWork Simulation workload flags
291 * \returns The event to synchronize the stream that consumes forces on device.
293 GpuEventSynchronizer* getLocalForcesReadyOnDeviceEvent(StepWorkload stepWork,
294 SimulationWorkload simulationWork);
296 /*! \brief Getter for the event synchronizer for when forces are reduced on the GPU.
298 * \param[in] atomLocality Locality of the particles to wait for.
299 * \returns The event to mark when forces are reduced on the GPU.
301 GpuEventSynchronizer* fReducedOnDevice(AtomLocality atomLocality);
303 //! \brief Consume the event for when the forces are reduced on device.
304 void consumeForcesReducedOnDeviceEvent(AtomLocality atomLocality);
306 /*! \brief Getter for the event synchronizer for the forces are ready for GPU update.
308 * \param[in] atomLocality Locality of the particles to wait for.
309 * \returns The event to mark when forces are ready for GPU update.
311 GpuEventSynchronizer* fReadyOnDevice(AtomLocality atomLocality);
313 /*! \brief Copy forces from the GPU memory.
315 * \param[in] h_f Forces buffer in the host memory.
316 * \param[in] atomLocality Locality of the particles to copy.
318 void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
320 /*! \brief Wait until forces are available on the host.
322 * \param[in] atomLocality Locality of the particles to wait for.
324 void waitForcesReadyOnHost(AtomLocality atomLocality);
326 /*! \brief Getter for the update stream.
328 * \todo This is temporary here, until the management of this stream is taken over.
330 * \returns The device command stream to use in update-constraints.
332 const DeviceStream* getUpdateStream();
334 /*! \brief Getter for the number of local atoms.
336 * \returns The number of local atoms.
338 int numAtomsLocal() const;
340 /*! \brief Getter for the total number of atoms.
342 * \returns The total number of atoms.
344 int numAtomsAll() const;
348 const DeviceStream* pmeStream_;
349 //! GPU NBNXM local stream.
350 const DeviceStream* localStream_;
351 //! GPU NBNXM non-local stream.
352 const DeviceStream* nonLocalStream_;
353 //! GPU Update-constraints stream.
354 const DeviceStream* updateStream_;
356 // Streams to use for coordinates H2D and D2H copies (one event for each atom locality)
357 EnumerationArray<AtomLocality, const DeviceStream*> xCopyStreams_ = { { nullptr } };
358 // Streams to use for velocities H2D and D2H copies (one event for each atom locality)
359 EnumerationArray<AtomLocality, const DeviceStream*> vCopyStreams_ = { { nullptr } };
360 // Streams to use for forces H2D and D2H copies (one event for each atom locality)
361 EnumerationArray<AtomLocality, const DeviceStream*> fCopyStreams_ = { { nullptr } };
362 // Streams internal to this module
363 std::unique_ptr<DeviceStream> copyInStream_;
364 std::unique_ptr<DeviceStream> memsetStream_;
366 /*! \brief An array of events that indicate H2D copy is complete (one event for each atom locality)
368 * \todo Reconsider naming. It should be xCopiedToDevice or xH2DCopyComplete, etc.
370 EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnDevice_;
371 //! A pointer to an event that the coordinates are ready after update-constraints execution
372 GpuEventSynchronizer* xUpdatedOnDeviceEvent_ = nullptr;
373 //! An array of events that indicate D2H copy of coordinates is complete (one event for each atom locality)
374 EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnHost_;
376 //! An array of events that indicate D2H copy of velocities is complete (one event for each atom locality)
377 EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnHost_;
379 //! An array of events that indicate H2D copy of forces is complete (one event for each atom locality)
380 EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnDevice_;
381 //! An array of events that indicate the forces were reduced on the GPU (one event for each atom locality)
382 EnumerationArray<AtomLocality, GpuEventSynchronizer> fReducedOnDevice_;
383 //! An array of events that indicate D2H copy of forces is complete (one event for each atom locality)
384 EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnHost_;
386 //! GPU context (for OpenCL builds)
387 const DeviceContext& deviceContext_;
388 //! Default GPU calls behavior
389 GpuApiCallBehavior transferKind_ = GpuApiCallBehavior::Async;
390 //! Required minimum divisor of the allocation size of the coordinates buffer
391 int allocationBlockSizeDivisor_ = 0;
393 //! Number of local atoms
394 int numAtomsLocal_ = -1;
395 //! Total number of atoms
396 int numAtomsAll_ = -1;
398 //! Device positions buffer
399 DeviceBuffer<RVec> d_x_;
400 //! Number of particles saved in the positions buffer
402 //! Allocation size for the positions buffer
403 int d_xCapacity_ = -1;
405 //! Device velocities buffer
406 DeviceBuffer<RVec> d_v_;
407 //! Number of particles saved in the velocities buffer
409 //! Allocation size for the velocities buffer
410 int d_vCapacity_ = -1;
412 //! Device force buffer
413 DeviceBuffer<RVec> d_f_;
414 //! Number of particles saved in the force buffer
416 //! Allocation size for the force buffer
417 int d_fCapacity_ = -1;
419 //! \brief Pointer to wallcycle structure.
420 gmx_wallcycle* wcycle_;
422 //! Whether this instance of the class is used on a PME-only rank
423 bool isPmeOnly_ = false;
425 /*! \brief Performs the copy of data from host to device buffer.
427 * \todo Template on locality.
429 * \param[out] d_data Device-side buffer.
430 * \param[in] h_data Host-side buffer.
431 * \param[in] dataSize Device-side data allocation size.
432 * \param[in] atomLocality If all, local or non-local ranges should be copied.
433 * \param[in] deviceStream GPU stream to execute copy in.
435 void copyToDevice(DeviceBuffer<RVec> d_data,
436 gmx::ArrayRef<const gmx::RVec> h_data,
438 AtomLocality atomLocality,
439 const DeviceStream& deviceStream);
441 /*! \brief Performs the copy of data from device to host buffer.
443 * \param[out] h_data Host-side buffer.
444 * \param[in] d_data Device-side buffer.
445 * \param[in] dataSize Device-side data allocation size.
446 * \param[in] atomLocality If all, local or non-local ranges should be copied.
447 * \param[in] deviceStream GPU stream to execute copy in.
449 void copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
450 DeviceBuffer<RVec> d_data,
452 AtomLocality atomLocality,
453 const DeviceStream& deviceStream);
455 /*! \brief Performs the clearing of data in device buffer.
457 * \todo Template on locality.
459 * \param[out] d_data Device-side buffer.
460 * \param[in] dataSize Device-side data allocation size.
461 * \param[in] atomLocality If all, local or non-local ranges should be cleared.
462 * \param[in] deviceStream GPU stream to execute copy in.
464 void clearOnDevice(DeviceBuffer<RVec> d_data,
466 AtomLocality atomLocality,
467 const DeviceStream& deviceStream) const;
472 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H