af073b6284dde35d568465192004b375be8b2dd0
[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, 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 == GMX_GPU_CUDA
52 #    include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
53 #elif GMX_GPU == GMX_GPU_OPENCL
54 #    include "gromacs/gpu_utils/gpueventsynchronizer_ocl.h"
55 #endif
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
58 #include "gromacs/utility/classhelpers.h"
59 #include "gromacs/utility/enumerationhelpers.h"
60
61 namespace gmx
62 {
63
64 class StatePropagatorDataGpu::Impl
65 {
66 public:
67     Impl();
68
69
70     /*! \brief Constructor
71      *
72      * The buffers are reallocated only at the reinit call, the padding is
73      * used there for the coordinates buffer. It is needed for PME and added at
74      * the end of the buffer. It is assumed that if the rank has PME duties on the
75      * GPU, all coordinates are copied to the GPU and hence, for this rank, the
76      * coordinates buffer is not split into local and non-local ranges. For other
77      * ranks, the padding size is zero. This works because only one rank ever does
78      * PME work on the GPU, and if that rank also does PP work that is the only
79      * rank. So all coordinates are always transferred.
80      *
81      * In OpenCL, only pmeStream is used since it is the only stream created in
82      * PME context. The local and non-local streams are only needed when buffer
83      * ops are offloaded. This feature is currently not available in OpenCL and
84      * hence these streams are not set in these builds.
85      *
86      * \note In CUDA, the update stream is created in the constructor as a temporary
87      *       solution, in place until the stream manager is introduced.
88      *       Note that this makes it impossible to construct this object in CUDA
89      *       builds executing on a host without any CUDA-capable device available.
90      *
91      * \note In CUDA, \p deviceContext is unused, hence always nullptr;
92      *       all stream arguments can also be nullptr in runs where the
93      *       respective streams are not required.
94      *       In OpenCL, \p deviceContext needs to be a valid device context.
95      *       In OpenCL runs StatePropagatorDataGpu is currently only used
96      *       with PME offload, and only on ranks with PME duty. Hence, the
97      *       \p pmeStream argument needs to be a valid OpenCL queue object
98      *       which must have been created in \p deviceContext.
99      *
100      * \todo Make a \p CommandStream visible in the CPU parts of the code so we
101      *       will not have to pass a void*.
102      * \todo Make a \p DeviceContext object visible in CPU parts of the code so we
103      *       will not have to pass a void*.
104      *
105      *  \param[in] pmeStream       Device PME stream, nullptr allowed.
106      *  \param[in] localStream     Device NBNXM local stream, nullptr allowed.
107      *  \param[in] nonLocalStream  Device NBNXM non-local stream, nullptr allowed.
108      *  \param[in] deviceContext   Device context, nullptr allowed.
109      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
110      *  \param[in] paddingSize     Padding size for coordinates buffer.
111      */
112     Impl(const void*        pmeStream,
113          const void*        localStream,
114          const void*        nonLocalStream,
115          const void*        deviceContext,
116          GpuApiCallBehavior transferKind,
117          int                paddingSize);
118
119     /*! \brief Constructor to use in PME-only rank and in tests.
120      *
121      *  This constructor should be used if only a coordinate device buffer should be managed
122      *  using a single stream. Any operation on force or velocity buffer as well as copy of
123      *  non-local coordinates will exit with assertion failure. Note, that the pmeStream can
124      *  not be a nullptr and the constructor will exit with an assertion failure.
125      *
126      *  \todo Currently, unsupported copy operations are blocked by assertion that the stream
127      *        not nullptr. This should be improved.
128      *
129      *  \param[in] pmeStream       Device PME stream, nullptr is not allowed.
130      *  \param[in] deviceContext   Device context, nullptr allowed for non-OpenCL builds.
131      *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
132      *  \param[in] paddingSize     Padding size for coordinates buffer.
133      */
134     Impl(const void* pmeStream, const void* deviceContext, GpuApiCallBehavior transferKind, int paddingSize);
135
136     ~Impl();
137
138
139     /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
140      *
141      * Reallocates coordinate, velocities and force buffers on the device.
142      *
143      * \note
144      * The coordinates buffer is (re)allocated, when required by PME, with a padding,
145      * the size of which is set by the constructor. The padding region clearing kernel
146      * is scheduled in the \p pmeStream_ (unlike the coordinates H2D) as only the PME
147      * task uses this padding area.
148      *
149      * \note
150      * The force buffer is cleared if its size increases, so that previously unused
151      * memory is cleared before forces are accumulated.
152      *
153      *  \param[in] numAtomsLocal  Number of atoms in local domain.
154      *  \param[in] numAtomsAll    Total number of atoms to handle.
155      */
156     void reinit(int numAtomsLocal, int numAtomsAll);
157
158     /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
159      *
160      * \todo There are at least three versions of the function with this functionality in the code:
161      *       this one and two more in NBNXM. These should be unified in a shape of a general function
162      *       in DD.
163      *
164      * \param[in]  atomLocality    If all, local or non-local ranges are needed.
165      *
166      * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
167      */
168     std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality atomLocality);
169
170
171     /*! \brief Get the positions buffer on the GPU.
172      *
173      *  \returns GPU positions buffer.
174      */
175     DeviceBuffer<float> getCoordinates();
176
177     /*! \brief Copy positions to the GPU memory.
178      *
179      *  \param[in] h_x           Positions in the host memory.
180      *  \param[in] atomLocality  Locality of the particles to copy.
181      */
182     void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec> h_x, AtomLocality atomLocality);
183
184     /*! \brief Get the event synchronizer of the coordinates ready for the consumption on the device.
185      *
186      * Returns the event synchronizer which indicates that the coordinates are ready for the
187      * consumption on the device. Takes into account that the producer may be different.
188      *
189      * If the update is offloaded, and the current step is not a DD/search step, the returned
190      * synchronizer indicates the completion of GPU update-constraint kernels. Otherwise, on search
191      * steps and if update is not offloaded, the coordinates are provided by the H2D copy and the
192      * returned synchronizer indicates that the copy is complete.
193      *
194      *  \param[in] atomLocality    Locality of the particles to wait for.
195      *  \param[in] simulationWork  The simulation lifetime flags.
196      *  \param[in] stepWork        The step lifetime flags.
197      *
198      *  \returns  The event to synchronize the stream that consumes coordinates on device.
199      */
200     GpuEventSynchronizer* getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
201                                                            const SimulationWorkload& simulationWork,
202                                                            const StepWorkload&       stepWork);
203
204     /*! \brief Blocking wait until coordinates are copied to the device.
205      *
206      * Synchronizes the stream in which the copy was executed.
207      *
208      *  \param[in] atomLocality  Locality of the particles to wait for.
209      */
210     void waitCoordinatesCopiedToDevice(AtomLocality atomLocality);
211
212     /*! \brief Getter for the event synchronizer for the update is done on th GPU
213      *
214      *  \returns  The event to synchronize the stream coordinates wre updated on device.
215      */
216     GpuEventSynchronizer* xUpdatedOnDevice();
217
218     /*! \brief Copy positions from the GPU memory.
219      *
220      *  \param[in] h_x           Positions buffer in the host memory.
221      *  \param[in] atomLocality  Locality of the particles to copy.
222      */
223     void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality);
224
225     /*! \brief Wait until coordinates are available on the host.
226      *
227      *  \param[in] atomLocality  Locality of the particles to wait for.
228      */
229     void waitCoordinatesReadyOnHost(AtomLocality atomLocality);
230
231
232     /*! \brief Get the velocities buffer on the GPU.
233      *
234      *  \returns GPU velocities buffer.
235      */
236     DeviceBuffer<float> getVelocities();
237
238     /*! \brief Copy velocities to the GPU memory.
239      *
240      *  \param[in] h_v           Velocities in the host memory.
241      *  \param[in] atomLocality  Locality of the particles to copy.
242      */
243     void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec> h_v, AtomLocality atomLocality);
244
245     /*! \brief Get the event synchronizer on the H2D velocities copy.
246      *
247      *  \param[in] atomLocality  Locality of the particles to wait for.
248      *
249      *  \returns  The event to synchronize the stream that consumes velocities on device.
250      */
251     GpuEventSynchronizer* getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality);
252
253     /*! \brief Copy velocities from the GPU memory.
254      *
255      *  \param[in] h_v           Velocities buffer in the host memory.
256      *  \param[in] atomLocality  Locality of the particles to copy.
257      */
258     void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality);
259
260     /*! \brief Wait until velocities are available on the host.
261      *
262      *  \param[in] atomLocality  Locality of the particles to wait for.
263      */
264     void waitVelocitiesReadyOnHost(AtomLocality atomLocality);
265
266
267     /*! \brief Get the force buffer on the GPU.
268      *
269      *  \returns GPU force buffer.
270      */
271     DeviceBuffer<float> getForces();
272
273     /*! \brief Copy forces to the GPU memory.
274      *
275      *  \param[in] h_f           Forces in the host memory.
276      *  \param[in] atomLocality  Locality of the particles to copy.
277      */
278     void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality);
279
280     /*! \brief Get the event synchronizer for the forces ready on device.
281      *
282      *  Returns either of the event synchronizers, depending on the offload scenario
283      *  for the current simulation timestep:
284      *  1. The forces are copied to the device (when GPU buffer ops are off)
285      *  2. The forces are reduced on the device (GPU buffer ops are on)
286      *
287      *  \todo Pass step workload instead of the useGpuFBufferOps boolean.
288      *
289      *  \param[in] atomLocality      Locality of the particles to wait for.
290      *  \param[in] useGpuFBufferOps  If the force buffer ops are offloaded to the GPU.
291      *
292      *  \returns  The event to synchronize the stream that consumes forces on device.
293      */
294     GpuEventSynchronizer* getForcesReadyOnDeviceEvent(AtomLocality atomLocality, bool useGpuFBufferOps);
295
296     /*! \brief Getter for the event synchronizer for the forces are reduced on the GPU.
297      *
298      *  \returns  The event to mark when forces are reduced on the GPU.
299      */
300     GpuEventSynchronizer* fReducedOnDevice();
301
302     /*! \brief Copy forces from the GPU memory.
303      *
304      *  \param[in] h_f           Forces buffer in the host memory.
305      *  \param[in] atomLocality  Locality of the particles to copy.
306      */
307     void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality);
308
309     /*! \brief Wait until forces are available on the host.
310      *
311      *  \param[in] atomLocality  Locality of the particles to wait for.
312      */
313     void waitForcesReadyOnHost(AtomLocality atomLocality);
314
315     /*! \brief Getter for the update stream.
316      *
317      *  \todo This is temporary here, until the management of this stream is taken over.
318      *
319      *  \returns The device command stream to use in update-constraints.
320      */
321     void* getUpdateStream();
322
323     /*! \brief Getter for the number of local atoms.
324      *
325      *  \returns The number of local atoms.
326      */
327     int numAtomsLocal();
328
329     /*! \brief Getter for the total number of atoms.
330      *
331      *  \returns The total number of atoms.
332      */
333     int numAtomsAll();
334
335 private:
336     //! GPU PME stream.
337     CommandStream pmeStream_ = nullptr;
338     //! GPU NBNXM local stream.
339     CommandStream localStream_ = nullptr;
340     //! GPU NBNXM non-local stream
341     CommandStream nonLocalStream_ = nullptr;
342     //! GPU Update-constreaints stream.
343     CommandStream updateStream_ = nullptr;
344
345     // Streams to use for coordinates H2D and D2H copies (one event for each atom locality)
346     EnumerationArray<AtomLocality, CommandStream> xCopyStreams_ = { { nullptr } };
347     // Streams to use for velocities H2D and D2H copies (one event for each atom locality)
348     EnumerationArray<AtomLocality, CommandStream> vCopyStreams_ = { { nullptr } };
349     // Streams to use for forces H2D and D2H copies (one event for each atom locality)
350     EnumerationArray<AtomLocality, CommandStream> fCopyStreams_ = { { nullptr } };
351
352     /*! \brief An array of events that indicate H2D copy is complete (one event for each atom locality)
353      *
354      * \todo Reconsider naming. It should be xCopiedToDevice or xH2DCopyComplete, etc.
355      */
356     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnDevice_;
357     //! An event that the coordinates are ready after update-constraints execution
358     GpuEventSynchronizer xUpdatedOnDevice_;
359     //! An array of events that indicate D2H copy of coordinates is complete (one event for each atom locality)
360     EnumerationArray<AtomLocality, GpuEventSynchronizer> xReadyOnHost_;
361
362     //! An array of events that indicate H2D copy of velocities is complete (one event for each atom locality)
363     EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnDevice_;
364     //! An array of events that indicate D2H copy of velocities is complete (one event for each atom locality)
365     EnumerationArray<AtomLocality, GpuEventSynchronizer> vReadyOnHost_;
366
367     //! An array of events that indicate H2D copy of forces is complete (one event for each atom locality)
368     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnDevice_;
369     //! An event that the forces were reduced on the GPU
370     GpuEventSynchronizer fReducedOnDevice_;
371     //! An array of events that indicate D2H copy of forces is complete (one event for each atom locality)
372     EnumerationArray<AtomLocality, GpuEventSynchronizer> fReadyOnHost_;
373
374     /*! \brief GPU context (for OpenCL builds)
375      * \todo Make a Context class usable in CPU code
376      */
377     DeviceContext deviceContext_ = nullptr;
378     //! Default GPU calls behavior
379     GpuApiCallBehavior transferKind_ = GpuApiCallBehavior::Async;
380     //! Padding size for the coordinates buffer
381     int paddingSize_ = 0;
382
383     //! Number of local atoms
384     int numAtomsLocal_ = -1;
385     //! Total number of atoms
386     int numAtomsAll_ = -1;
387
388     //! Device positions buffer
389     DeviceBuffer<float> d_x_;
390     //! Number of particles saved in the positions buffer
391     int d_xSize_ = -1;
392     //! Allocation size for the positions buffer
393     int d_xCapacity_ = -1;
394
395     //! Device velocities buffer
396     DeviceBuffer<float> d_v_;
397     //! Number of particles saved in the velocities buffer
398     int d_vSize_ = -1;
399     //! Allocation size for the velocities buffer
400     int d_vCapacity_ = -1;
401
402     //! Device force buffer
403     DeviceBuffer<float> d_f_;
404     //! Number of particles saved in the force buffer
405     int d_fSize_ = -1;
406     //! Allocation size for the force buffer
407     int d_fCapacity_ = -1;
408
409     /*! \brief Performs the copy of data from host to device buffer.
410      *
411      * \todo Template on locality.
412      *
413      *  \param[out] d_data         Device-side buffer.
414      *  \param[in]  h_data         Host-side buffer.
415      *  \param[in]  dataSize       Device-side data allocation size.
416      *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
417      *  \param[in]  commandStream  GPU stream to execute copy in.
418      */
419     void copyToDevice(DeviceBuffer<float>            d_data,
420                       gmx::ArrayRef<const gmx::RVec> h_data,
421                       int                            dataSize,
422                       AtomLocality                   atomLocality,
423                       CommandStream                  commandStream);
424
425     /*! \brief Performs the copy of data from device to host buffer.
426      *
427      *  \param[out] h_data         Host-side buffer.
428      *  \param[in]  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 copied.
431      *  \param[in]  commandStream  GPU stream to execute copy in.
432      */
433     void copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
434                         DeviceBuffer<float>      d_data,
435                         int                      dataSize,
436                         AtomLocality             atomLocality,
437                         CommandStream            commandStream);
438 };
439
440 } // namespace gmx
441
442 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H