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