dbe435d508ddd13c54aeec0eeb7c21d26cae7784
[alexxy/gromacs.git] / src / gromacs / mdtypes / state_propagator_data_gpu_impl.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Declaration of low-level functions and fields of GPU state propagator object.
38  *
39  * \author Artem Zhmurov <zhmurov@gmail.com>
40  *
41  * \ingroup module_mdtypes
42  */
43 #ifndef GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H
44 #define GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H
45
46 #include "gmxpre.h"
47
48 #include "config.h"
49
50 #include <memory>
51
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"
57
58 struct gmx_wallcycle;
59
60 namespace gmx
61 {
62
63 class StatePropagatorDataGpu::Impl
64 {
65 public:
66     Impl();
67
68
69     /*! \brief Constructor
70      *
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.
79      *
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.
84      *
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.
89      */
90     Impl(const DeviceStreamManager& deviceStreamManager,
91          GpuApiCallBehavior         transferKind,
92          int                        allocationBlockSizeDivisor,
93          gmx_wallcycle*             wcycle);
94
95     /*! \brief Constructor to use in PME-only rank and in tests.
96      *
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.
101      *
102      *  \todo Currently, unsupported copy operations are blocked by assertion that the stream
103      *        not nullptr. This should be improved.
104      *
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.
110      */
111     Impl(const DeviceStream*  pmeStream,
112          const DeviceContext& deviceContext,
113          GpuApiCallBehavior   transferKind,
114          int                  allocationBlockSizeDivisor,
115          gmx_wallcycle*       wcycle);
116
117     ~Impl();
118
119
120     /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
121      *
122      * Reallocates coordinate, velocities and force buffers on the device.
123      *
124      * \note
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.
129      *
130      * \note
131      * The force buffer is cleared if its size increases, so that previously unused
132      * memory is cleared before forces are accumulated.
133      *
134      *  \param[in] numAtomsLocal  Number of atoms in local domain.
135      *  \param[in] numAtomsAll    Total number of atoms to handle.
136      */
137     void reinit(int numAtomsLocal, int numAtomsAll);
138
139     /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
140      *
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
143      *       in DD.
144      *
145      * \param[in]  atomLocality    If all, local or non-local ranges are needed.
146      *
147      * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
148      */
149     std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality) const;
150
151
152     /*! \brief Get the positions buffer on the GPU.
153      *
154      *  \returns GPU positions buffer.
155      */
156     DeviceBuffer<RVec> getCoordinates();
157
158     /*! \brief Copy positions to the GPU memory.
159      *
160      *  \param[in] h_x           Positions in the host memory.
161      *  \param[in] atomLocality  Locality of the particles to copy.
162      */
163     void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
164
165     /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
166      *
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.
169      *
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.
174      *
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.
179      *
180      *  \returns  The event to synchronize the stream that consumes coordinates on device.
181      */
182     GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
183                                                            const SimulationWorkload& simulationWork,
184                                                            const StepWorkload&       stepWork,
185                                                            GpuEventSynchronizer* gpuCoordinateHaloLaunched = nullptr);
186
187     /*! \brief Blocking wait until coordinates are copied to the device.
188      *
189      * Synchronizes the stream in which the copy was executed.
190      *
191      *  \param[in] atomLocality  Locality of the particles to wait for.
192      */
193     void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
194
195     /*! \brief Setter for the event synchronizer for the update is done on th GPU
196      *
197      *  \param[in] xUpdatedOnDeviceEvent  The event to synchronize the stream coordinates wre updated on device.
198      */
199     void setXUpdatedOnDeviceEvent(GpuEventSynchronizer* xUpdatedOnDeviceEvent);
200
201     /*! \brief Copy positions from the GPU memory, with an optional explicit dependency.
202      *
203      *  \param[in] h_x           Positions buffer in the host memory.
204      *  \param[in] atomLocality  Locality of the particles to copy.
205      *  \param[in] dependency    Dependency event for this operation.
206      */
207     void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x,
208                                 AtomLocality             atomLocality,
209                                 GpuEventSynchronizer*    dependency = nullptr);
210
211     /*! \brief Wait until coordinates are available on the host.
212      *
213      *  \param[in] atomLocality  Locality of the particles to wait for.
214      */
215     void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
216
217
218     /*! \brief Get the velocities buffer on the GPU.
219      *
220      *  \returns GPU velocities buffer.
221      */
222     DeviceBuffer<RVec> getVelocities();
223
224     /*! \brief Copy velocities to the GPU memory.
225      *
226      *  \param[in] h_v           Velocities in the host memory.
227      *  \param[in] atomLocality  Locality of the particles to copy.
228      */
229     void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
230
231     /*! \brief Copy velocities from the GPU memory.
232      *
233      *  \param[in] h_v           Velocities buffer in the host memory.
234      *  \param[in] atomLocality  Locality of the particles to copy.
235      */
236     void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
237
238     /*! \brief Wait until velocities are available on the host.
239      *
240      *  \param[in] atomLocality  Locality of the particles to wait for.
241      */
242     void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
243
244
245     /*! \brief Get the force buffer on the GPU.
246      *
247      *  \returns GPU force buffer.
248      */
249     DeviceBuffer<RVec> getForces();
250
251     /*! \brief Copy forces to the GPU memory.
252      *
253      *  \param[in] h_f           Forces in the host memory.
254      *  \param[in] atomLocality  Locality of the particles to copy.
255      */
256     void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
257
258     /*! \brief Clear forces in the GPU memory.
259      *
260      *  \param[in] atomLocality  Locality of the particles to clear.
261      *  \param[in] dependency    Dependency event for this operation.
262      */
263     void clearForcesOnGpu(AtomLocality atomLocality, GpuEventSynchronizer* dependency);
264
265     /*! \brief Get the event synchronizer for the forces ready on device.
266      *
267      *  Returns either of the event synchronizers, depending on the offload scenario
268      *  for the current simulation timestep:
269      *  1. The forces are copied to the device (when GPU buffer ops are off)
270      *  2. The forces are reduced on the device (GPU buffer ops are on)
271      *
272      *  \param[in] stepWork        Step workload flags
273      *  \param[in] simulationWork  Simulation workload flags
274      *
275      *  \returns  The event to synchronize the stream that consumes forces on device.
276      */
277     GpuEventSynchronizer* getLocalForcesReadyOnDeviceEvent(StepWorkload       stepWork,
278                                                            SimulationWorkload simulationWork);
279
280     /*! \brief Getter for the event synchronizer for when forces are reduced on the GPU.
281      *
282      *  \param[in] atomLocality      Locality of the particles to wait for.
283      *  \returns                     The event to mark when forces are reduced on the GPU.
284      */
285     GpuEventSynchronizer* fReducedOnDevice(AtomLocality atomLocality);
286
287     /*! \brief Getter for the event synchronizer for the forces are ready for GPU update.
288      *
289      *  \param[in] atomLocality      Locality of the particles to wait for.
290      *  \returns                     The event to mark when forces are ready for GPU update.
291      */
292     GpuEventSynchronizer* fReadyOnDevice(AtomLocality atomLocality);
293
294     /*! \brief Copy forces from the GPU memory.
295      *
296      *  \param[in] h_f           Forces buffer in the host memory.
297      *  \param[in] atomLocality  Locality of the particles to copy.
298      */
299     void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
300
301     /*! \brief Wait until forces are available on the host.
302      *
303      *  \param[in] atomLocality  Locality of the particles to wait for.
304      */
305     void waitForcesReadyOnHost(AtomLocality atomLocality);
306
307     /*! \brief Getter for the update stream.
308      *
309      *  \todo This is temporary here, until the management of this stream is taken over.
310      *
311      *  \returns The device command stream to use in update-constraints.
312      */
313     const DeviceStream* getUpdateStream();
314
315     /*! \brief Getter for the number of local atoms.
316      *
317      *  \returns The number of local atoms.
318      */
319     int numAtomsLocal() const;
320
321     /*! \brief Getter for the total number of atoms.
322      *
323      *  \returns The total number of atoms.
324      */
325     int numAtomsAll() const;
326
327 private:
328     //! GPU PME stream.
329     const DeviceStream* pmeStream_;
330     //! GPU NBNXM local stream.
331     const DeviceStream* localStream_;
332     //! GPU NBNXM non-local stream.
333     const DeviceStream* nonLocalStream_;
334     //! GPU Update-constraints stream.
335     const DeviceStream* updateStream_;
336
337     // Streams to use for coordinates H2D and D2H copies (one event for each atom locality)
338     EnumerationArray<AtomLocality, const DeviceStream*> xCopyStreams_ = { { nullptr } };
339     // Streams to use for velocities H2D and D2H copies (one event for each atom locality)
340     EnumerationArray<AtomLocality, const DeviceStream*> vCopyStreams_ = { { nullptr } };
341     // Streams to use for forces H2D and D2H copies (one event for each atom locality)
342     EnumerationArray<AtomLocality, const DeviceStream*> fCopyStreams_ = { { nullptr } };
343     // Streams internal to this module
344     std::unique_ptr<DeviceStream> copyInStream_;
345     std::unique_ptr<DeviceStream> memsetStream_;
346
347     /*! \brief An array of events that indicate H2D copy is complete (one event for each atom locality)
348      *
349      * \todo Reconsider naming. It should be xCopiedToDevice or xH2DCopyComplete, etc.
350      */
351     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnDevice_;
352     //! A pointer to an event that the coordinates are ready after update-constraints execution
353     GpuEventSynchronizer* xUpdatedOnDeviceEvent_ = nullptr;
354     //! An array of events that indicate D2H copy of coordinates is complete (one event for each atom locality)
355     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnHost_;
356
357     //! An array of events that indicate D2H copy of velocities is complete (one event for each atom locality)
358     EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnHost_;
359
360     //! An array of events that indicate H2D copy of forces is complete (one event for each atom locality)
361     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnDevice_;
362     //! An array of events that indicate the forces were reduced on the GPU (one event for each atom locality)
363     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReducedOnDevice_;
364     //! An array of events that indicate D2H copy of forces is complete (one event for each atom locality)
365     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnHost_;
366
367     //! GPU context (for OpenCL builds)
368     const DeviceContext& deviceContext_;
369     //! Default GPU calls behavior
370     GpuApiCallBehavior transferKind_ = GpuApiCallBehavior::Async;
371     //! Required minimum divisor of the allocation size of the coordinates buffer
372     int allocationBlockSizeDivisor_ = 0;
373
374     //! Number of local atoms
375     int numAtomsLocal_ = -1;
376     //! Total number of atoms
377     int numAtomsAll_ = -1;
378
379     //! Device positions buffer
380     DeviceBuffer<RVec> d_x_;
381     //! Number of particles saved in the positions buffer
382     int d_xSize_ = -1;
383     //! Allocation size for the positions buffer
384     int d_xCapacity_ = -1;
385
386     //! Device velocities buffer
387     DeviceBuffer<RVec> d_v_;
388     //! Number of particles saved in the velocities buffer
389     int d_vSize_ = -1;
390     //! Allocation size for the velocities buffer
391     int d_vCapacity_ = -1;
392
393     //! Device force buffer
394     DeviceBuffer<RVec> d_f_;
395     //! Number of particles saved in the force buffer
396     int d_fSize_ = -1;
397     //! Allocation size for the force buffer
398     int d_fCapacity_ = -1;
399
400     //! \brief Pointer to wallcycle structure.
401     gmx_wallcycle* wcycle_;
402
403     /*! \brief Performs the copy of data from host to device buffer.
404      *
405      * \todo Template on locality.
406      *
407      *  \param[out] d_data         Device-side buffer.
408      *  \param[in]  h_data         Host-side buffer.
409      *  \param[in]  dataSize       Device-side data allocation size.
410      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
411      *  \param[in]  deviceStream   GPU stream to execute copy in.
412      */
413     void copyToDevice(DeviceBuffer<RVec>             d_data,
414                       gmx::ArrayRef<const gmx::RVec> h_data,
415                       int                            dataSize,
416                       AtomLocality                   atomLocality,
417                       const DeviceStream&            deviceStream);
418
419     /*! \brief Performs the copy of data from device to host buffer.
420      *
421      *  \param[out] h_data         Host-side buffer.
422      *  \param[in]  d_data         Device-side buffer.
423      *  \param[in]  dataSize       Device-side data allocation size.
424      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
425      *  \param[in]  deviceStream   GPU stream to execute copy in.
426      */
427     void copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
428                         DeviceBuffer<RVec>       d_data,
429                         int                      dataSize,
430                         AtomLocality             atomLocality,
431                         const DeviceStream&      deviceStream);
432
433     /*! \brief Performs the clearing of data in device buffer.
434      *
435      * \todo Template on locality.
436      *
437      *  \param[out] d_data         Device-side buffer.
438      *  \param[in]  dataSize       Device-side data allocation size.
439      *  \param[in]  atomLocality   If all, local or non-local ranges should be cleared.
440      *  \param[in]  deviceStream   GPU stream to execute copy in.
441      */
442     void clearOnDevice(DeviceBuffer<RVec>  d_data,
443                        int                 dataSize,
444                        AtomLocality        atomLocality,
445                        const DeviceStream& deviceStream) const;
446 };
447
448 } // namespace gmx
449
450 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H