ad07df2c132236c3763e2b7ce22d14f76506e454
[alexxy/gromacs.git] / src / gromacs / mdtypes / state_propagator_data_gpu.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 interfaces for GPU state data propagator object.
38  *
39  * This object stores and manages positions, velocities and forces for
40  * all particles in the system on the GPU.
41  *
42  * \todo Add cycle counters.
43  * \todo Add synchronization points.
44  *
45  * \author Artem Zhmurov <zhmurov@gmail.com>
46  *
47  * \inlibraryapi
48  * \ingroup module_mdtypes
49  */
50 #ifndef GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_H
51 #define GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_H
52
53 #include <memory>
54 #include <tuple>
55
56 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
57 #include "gromacs/gpu_utils/gpu_utils.h"
58 #include "gromacs/math/vectypes.h"
59 #include "gromacs/mdtypes/simulation_workload.h"
60 #include "gromacs/utility/arrayref.h"
61 #include "gromacs/utility/classhelpers.h"
62
63 #include "locality.h"
64
65 class DeviceContext;
66 class DeviceStream;
67 class GpuEventSynchronizer;
68 struct gmx_wallcycle;
69
70 namespace gmx
71 {
72 class DeviceStreamManager;
73
74 class StatePropagatorDataGpu
75 {
76 public:
77     /*! \brief Constructor
78      *
79      * The buffers are reallocated only at the reinit call, the padding is
80      * used there for the coordinates buffer. It is needed for PME and added at
81      * the end of the buffer. It is assumed that if the rank has PME duties on the
82      * GPU, all coordinates are copied to the GPU and hence, for this rank, the
83      * coordinates buffer is not split into local and non-local ranges. For other
84      * ranks, the padding size is zero. This works because only one rank ever does
85      * PME work on the GPU, and if that rank also does PP work that is the only
86      * rank. So all coordinates are always transferred.
87      *
88      * In OpenCL, only pmeStream is used since it is the only stream created in
89      * PME context. The local and non-local streams are only needed when buffer
90      * ops are offloaded. This feature is currently not available in OpenCL and
91      * hence these streams are not set in these builds.
92      *
93      *  \param[in] deviceStreamManager         Object that owns the DeviceContext and DeviceStreams.
94      *  \param[in] transferKind                H2D/D2H transfer call behavior (synchronous or not).
95      *  \param[in] allocationBlockSizeDivisor  Deterines padding size for coordinates buffer.
96      *  \param[in] wcycle                      Wall cycle counter data.
97      */
98     StatePropagatorDataGpu(const DeviceStreamManager& deviceStreamManager,
99                            GpuApiCallBehavior         transferKind,
100                            int                        allocationBlockSizeDivisor,
101                            gmx_wallcycle*             wcycle);
102
103     /*! \brief Constructor to use in PME-only rank and in tests.
104      *
105      *  This constructor should be used if only a coordinate device buffer should be managed
106      *  using a single stream. Any operation on force or velocity buffer as well as copy of
107      *  non-local coordinates will exit with assertion failure. Note, that the pmeStream can
108      *  not be a nullptr and the constructor will exit with an assertion failure.
109      *
110      *  \todo Currently, unsupported copy operations are blocked by assertion that the stream
111      *        not nullptr. This should be improved.
112      *
113      *  \param[in] pmeStream       Device PME stream, nullptr is not allowed.
114      *  \param[in] deviceContext   Device context, nullptr allowed for non-OpenCL builds.
115      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
116      *  \param[in] allocationBlockSizeDivisor Determines padding size for coordinates buffer.
117      *  \param[in] wcycle          Wall cycle counter data.
118      */
119     StatePropagatorDataGpu(const DeviceStream*  pmeStream,
120                            const DeviceContext& deviceContext,
121                            GpuApiCallBehavior   transferKind,
122                            int                  allocationBlockSizeDivisor,
123                            gmx_wallcycle*       wcycle);
124
125     //! Move constructor
126     StatePropagatorDataGpu(StatePropagatorDataGpu&& other) noexcept;
127     //! Move assignment
128     StatePropagatorDataGpu& operator=(StatePropagatorDataGpu&& other) noexcept;
129     //! Destructor
130     ~StatePropagatorDataGpu();
131
132     /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
133      *
134      * Reallocates coordinate, velocities and force buffers on the device.
135      *
136      * \note
137      * The coordinates buffer is (re)allocated, when required by PME, with a padding,
138      * the size of which is set by the constructor. The padding region clearing kernel
139      * is scheduled in the \p pmeStream_ (unlike the coordinates H2D) as only the PME
140      * task uses this padding area.
141      *
142      * \note
143      * The force buffer is cleared if its size increases, so that previously unused
144      * memory is cleared before forces are accumulated.
145      *
146      *  \param[in] numAtomsLocal  Number of atoms in local domain.
147      *  \param[in] numAtomsAll    Total number of atoms to handle.
148      */
149     void reinit(int numAtomsLocal, int numAtomsAll);
150
151     /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
152      *
153      * \todo There are at least three versions of the function with this functionality in the code:
154      *       this one and two more in NBNXM. These should be unified in a shape of a general function
155      *       in DD.
156      *
157      * \param[in]  atomLocality    If all, local or non-local ranges are needed.
158      *
159      * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
160      */
161     std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality) const;
162
163
164     /*! \brief Get the positions buffer on the GPU.
165      *
166      *  \returns GPU positions buffer.
167      */
168     DeviceBuffer<RVec> getCoordinates();
169
170     /*! \brief Copy positions to the GPU memory.
171      *
172      * Use \ref getCoordinatesReadyOnDeviceEvent to get the associated event synchronizer or
173      * \ref waitCoordinatesCopiedToDevice to wait for the copy completion.
174      * Note: the event is not marked in OpenCL, because it is not used.
175      *
176      *  \param[in] h_x           Positions in the host memory.
177      *  \param[in] atomLocality  Locality of the particles to copy.
178      */
179     void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
180
181     /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
182      *
183      * Returns the event synchronizer which indicates that the coordinates are ready for the
184      * consumption on the device. Takes into account that the producer may be different.
185      *
186      * If the update is offloaded, and the current step is not a DD/search step, the returned
187      * synchronizer indicates the completion of GPU update-constraint kernels. Otherwise, on search
188      * steps and if update is not offloaded, the coordinates are provided by the H2D copy and the
189      * returned synchronizer indicates that the copy is complete.
190      *
191      *  \param[in] atomLocality              Locality of the particles to wait for.
192      *  \param[in] simulationWork            The simulation lifetime flags.
193      *  \param[in] stepWork                  The step lifetime flags.
194      *  \param[in] gpuCoordinateHaloLaunched Event recorded when GPU coordinate halo has been launched.
195      *
196      *  \returns  The event to synchronize the stream that consumes coordinates on device.
197      */
198     GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
199                                                            const SimulationWorkload& simulationWork,
200                                                            const StepWorkload&       stepWork,
201                                                            GpuEventSynchronizer* gpuCoordinateHaloLaunched = nullptr);
202
203     /*! \brief Blocking wait until coordinates are copied to the device.
204      *
205      * Synchronizes the stream in which the copy was executed.
206      *
207      *  \param[in] atomLocality  Locality of the particles to wait for.
208      */
209     void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
210
211     /*! \brief Setter for the event synchronizer for the update is done on th GPU
212      *
213      *  \param[in] xUpdatedOnDeviceEvent  The event to synchronize the stream coordinates wre updated on device.
214      */
215     void setXUpdatedOnDeviceEvent(GpuEventSynchronizer* xUpdatedOnDeviceEvent);
216
217     /*! \brief Copy positions from the GPU memory, with an optional explicit dependency.
218      *
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.
222      */
223     void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x,
224                                 AtomLocality             atomLocality,
225                                 GpuEventSynchronizer*    dependency = nullptr);
226
227     /*! \brief Wait until coordinates are available on the host.
228      *
229      *  \param[in] atomLocality  Locality of the particles to wait for.
230      */
231     void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
232
233
234     /*! \brief Get the velocities buffer on the GPU.
235      *
236      *  \returns GPU velocities buffer.
237      */
238     DeviceBuffer<RVec> getVelocities();
239
240     /*! \brief Copy velocities to the GPU memory.
241      *
242      * Does not mark any event, because we don't use it anywhere at the moment.
243      *
244      *  \param[in] h_v           Velocities in the host memory.
245      *  \param[in] atomLocality  Locality of the particles to copy.
246      */
247     void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
248
249     /*! \brief Copy velocities from the GPU memory.
250      *
251      *  \param[in] h_v           Velocities buffer in the host memory.
252      *  \param[in] atomLocality  Locality of the particles to copy.
253      */
254     void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
255
256     /*! \brief Wait until velocities are available on the host.
257      *
258      *  \param[in] atomLocality  Locality of the particles to wait for.
259      */
260     void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
261
262
263     /*! \brief Get the force buffer on the GPU.
264      *
265      *  \returns GPU force buffer.
266      */
267     DeviceBuffer<RVec> getForces();
268
269     /*! \brief Copy forces to the GPU memory.
270      *
271      *  \param[in] h_f           Forces in the host memory.
272      *  \param[in] atomLocality  Locality of the particles to copy.
273      */
274     void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
275
276     /*! \brief Clear forces in the GPU memory.
277      *
278      *  \param[in] atomLocality  Locality of the particles to clear.
279      *  \param[in] dependency    Dependency event for this operation.
280      */
281     void clearForcesOnGpu(AtomLocality atomLocality, GpuEventSynchronizer* dependency);
282
283     /*! \brief Get the event synchronizer for the forces ready on device.
284      *
285      *  Returns either of the event synchronizers, depending on the offload scenario
286      *  for the current simulation timestep:
287      *  1. The forces are copied to the device (when GPU buffer ops are off)
288      *  2. The forces are reduced on the device (GPU buffer ops are on)
289      *
290      *  \param[in] stepWork        Step workload flags
291      *  \param[in] simulationWork  Simulation workload flags
292      *
293      *  \returns  The event to synchronize the stream that consumes forces on device.
294      */
295     GpuEventSynchronizer* getLocalForcesReadyOnDeviceEvent(StepWorkload       stepWork,
296                                                            SimulationWorkload simulationWork);
297
298     /*! \brief Getter for the event synchronizer for the forces are reduced on the GPU.
299      *
300      *  \param[in] atomLocality      Locality of the particles to wait for.
301      *  \returns                     The event to mark when forces are reduced on the GPU.
302      */
303     GpuEventSynchronizer* fReducedOnDevice(AtomLocality atomLocality);
304
305     /*! \brief Getter for the event synchronizer for the forces are ready on the GPU.
306      *
307      *  \param[in] atomLocality      Locality of the particles to wait for.
308      *  \returns                     The event to mark when forces are ready on the GPU.
309      */
310     GpuEventSynchronizer* fReadyOnDevice(AtomLocality atomLocality);
311
312     /*! \brief Copy forces from the GPU memory.
313      *
314      *  \param[in] h_f           Forces buffer in the host memory.
315      *  \param[in] atomLocality  Locality of the particles to copy.
316      */
317     void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
318
319     /*! \brief Wait until forces are available on the host.
320      *
321      *  \param[in] atomLocality  Locality of the particles to wait for.
322      */
323     void waitForcesReadyOnHost(AtomLocality atomLocality);
324
325     /*! \brief Getter for the update stream.
326      *
327      *  \todo This is temporary here, until the management of this stream is taken over.
328      *
329      *  \returns The device command stream to use in update-constraints.
330      */
331     const DeviceStream* getUpdateStream();
332
333     /*! \brief Getter for the number of local atoms.
334      *
335      *  \returns The number of local atoms.
336      */
337     int numAtomsLocal() const;
338
339     /*! \brief Getter for the total number of atoms.
340      *
341      *  \returns The total number of atoms.
342      */
343     int numAtomsAll() const;
344
345 private:
346     class Impl;
347     std::unique_ptr<Impl> impl_;
348     GMX_DISALLOW_COPY_AND_ASSIGN(StatePropagatorDataGpu);
349 };
350
351 } // namespace gmx
352
353 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_H