0b73395af8939a4e8d50bbe409be8aab77cb9d2d
[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 "gromacs/gpu_utils/devicebuffer.h"
49 #include "gromacs/math/vectypes.h"
50 #include "gromacs/mdtypes/state_propagator_data_gpu.h"
51 #include "gromacs/utility/classhelpers.h"
52
53 namespace gmx
54 {
55
56 class StatePropagatorDataGpu::Impl
57 {
58     public:
59
60         Impl();
61
62
63         /*! \brief Constructor
64          *
65          * The buffers are reallocated only at the reinit call, the padding is
66          * used there for the coordinates buffer. It is needed for PME and added at
67          * the end of the buffer. It is assumed that if the rank has PME duties on the
68          * GPU, all coordinates are copied to the GPU and hence, for this rank, the
69          * coordinates buffer is not split into local and non-local ranges. For other
70          * ranks, the padding size is zero. This works because only one rank ever does
71          * PME work on the GPU, and if that rank also does PP work that is the only
72          * rank. So all coordinates are always transferred.
73          *
74          * \note \p commandStream and \p deviceContext are allowed to be nullptr if
75          *       StatePropagatorDataGpu is not used in the OpenCL run (e.g. if PME
76          *       does not run on the GPU).
77          *
78          * \todo A CommandStream is now visible in the CPU parts of the code so we
79          *       can stop passing a void*.
80          * \todo A DeviceContext object is visible in CPU parts of the code so we
81          *       can stop passing a void*.
82          *
83          *  \param[in] pmeStream       Device PME stream, nullptr allowed.
84          *  \param[in] localStream     Device NBNXM local stream, nullptr allowed.
85          *  \param[in] nonLocalStream  Device NBNXM non-local stream, nullptr allowed.
86          *  \param[in] deviceContext   Device context, nullptr allowed.
87          *  \param[in] transferKind    H2D/D2H transfer call behavior (synchronous or not).
88          *  \param[in] paddingSize     Padding size for coordinates buffer.
89          */
90         Impl(const void        *pmeStream,
91              const void        *localStream,
92              const void        *nonLocalStream,
93              const void        *deviceContext,
94              GpuApiCallBehavior transferKind,
95              int                paddingSize);
96
97         ~Impl();
98
99
100         /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
101          *
102          * The coordinates buffer is reallocated with the padding added at the end. The
103          * size of padding is set by the constructor.
104          *
105          * \note The PME requires clearing of the padding, which is done in the pmeStream_.
106          *       Hence the pmeStream_ should be created in the gpuContext_.
107          *
108          *  \param[in] numAtomsLocal  Number of atoms in local domain.
109          *  \param[in] numAtomsAll    Total number of atoms to handle.
110          */
111         void reinit(int numAtomsLocal, int numAtomsAll);
112
113         /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
114          *
115          * \todo There are at least three versions of the function with this functionality in the code:
116          *       this one and two more in NBNXM. These should be unified in a shape of a general function
117          *       in DD.
118          *
119          * \param[in]  atomLocality    If all, local or non-local ranges are needed.
120          *
121          * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
122          */
123         std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality  atomLocality);
124
125
126         /*! \brief Get the positions buffer on the GPU.
127          *
128          *  \returns GPU positions buffer.
129          */
130         DeviceBuffer<float> getCoordinates();
131
132         /*! \brief Copy positions to the GPU memory.
133          *
134          *  \param[in] h_x           Positions in the host memory.
135          *  \param[in] atomLocality  Locality of the particles to copy.
136          */
137         void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec>  h_x,
138                                   AtomLocality                    atomLocality);
139
140         /*! \brief Copy positions from the GPU memory.
141          *
142          *  \param[in] h_x           Positions buffer in the host memory.
143          *  \param[in] atomLocality  Locality of the particles to copy.
144          */
145         void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec>  h_x,
146                                     AtomLocality              atomLocality);
147
148
149         /*! \brief Get the velocities buffer on the GPU.
150          *
151          *  \returns GPU velocities buffer.
152          */
153         DeviceBuffer<float> getVelocities();
154
155         /*! \brief Copy velocities to the GPU memory.
156          *
157          *  \param[in] h_v           Velocities in the host memory.
158          *  \param[in] atomLocality  Locality of the particles to copy.
159          */
160         void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec>  h_v,
161                                  AtomLocality                    atomLocality);
162
163         /*! \brief Copy velocities from the GPU memory.
164          *
165          *  \param[in] h_v           Velocities buffer in the host memory.
166          *  \param[in] atomLocality  Locality of the particles to copy.
167          */
168         void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec>  h_v,
169                                    AtomLocality              atomLocality);
170
171
172         /*! \brief Get the force buffer on the GPU.
173          *
174          *  \returns GPU force buffer.
175          */
176         DeviceBuffer<float> getForces();
177
178         /*! \brief Copy forces to the GPU memory.
179          *
180          *  \param[in] h_f           Forces in the host memory.
181          *  \param[in] atomLocality  Locality of the particles to copy.
182          */
183         void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec>  h_f,
184                              AtomLocality                    atomLocality);
185
186         /*! \brief Copy forces from the GPU memory.
187          *
188          *  \param[in] h_f           Forces buffer in the host memory.
189          *  \param[in] atomLocality  Locality of the particles to copy.
190          */
191         void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec>  h_f,
192                                AtomLocality              atomLocality);
193
194         /*! \brief Getter for the update stream.
195          *
196          *  \todo This is temporary here, until the management of this stream is taken over.
197          *
198          *  \returns The device command stream to use in update-constraints.
199          */
200         void* getUpdateStream();
201
202         /*! \brief Getter for the number of local atoms.
203          *
204          *  \returns The number of local atoms.
205          */
206         int numAtomsLocal();
207
208         /*! \brief Getter for the total number of atoms.
209          *
210          *  \returns The total number of atoms.
211          */
212         int numAtomsAll();
213
214     private:
215
216         //! GPU PME stream.
217         CommandStream        pmeStream_                  = nullptr;
218         //! GPU NBNXM local stream.
219         CommandStream        localStream_                = nullptr;
220         //! GPU NBNXM non-local stream
221         CommandStream        nonLocalStream_             = nullptr;
222         //! GPU Update-constreaints stream.
223         CommandStream        updateStream_               = nullptr;
224         /*! \brief GPU context (for OpenCL builds)
225          * \todo Make a Context class usable in CPU code
226          */
227         DeviceContext        deviceContext_              = nullptr;
228         //! Default GPU calls behavior
229         GpuApiCallBehavior   transferKind_               = GpuApiCallBehavior::Async;
230         //! Padding size for the coordinates buffer
231         int                  paddingSize_                = 0;
232
233         //! Number of local atoms
234         int                  numAtomsLocal_              = -1;
235         //! Total number of atoms
236         int                  numAtomsAll_                = -1;
237
238         //! Device positions buffer
239         DeviceBuffer<float>  d_x_;
240         //! Number of particles saved in the positions buffer
241         int                  d_xSize_                    = -1;
242         //! Allocation size for the positions buffer
243         int                  d_xCapacity_                = -1;
244
245         //! Device velocities buffer
246         DeviceBuffer<float>  d_v_;
247         //! Number of particles saved in the velocities buffer
248         int                  d_vSize_                    = -1;
249         //! Allocation size for the velocities buffer
250         int                  d_vCapacity_                = -1;
251
252         //! Device force buffer
253         DeviceBuffer<float>  d_f_;
254         //! Number of particles saved in the force buffer
255         int                  d_fSize_                    = -1;
256         //! Allocation size for the force buffer
257         int                  d_fCapacity_                = -1;
258
259         /*! \brief Performs the copy of data from host to device buffer.
260          *
261          * \todo Template on locality.
262          *
263          *  \param[out] d_data         Device-side buffer.
264          *  \param[in]  h_data         Host-side buffer.
265          *  \param[in]  dataSize       Device-side data allocation size.
266          *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
267          *  \param[in]  commandStream  GPU stream to execute copy in.
268          */
269         void copyToDevice(DeviceBuffer<float>                   d_data,
270                           gmx::ArrayRef<const gmx::RVec>        h_data,
271                           int                                   dataSize,
272                           AtomLocality                          atomLocality,
273                           CommandStream                         commandStream);
274
275         /*! \brief Performs the copy of data from device to host buffer.
276          *
277          *  \param[out] h_data         Host-side buffer.
278          *  \param[in]  d_data         Device-side buffer.
279          *  \param[in]  dataSize       Device-side data allocation size.
280          *  \param[in]  atomLocality   If all, local or non-local ranges should be copied.
281          *  \param[in]  commandStream  GPU stream to execute copy in.
282          */
283         void copyFromDevice(gmx::ArrayRef<gmx::RVec>  h_data,
284                             DeviceBuffer<float>       d_data,
285                             int                       dataSize,
286                             AtomLocality              atomLocality,
287                             CommandStream             commandStream);
288 };
289
290 }      // namespace gmx
291
292 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H