cfceab3bff26cb44fe19ac97a9de4026b3c66529
[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 #if GMX_GPU_CUDA
54 #    include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
55 #elif GMX_GPU_OPENCL
56 #    include "gromacs/gpu_utils/gpueventsynchronizer_ocl.h"
57 #elif GMX_GPU_SYCL
58 #    include "gromacs/gpu_utils/gpueventsynchronizer_sycl.h"
59 #endif
60
61 #include "gromacs/math/vectypes.h"
62 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
63 #include "gromacs/utility/enumerationhelpers.h"
64
65 struct gmx_wallcycle;
66
67 namespace gmx
68 {
69
70 class StatePropagatorDataGpu::Impl
71 {
72 public:
73     Impl();
74
75
76     /*! \brief Constructor
77      *
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.
86      *
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.
91      *
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.
96      */
97     Impl(const DeviceStreamManager& deviceStreamManager,
98          GpuApiCallBehavior         transferKind,
99          int                        allocationBlockSizeDivisor,
100          gmx_wallcycle*             wcycle);
101
102     /*! \brief Constructor to use in PME-only rank and in tests.
103      *
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.
108      *
109      *  \todo Currently, unsupported copy operations are blocked by assertion that the stream
110      *        not nullptr. This should be improved.
111      *
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.
117      */
118     Impl(const DeviceStream*  pmeStream,
119          const DeviceContext& deviceContext,
120          GpuApiCallBehavior   transferKind,
121          int                  allocationBlockSizeDivisor,
122          gmx_wallcycle*       wcycle);
123
124     ~Impl();
125
126
127     /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
128      *
129      * Reallocates coordinate, velocities and force buffers on the device.
130      *
131      * \note
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.
136      *
137      * \note
138      * The force buffer is cleared if its size increases, so that previously unused
139      * memory is cleared before forces are accumulated.
140      *
141      *  \param[in] numAtomsLocal  Number of atoms in local domain.
142      *  \param[in] numAtomsAll    Total number of atoms to handle.
143      */
144     void reinit(int numAtomsLocal, int numAtomsAll);
145
146     /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
147      *
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
150      *       in DD.
151      *
152      * \param[in]  atomLocality    If all, local or non-local ranges are needed.
153      *
154      * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
155      */
156     std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality) const;
157
158
159     /*! \brief Get the positions buffer on the GPU.
160      *
161      *  \returns GPU positions buffer.
162      */
163     DeviceBuffer<RVec> getCoordinates();
164
165     /*! \brief Copy positions to the GPU memory.
166      *
167      *  \param[in] h_x           Positions in the host memory.
168      *  \param[in] atomLocality  Locality of the particles to copy.
169      */
170     void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
171
172     /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
173      *
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.
176      *
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.
181      *
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.
185      *
186      *  \returns  The event to synchronize the stream that consumes coordinates on device.
187      */
188     GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
189                                                            const SimulationWorkload& simulationWork,
190                                                            const StepWorkload&       stepWork);
191
192     /*! \brief Blocking wait until coordinates are copied to the device.
193      *
194      * Synchronizes the stream in which the copy was executed.
195      *
196      *  \param[in] atomLocality  Locality of the particles to wait for.
197      */
198     void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
199
200     /*! \brief Getter for the event synchronizer for the update is done on th GPU
201      *
202      *  \returns  The event to synchronize the stream coordinates wre updated on device.
203      */
204     GpuEventSynchronizer* xUpdatedOnDevice();
205
206     /*! \brief Copy positions from the GPU memory.
207      *
208      *  \param[in] h_x           Positions buffer in the host memory.
209      *  \param[in] atomLocality  Locality of the particles to copy.
210      */
211     void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality);
212
213     /*! \brief Wait until coordinates are available on the host.
214      *
215      *  \param[in] atomLocality  Locality of the particles to wait for.
216      */
217     void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
218
219
220     /*! \brief Get the velocities buffer on the GPU.
221      *
222      *  \returns GPU velocities buffer.
223      */
224     DeviceBuffer<RVec> getVelocities();
225
226     /*! \brief Copy velocities to the GPU memory.
227      *
228      *  \param[in] h_v           Velocities in the host memory.
229      *  \param[in] atomLocality  Locality of the particles to copy.
230      */
231     void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
232
233     /*! \brief Copy velocities from the GPU memory.
234      *
235      *  \param[in] h_v           Velocities buffer in the host memory.
236      *  \param[in] atomLocality  Locality of the particles to copy.
237      */
238     void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
239
240     /*! \brief Wait until velocities are available on the host.
241      *
242      *  \param[in] atomLocality  Locality of the particles to wait for.
243      */
244     void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
245
246
247     /*! \brief Get the force buffer on the GPU.
248      *
249      *  \returns GPU force buffer.
250      */
251     DeviceBuffer<RVec> getForces();
252
253     /*! \brief Copy forces to the GPU memory.
254      *
255      *  \param[in] h_f           Forces in the host memory.
256      *  \param[in] atomLocality  Locality of the particles to copy.
257      */
258     void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
259
260     /*! \brief Clear forces in the GPU memory.
261      *
262      *  \param[in] atomLocality  Locality of the particles to clear.
263      */
264     void clearForcesOnGpu(AtomLocality atomLocality);
265
266     /*! \brief Get the event synchronizer for the forces ready on device.
267      *
268      *  Returns either of the event synchronizers, depending on the offload scenario
269      *  for the current simulation timestep:
270      *  1. The forces are copied to the device (when GPU buffer ops are off)
271      *  2. The forces are reduced on the device (GPU buffer ops are on)
272      *
273      *  \todo Pass step workload instead of the useGpuFBufferOps boolean.
274      *
275      *  \param[in] atomLocality      Locality of the particles to wait for.
276      *  \param[in] useGpuFBufferOps  If the force buffer ops are offloaded to the GPU.
277      *
278      *  \returns  The event to synchronize the stream that consumes forces on device.
279      */
280     GpuEventSynchronizer* getForcesReadyOnDeviceEvent(AtomLocality atomLocality, bool useGpuFBufferOps);
281
282     /*! \brief Getter for the event synchronizer for the forces are reduced on the GPU.
283      *
284      *  \returns  The event to mark when forces are reduced on the GPU.
285      */
286     GpuEventSynchronizer* fReducedOnDevice();
287
288     /*! \brief Copy forces from the GPU memory.
289      *
290      *  \param[in] h_f           Forces buffer in the host memory.
291      *  \param[in] atomLocality  Locality of the particles to copy.
292      */
293     void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
294
295     /*! \brief Wait until forces are available on the host.
296      *
297      *  \param[in] atomLocality  Locality of the particles to wait for.
298      */
299     void waitForcesReadyOnHost(AtomLocality atomLocality);
300
301     /*! \brief Getter for the update stream.
302      *
303      *  \todo This is temporary here, until the management of this stream is taken over.
304      *
305      *  \returns The device command stream to use in update-constraints.
306      */
307     const DeviceStream* getUpdateStream();
308
309     /*! \brief Getter for the number of local atoms.
310      *
311      *  \returns The number of local atoms.
312      */
313     int numAtomsLocal() const;
314
315     /*! \brief Getter for the total number of atoms.
316      *
317      *  \returns The total number of atoms.
318      */
319     int numAtomsAll() const;
320
321 private:
322     //! GPU PME stream.
323     const DeviceStream* pmeStream_;
324     //! GPU NBNXM local stream.
325     const DeviceStream* localStream_;
326     //! GPU NBNXM non-local stream.
327     const DeviceStream* nonLocalStream_;
328     //! GPU Update-constreaints stream.
329     const DeviceStream* updateStream_;
330
331     // Streams to use for coordinates H2D and D2H copies (one event for each atom locality)
332     EnumerationArray<AtomLocality, const DeviceStream*> xCopyStreams_ = { { nullptr } };
333     // Streams to use for velocities H2D and D2H copies (one event for each atom locality)
334     EnumerationArray<AtomLocality, const DeviceStream*> vCopyStreams_ = { { nullptr } };
335     // Streams to use for forces H2D and D2H copies (one event for each atom locality)
336     EnumerationArray<AtomLocality, const DeviceStream*> fCopyStreams_ = { { nullptr } };
337
338     /*! \brief An array of events that indicate H2D copy is complete (one event for each atom locality)
339      *
340      * \todo Reconsider naming. It should be xCopiedToDevice or xH2DCopyComplete, etc.
341      */
342     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnDevice_;
343     //! An event that the coordinates are ready after update-constraints execution
344     GpuEventSynchronizer xUpdatedOnDevice_;
345     //! An array of events that indicate D2H copy of coordinates is complete (one event for each atom locality)
346     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnHost_;
347
348     //! An array of events that indicate D2H copy of velocities is complete (one event for each atom locality)
349     EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnHost_;
350
351     //! An array of events that indicate H2D copy of forces is complete (one event for each atom locality)
352     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnDevice_;
353     //! An event that the forces were reduced on the GPU
354     GpuEventSynchronizer fReducedOnDevice_;
355     //! An array of events that indicate D2H copy of forces is complete (one event for each atom locality)
356     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnHost_;
357
358     //! GPU context (for OpenCL builds)
359     const DeviceContext& deviceContext_;
360     //! Default GPU calls behavior
361     GpuApiCallBehavior transferKind_ = GpuApiCallBehavior::Async;
362     //! Required minimum divisor of the allocation size of the coordinates buffer
363     int allocationBlockSizeDivisor_ = 0;
364
365     //! Number of local atoms
366     int numAtomsLocal_ = -1;
367     //! Total number of atoms
368     int numAtomsAll_ = -1;
369
370     //! Device positions buffer
371     DeviceBuffer<RVec> d_x_;
372     //! Number of particles saved in the positions buffer
373     int d_xSize_ = -1;
374     //! Allocation size for the positions buffer
375     int d_xCapacity_ = -1;
376
377     //! Device velocities buffer
378     DeviceBuffer<RVec> d_v_;
379     //! Number of particles saved in the velocities buffer
380     int d_vSize_ = -1;
381     //! Allocation size for the velocities buffer
382     int d_vCapacity_ = -1;
383
384     //! Device force buffer
385     DeviceBuffer<RVec> d_f_;
386     //! Number of particles saved in the force buffer
387     int d_fSize_ = -1;
388     //! Allocation size for the force buffer
389     int d_fCapacity_ = -1;
390
391     //! \brief Pointer to wallcycle structure.
392     gmx_wallcycle* wcycle_;
393
394     /*! \brief Performs the copy of data from host to device buffer.
395      *
396      * \todo Template on locality.
397      *
398      *  \param[out] d_data         Device-side buffer.
399      *  \param[in]  h_data         Host-side buffer.
400      *  \param[in]  dataSize       Device-side data allocation size.
401      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
402      *  \param[in]  deviceStream   GPU stream to execute copy in.
403      */
404     void copyToDevice(DeviceBuffer<RVec>             d_data,
405                       gmx::ArrayRef<const gmx::RVec> h_data,
406                       int                            dataSize,
407                       AtomLocality                   atomLocality,
408                       const DeviceStream&            deviceStream);
409
410     /*! \brief Performs the copy of data from device to host buffer.
411      *
412      *  \param[out] h_data         Host-side buffer.
413      *  \param[in]  d_data         Device-side buffer.
414      *  \param[in]  dataSize       Device-side data allocation size.
415      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
416      *  \param[in]  deviceStream   GPU stream to execute copy in.
417      */
418     void copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
419                         DeviceBuffer<RVec>       d_data,
420                         int                      dataSize,
421                         AtomLocality             atomLocality,
422                         const DeviceStream&      deviceStream);
423
424     /*! \brief Performs the clearing of data in device buffer.
425      *
426      * \todo Template on locality.
427      *
428      *  \param[out] d_data         Device-side buffer.
429      *  \param[in]  dataSize       Device-side data allocation size.
430      *  \param[in]  atomLocality   If all, local or non-local ranges should be cleared.
431      *  \param[in]  deviceStream   GPU stream to execute copy in.
432      */
433     void clearOnDevice(DeviceBuffer<RVec>  d_data,
434                        int                 dataSize,
435                        AtomLocality        atomLocality,
436                        const DeviceStream& deviceStream) const;
437 };
438
439 } // namespace gmx
440
441 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H