Define stubs for GPU data structures in SYCL to make SYCL compilation possible
[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, 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 "gromacs/gpu_utils/devicebuffer.h"
51 #if GMX_GPU_CUDA
52 #    include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
53 #elif GMX_GPU_OPENCL
54 #    include "gromacs/gpu_utils/gpueventsynchronizer_ocl.h"
55 #elif GMX_GPU_SYCL
56 #    include "gromacs/gpu_utils/gpueventsynchronizer_sycl.h"
57 #endif
58 #include "gromacs/math/vectypes.h"
59 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
60 #include "gromacs/utility/classhelpers.h"
61 #include "gromacs/utility/enumerationhelpers.h"
62
63 struct gmx_wallcycle;
64
65 namespace gmx
66 {
67
68 class StatePropagatorDataGpu::Impl
69 {
70 public:
71     Impl();
72
73
74     /*! \brief Constructor
75      *
76      * The buffers are reallocated only at the reinit call, the padding is
77      * used there for the coordinates buffer. It is needed for PME and added at
78      * the end of the buffer. It is assumed that if the rank has PME duties on the
79      * GPU, all coordinates are copied to the GPU and hence, for this rank, the
80      * coordinates buffer is not split into local and non-local ranges. For other
81      * ranks, the padding size is zero. This works because only one rank ever does
82      * PME work on the GPU, and if that rank also does PP work that is the only
83      * rank. So all coordinates are always transferred.
84      *
85      * In OpenCL, only pmeStream is used since it is the only stream created in
86      * PME context. The local and non-local streams are only needed when buffer
87      * ops are offloaded. This feature is currently not available in OpenCL and
88      * hence these streams are not set in these builds.
89      *
90      *  \param[in] deviceStreamManager         Object that owns the DeviceContext and DeviceStreams.
91      *  \param[in] transferKind                H2D/D2H transfer call behavior (synchronous or not).
92      *  \param[in] allocationBlockSizeDivisor  Determines the padding size for coordinates buffer.
93      *  \param[in] wcycle                      Wall cycle counter data.
94      */
95     Impl(const DeviceStreamManager& deviceStreamManager,
96          GpuApiCallBehavior         transferKind,
97          int                        allocationBlockSizeDivisor,
98          gmx_wallcycle*             wcycle);
99
100     /*! \brief Constructor to use in PME-only rank and in tests.
101      *
102      *  This constructor should be used if only a coordinate device buffer should be managed
103      *  using a single stream. Any operation on force or velocity buffer as well as copy of
104      *  non-local coordinates will exit with assertion failure. Note, that the pmeStream can
105      *  not be a nullptr and the constructor will exit with an assertion failure.
106      *
107      *  \todo Currently, unsupported copy operations are blocked by assertion that the stream
108      *        not nullptr. This should be improved.
109      *
110      *  \param[in] pmeStream       Device PME stream, nullptr is not allowed.
111      *  \param[in] deviceContext   Device context, nullptr allowed for non-OpenCL builds.
112      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
113      *  \param[in] allocationBlockSizeDivisor  Determines the padding size for coordinates buffer.
114      *  \param[in] wcycle          Wall cycle counter data.
115      */
116     Impl(const DeviceStream*  pmeStream,
117          const DeviceContext& deviceContext,
118          GpuApiCallBehavior   transferKind,
119          int                  allocationBlockSizeDivisor,
120          gmx_wallcycle*       wcycle);
121
122     ~Impl();
123
124
125     /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
126      *
127      * Reallocates coordinate, velocities and force buffers on the device.
128      *
129      * \note
130      * The coordinates buffer is (re)allocated, when required by PME, with a padding,
131      * the size of which is set by the constructor. The padding region clearing kernel
132      * is scheduled in the \p pmeStream_ (unlike the coordinates H2D) as only the PME
133      * task uses this padding area.
134      *
135      * \note
136      * The force buffer is cleared if its size increases, so that previously unused
137      * memory is cleared before forces are accumulated.
138      *
139      *  \param[in] numAtomsLocal  Number of atoms in local domain.
140      *  \param[in] numAtomsAll    Total number of atoms to handle.
141      */
142     void reinit(int numAtomsLocal, int numAtomsAll);
143
144     /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
145      *
146      * \todo There are at least three versions of the function with this functionality in the code:
147      *       this one and two more in NBNXM. These should be unified in a shape of a general function
148      *       in DD.
149      *
150      * \param[in]  atomLocality    If all, local or non-local ranges are needed.
151      *
152      * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
153      */
154     std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality);
155
156
157     /*! \brief Get the positions buffer on the GPU.
158      *
159      *  \returns GPU positions buffer.
160      */
161     DeviceBuffer<RVec> getCoordinates();
162
163     /*! \brief Copy positions to the GPU memory.
164      *
165      *  \param[in] h_x           Positions in the host memory.
166      *  \param[in] atomLocality  Locality of the particles to copy.
167      */
168     void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
169
170     /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
171      *
172      * Returns the event synchronizer which indicates that the coordinates are ready for the
173      * consumption on the device. Takes into account that the producer may be different.
174      *
175      * If the update is offloaded, and the current step is not a DD/search step, the returned
176      * synchronizer indicates the completion of GPU update-constraint kernels. Otherwise, on search
177      * steps and if update is not offloaded, the coordinates are provided by the H2D copy and the
178      * returned synchronizer indicates that the copy is complete.
179      *
180      *  \param[in] atomLocality    Locality of the particles to wait for.
181      *  \param[in] simulationWork  The simulation lifetime flags.
182      *  \param[in] stepWork        The step lifetime flags.
183      *
184      *  \returns  The event to synchronize the stream that consumes coordinates on device.
185      */
186     GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
187                                                            const SimulationWorkload& simulationWork,
188                                                            const StepWorkload&       stepWork);
189
190     /*! \brief Blocking wait until coordinates are copied to the device.
191      *
192      * Synchronizes the stream in which the copy was executed.
193      *
194      *  \param[in] atomLocality  Locality of the particles to wait for.
195      */
196     void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
197
198     /*! \brief Getter for the event synchronizer for the update is done on th GPU
199      *
200      *  \returns  The event to synchronize the stream coordinates wre updated on device.
201      */
202     GpuEventSynchronizer* xUpdatedOnDevice();
203
204     /*! \brief Copy positions from the GPU memory.
205      *
206      *  \param[in] h_x           Positions buffer in the host memory.
207      *  \param[in] atomLocality  Locality of the particles to copy.
208      */
209     void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality);
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 Get the event synchronizer on the H2D velocities copy.
232      *
233      *  \param[in] atomLocality  Locality of the particles to wait for.
234      *
235      *  \returns  The event to synchronize the stream that consumes velocities on device.
236      */
237     GpuEventSynchronizer* getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality);
238
239     /*! \brief Copy velocities from the GPU memory.
240      *
241      *  \param[in] h_v           Velocities buffer in the host memory.
242      *  \param[in] atomLocality  Locality of the particles to copy.
243      */
244     void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
245
246     /*! \brief Wait until velocities are available on the host.
247      *
248      *  \param[in] atomLocality  Locality of the particles to wait for.
249      */
250     void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
251
252
253     /*! \brief Get the force buffer on the GPU.
254      *
255      *  \returns GPU force buffer.
256      */
257     DeviceBuffer<RVec> getForces();
258
259     /*! \brief Copy forces to the GPU memory.
260      *
261      *  \param[in] h_f           Forces in the host memory.
262      *  \param[in] atomLocality  Locality of the particles to copy.
263      */
264     void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, 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();
314
315     /*! \brief Getter for the total number of atoms.
316      *
317      *  \returns The total number of atoms.
318      */
319     int numAtomsAll();
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 H2D copy of velocities is complete (one event for each atom locality)
349     EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnDevice_;
350     //! An array of events that indicate D2H copy of velocities is complete (one event for each atom locality)
351     EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnHost_;
352
353     //! An array of events that indicate H2D copy of forces is complete (one event for each atom locality)
354     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnDevice_;
355     //! An event that the forces were reduced on the GPU
356     GpuEventSynchronizer fReducedOnDevice_;
357     //! An array of events that indicate D2H copy of forces is complete (one event for each atom locality)
358     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnHost_;
359
360     //! GPU context (for OpenCL builds)
361     const DeviceContext& deviceContext_;
362     //! Default GPU calls behavior
363     GpuApiCallBehavior transferKind_ = GpuApiCallBehavior::Async;
364     //! Required minimum divisor of the allocation size of the coordinates buffer
365     int allocationBlockSizeDivisor_ = 0;
366
367     //! Number of local atoms
368     int numAtomsLocal_ = -1;
369     //! Total number of atoms
370     int numAtomsAll_ = -1;
371
372     //! Device positions buffer
373     DeviceBuffer<RVec> d_x_;
374     //! Number of particles saved in the positions buffer
375     int d_xSize_ = -1;
376     //! Allocation size for the positions buffer
377     int d_xCapacity_ = -1;
378
379     //! Device velocities buffer
380     DeviceBuffer<RVec> d_v_;
381     //! Number of particles saved in the velocities buffer
382     int d_vSize_ = -1;
383     //! Allocation size for the velocities buffer
384     int d_vCapacity_ = -1;
385
386     //! Device force buffer
387     DeviceBuffer<RVec> d_f_;
388     //! Number of particles saved in the force buffer
389     int d_fSize_ = -1;
390     //! Allocation size for the force buffer
391     int d_fCapacity_ = -1;
392
393     //! \brief Pointer to wallcycle structure.
394     gmx_wallcycle* wcycle_;
395
396     /*! \brief Performs the copy of data from host to device buffer.
397      *
398      * \todo Template on locality.
399      *
400      *  \param[out] d_data         Device-side buffer.
401      *  \param[in]  h_data         Host-side buffer.
402      *  \param[in]  dataSize       Device-side data allocation size.
403      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
404      *  \param[in]  deviceStream   GPU stream to execute copy in.
405      */
406     void copyToDevice(DeviceBuffer<RVec>             d_data,
407                       gmx::ArrayRef<const gmx::RVec> h_data,
408                       int                            dataSize,
409                       AtomLocality                   atomLocality,
410                       const DeviceStream&            deviceStream);
411
412     /*! \brief Performs the copy of data from device to host buffer.
413      *
414      *  \param[out] h_data         Host-side buffer.
415      *  \param[in]  d_data         Device-side buffer.
416      *  \param[in]  dataSize       Device-side data allocation size.
417      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
418      *  \param[in]  deviceStream   GPU stream to execute copy in.
419      */
420     void copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
421                         DeviceBuffer<RVec>       d_data,
422                         int                      dataSize,
423                         AtomLocality             atomLocality,
424                         const DeviceStream&      deviceStream);
425 };
426
427 } // namespace gmx
428
429 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H