4353eecf4ecf278d25ef632b8230b8cd45b03715
[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      *
179      *  \returns  The event to synchronize the stream that consumes coordinates on device.
180      */
181     GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
182                                                            const SimulationWorkload& simulationWork,
183                                                            const StepWorkload&       stepWork);
184
185     /*! \brief Blocking wait until coordinates are copied to the device.
186      *
187      * Synchronizes the stream in which the copy was executed.
188      *
189      *  \param[in] atomLocality  Locality of the particles to wait for.
190      */
191     void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
192
193     /*! \brief Setter for the event synchronizer for the update is done on th GPU
194      *
195      *  \param[in] xUpdatedOnDeviceEvent  The event to synchronize the stream coordinates wre updated on device.
196      */
197     void setXUpdatedOnDeviceEvent(GpuEventSynchronizer* xUpdatedOnDeviceEvent);
198
199     /*! \brief Copy positions from the GPU memory.
200      *
201      *  \param[in] h_x           Positions buffer in the host memory.
202      *  \param[in] atomLocality  Locality of the particles to copy.
203      */
204     void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality);
205
206     /*! \brief Wait until coordinates are available on the host.
207      *
208      *  \param[in] atomLocality  Locality of the particles to wait for.
209      */
210     void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
211
212
213     /*! \brief Get the velocities buffer on the GPU.
214      *
215      *  \returns GPU velocities buffer.
216      */
217     DeviceBuffer<RVec> getVelocities();
218
219     /*! \brief Copy velocities to the GPU memory.
220      *
221      *  \param[in] h_v           Velocities in the host memory.
222      *  \param[in] atomLocality  Locality of the particles to copy.
223      */
224     void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
225
226     /*! \brief Copy velocities from the GPU memory.
227      *
228      *  \param[in] h_v           Velocities buffer in the host memory.
229      *  \param[in] atomLocality  Locality of the particles to copy.
230      */
231     void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
232
233     /*! \brief Wait until velocities are available on the host.
234      *
235      *  \param[in] atomLocality  Locality of the particles to wait for.
236      */
237     void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
238
239
240     /*! \brief Get the force buffer on the GPU.
241      *
242      *  \returns GPU force buffer.
243      */
244     DeviceBuffer<RVec> getForces();
245
246     /*! \brief Copy forces to the GPU memory.
247      *
248      *  \param[in] h_f           Forces in the host memory.
249      *  \param[in] atomLocality  Locality of the particles to copy.
250      */
251     void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
252
253     /*! \brief Clear forces in the GPU memory.
254      *
255      *  \param[in] atomLocality  Locality of the particles to clear.
256      */
257     void clearForcesOnGpu(AtomLocality atomLocality);
258
259     /*! \brief Get the event synchronizer for the forces ready on device.
260      *
261      *  Returns either of the event synchronizers, depending on the offload scenario
262      *  for the current simulation timestep:
263      *  1. The forces are copied to the device (when GPU buffer ops are off)
264      *  2. The forces are reduced on the device (GPU buffer ops are on)
265      *
266      *  \todo Pass step workload instead of the useGpuFBufferOps boolean.
267      *
268      *  \param[in] atomLocality      Locality of the particles to wait for.
269      *  \param[in] useGpuFBufferOps  If the force buffer ops are offloaded to the GPU.
270      *
271      *  \returns  The event to synchronize the stream that consumes forces on device.
272      */
273     GpuEventSynchronizer* getForcesReadyOnDeviceEvent(AtomLocality atomLocality, bool useGpuFBufferOps);
274
275     /*! \brief Getter for the event synchronizer for the forces are reduced on the GPU.
276      *
277      *  \returns  The event to mark when forces are reduced on the GPU.
278      */
279     GpuEventSynchronizer* fReducedOnDevice();
280
281     /*! \brief Copy forces from the GPU memory.
282      *
283      *  \param[in] h_f           Forces buffer in the host memory.
284      *  \param[in] atomLocality  Locality of the particles to copy.
285      */
286     void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
287
288     /*! \brief Wait until forces are available on the host.
289      *
290      *  \param[in] atomLocality  Locality of the particles to wait for.
291      */
292     void waitForcesReadyOnHost(AtomLocality atomLocality);
293
294     /*! \brief Getter for the update stream.
295      *
296      *  \todo This is temporary here, until the management of this stream is taken over.
297      *
298      *  \returns The device command stream to use in update-constraints.
299      */
300     const DeviceStream* getUpdateStream();
301
302     /*! \brief Getter for the number of local atoms.
303      *
304      *  \returns The number of local atoms.
305      */
306     int numAtomsLocal() const;
307
308     /*! \brief Getter for the total number of atoms.
309      *
310      *  \returns The total number of atoms.
311      */
312     int numAtomsAll() const;
313
314 private:
315     //! GPU PME stream.
316     const DeviceStream* pmeStream_;
317     //! GPU NBNXM local stream.
318     const DeviceStream* localStream_;
319     //! GPU NBNXM non-local stream.
320     const DeviceStream* nonLocalStream_;
321     //! GPU Update-constraints stream.
322     const DeviceStream* updateStream_;
323
324     // Streams to use for coordinates H2D and D2H copies (one event for each atom locality)
325     EnumerationArray<AtomLocality, const DeviceStream*> xCopyStreams_ = { { nullptr } };
326     // Streams to use for velocities H2D and D2H copies (one event for each atom locality)
327     EnumerationArray<AtomLocality, const DeviceStream*> vCopyStreams_ = { { nullptr } };
328     // Streams to use for forces H2D and D2H copies (one event for each atom locality)
329     EnumerationArray<AtomLocality, const DeviceStream*> fCopyStreams_ = { { nullptr } };
330
331     /*! \brief An array of events that indicate H2D copy is complete (one event for each atom locality)
332      *
333      * \todo Reconsider naming. It should be xCopiedToDevice or xH2DCopyComplete, etc.
334      */
335     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnDevice_;
336     //! A pointer to an event that the coordinates are ready after update-constraints execution
337     GpuEventSynchronizer* xUpdatedOnDeviceEvent_ = nullptr;
338     //! An array of events that indicate D2H copy of coordinates is complete (one event for each atom locality)
339     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnHost_;
340
341     //! An array of events that indicate D2H copy of velocities is complete (one event for each atom locality)
342     EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnHost_;
343
344     //! An array of events that indicate H2D copy of forces is complete (one event for each atom locality)
345     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnDevice_;
346     //! An event that the forces were reduced on the GPU
347     GpuEventSynchronizer fReducedOnDevice_;
348     //! An array of events that indicate D2H copy of forces is complete (one event for each atom locality)
349     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnHost_;
350
351     //! GPU context (for OpenCL builds)
352     const DeviceContext& deviceContext_;
353     //! Default GPU calls behavior
354     GpuApiCallBehavior transferKind_ = GpuApiCallBehavior::Async;
355     //! Required minimum divisor of the allocation size of the coordinates buffer
356     int allocationBlockSizeDivisor_ = 0;
357
358     //! Number of local atoms
359     int numAtomsLocal_ = -1;
360     //! Total number of atoms
361     int numAtomsAll_ = -1;
362
363     //! Device positions buffer
364     DeviceBuffer<RVec> d_x_;
365     //! Number of particles saved in the positions buffer
366     int d_xSize_ = -1;
367     //! Allocation size for the positions buffer
368     int d_xCapacity_ = -1;
369
370     //! Device velocities buffer
371     DeviceBuffer<RVec> d_v_;
372     //! Number of particles saved in the velocities buffer
373     int d_vSize_ = -1;
374     //! Allocation size for the velocities buffer
375     int d_vCapacity_ = -1;
376
377     //! Device force buffer
378     DeviceBuffer<RVec> d_f_;
379     //! Number of particles saved in the force buffer
380     int d_fSize_ = -1;
381     //! Allocation size for the force buffer
382     int d_fCapacity_ = -1;
383
384     //! \brief Pointer to wallcycle structure.
385     gmx_wallcycle* wcycle_;
386
387     /*! \brief Performs the copy of data from host to device buffer.
388      *
389      * \todo Template on locality.
390      *
391      *  \param[out] d_data         Device-side buffer.
392      *  \param[in]  h_data         Host-side buffer.
393      *  \param[in]  dataSize       Device-side data allocation size.
394      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
395      *  \param[in]  deviceStream   GPU stream to execute copy in.
396      */
397     void copyToDevice(DeviceBuffer<RVec>             d_data,
398                       gmx::ArrayRef<const gmx::RVec> h_data,
399                       int                            dataSize,
400                       AtomLocality                   atomLocality,
401                       const DeviceStream&            deviceStream);
402
403     /*! \brief Performs the copy of data from device to host buffer.
404      *
405      *  \param[out] h_data         Host-side buffer.
406      *  \param[in]  d_data         Device-side buffer.
407      *  \param[in]  dataSize       Device-side data allocation size.
408      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
409      *  \param[in]  deviceStream   GPU stream to execute copy in.
410      */
411     void copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
412                         DeviceBuffer<RVec>       d_data,
413                         int                      dataSize,
414                         AtomLocality             atomLocality,
415                         const DeviceStream&      deviceStream);
416
417     /*! \brief Performs the clearing of data in device buffer.
418      *
419      * \todo Template on locality.
420      *
421      *  \param[out] d_data         Device-side buffer.
422      *  \param[in]  dataSize       Device-side data allocation size.
423      *  \param[in]  atomLocality   If all, local or non-local ranges should be cleared.
424      *  \param[in]  deviceStream   GPU stream to execute copy in.
425      */
426     void clearOnDevice(DeviceBuffer<RVec>  d_data,
427                        int                 dataSize,
428                        AtomLocality        atomLocality,
429                        const DeviceStream& deviceStream) const;
430 };
431
432 } // namespace gmx
433
434 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H