6ba40d987e93154e06907dfd940969f57f46b4ae
[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 gpuContext 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 Make CommandStream visible in the CPU parts of the code so we
79          *       will not have to pass a void*.
80          * \todo Make a Context object visible in CPU parts of the code so we
81          *       will not have to pass a void*.
82          *
83          *  \param[in] commandStream  GPU stream, nullptr allowed.
84          *  \param[in] gpuContext     GPU context, nullptr allowed.
85          *  \param[in] transferKind   H2D/D2H transfer call behavior (synchronous or not).
86          *  \param[in] paddingSize    Padding size for coordinates buffer.
87          */
88         Impl(const void        *commandStream,
89              const void        *gpuContext,
90              GpuApiCallBehavior transferKind,
91              int                paddingSize);
92
93         ~Impl();
94
95
96         /*! \brief Set the ranges for local and non-local atoms and reallocates buffers.
97          *
98          * The coordinates buffer is reallocated with the padding added at the end. The
99          * size of padding is set by the constructor.
100          *
101          *  \param[in] numAtomsLocal  Number of atoms in local domain.
102          *  \param[in] numAtomsAll    Total number of atoms to handle.
103          */
104         void reinit(int numAtomsLocal, int numAtomsAll);
105
106         /*! \brief Returns the range of atoms to be copied based on the copy type (all, local or non-local).
107          *
108          * \todo There are at least three versions of the function with this functionality in the code:
109          *       this one and two more in NBNXM. These should be unified in a shape of a general function
110          *       in DD.
111          *
112          * \param[in]  atomLocality    If all, local or non-local ranges are needed.
113          *
114          * \returns Tuple, containing the index of the first atom in the range and the total number of atoms in the range.
115          */
116         std::tuple<int, int> getAtomRangesFromAtomLocality(AtomLocality  atomLocality);
117
118
119         /*! \brief Get the positions buffer on the GPU.
120          *
121          *  \returns GPU positions buffer.
122          */
123         DeviceBuffer<float> getCoordinates();
124
125         /*! \brief Copy positions to the GPU memory.
126          *
127          *  \param[in] h_x           Positions in the host memory.
128          *  \param[in] atomLocality  Locality of the particles to copy.
129          */
130         void copyCoordinatesToGpu(gmx::ArrayRef<const gmx::RVec>  h_x,
131                                   AtomLocality                    atomLocality);
132
133         /*! \brief Copy positions from the GPU memory.
134          *
135          *  \param[in] h_x           Positions buffer in the host memory.
136          *  \param[in] atomLocality  Locality of the particles to copy.
137          */
138         void copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec>  h_x,
139                                     AtomLocality              atomLocality);
140
141
142         /*! \brief Get the velocities buffer on the GPU.
143          *
144          *  \returns GPU velocities buffer.
145          */
146         DeviceBuffer<float> getVelocities();
147
148         /*! \brief Copy velocities to the GPU memory.
149          *
150          *  \param[in] h_v           Velocities in the host memory.
151          *  \param[in] atomLocality  Locality of the particles to copy.
152          */
153         void copyVelocitiesToGpu(gmx::ArrayRef<const gmx::RVec>  h_v,
154                                  AtomLocality                    atomLocality);
155
156         /*! \brief Copy velocities from the GPU memory.
157          *
158          *  \param[in] h_v           Velocities buffer in the host memory.
159          *  \param[in] atomLocality  Locality of the particles to copy.
160          */
161         void copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec>  h_v,
162                                    AtomLocality              atomLocality);
163
164
165         /*! \brief Get the force buffer on the GPU.
166          *
167          *  \returns GPU force buffer.
168          */
169         DeviceBuffer<float> getForces();
170
171         /*! \brief Copy forces to the GPU memory.
172          *
173          *  \param[in] h_f           Forces in the host memory.
174          *  \param[in] atomLocality  Locality of the particles to copy.
175          */
176         void copyForcesToGpu(gmx::ArrayRef<const gmx::RVec>  h_f,
177                              AtomLocality                    atomLocality);
178
179         /*! \brief Copy forces from the GPU memory.
180          *
181          *  \param[in] h_f           Forces buffer in the host memory.
182          *  \param[in] atomLocality  Locality of the particles to copy.
183          */
184         void copyForcesFromGpu(gmx::ArrayRef<gmx::RVec>  h_f,
185                                AtomLocality              atomLocality);
186
187         /*! \brief Synchronize the underlying GPU stream
188          */
189         void synchronizeStream();
190
191         /*! \brief Getter for the number of local atoms.
192          *
193          *  \returns The number of local atoms.
194          */
195         int numAtomsLocal();
196
197         /*! \brief Getter for the total number of atoms.
198          *
199          *  \returns The total number of atoms.
200          */
201         int numAtomsAll();
202
203     private:
204
205         /*! \brief GPU stream.
206          * \todo The stream should be set to non-nullptr once the synchronization points are restored
207          */
208         CommandStream        commandStream_              = nullptr;
209         /*! \brief GPU context (for OpenCL builds)
210          * \todo Make a Context class usable in CPU code
211          */
212         Context              gpuContext_                 = nullptr;
213         //! Default GPU calls behavior
214         GpuApiCallBehavior   transferKind_               = GpuApiCallBehavior::Async;
215         //! Padding size for the coordinates buffer
216         int                  paddingSize_                = 0;
217
218         //! Number of local atoms
219         int                  numAtomsLocal_              = -1;
220         //! Total number of atoms
221         int                  numAtomsAll_                = -1;
222
223         //! Device positions buffer
224         DeviceBuffer<float>  d_x_;
225         //! Number of particles saved in the positions buffer
226         int                  d_xSize_                    = -1;
227         //! Allocation size for the positions buffer
228         int                  d_xCapacity_                = -1;
229
230         //! Device velocities buffer
231         DeviceBuffer<float>  d_v_;
232         //! Number of particles saved in the velocities buffer
233         int                  d_vSize_                    = -1;
234         //! Allocation size for the velocities buffer
235         int                  d_vCapacity_                = -1;
236
237         //! Device force buffer
238         DeviceBuffer<float>  d_f_;
239         //! Number of particles saved in the force buffer
240         int                  d_fSize_                    = -1;
241         //! Allocation size for the force buffer
242         int                  d_fCapacity_                = -1;
243
244         /*! \brief Performs the copy of data from host to device buffer.
245          *
246          * \todo Template on locality.
247          *
248          * \param[in,out]  d_data        Device-side buffer.
249          * \param[in,out]  h_data        Host-side buffer.
250          * \param[in]      dataSize      Device-side data allocation size.
251          * \param[in]      atomLocality  If all, local or non-local ranges should be copied.
252          */
253         void copyToDevice(DeviceBuffer<float>            d_data,
254                           gmx::ArrayRef<const gmx::RVec> h_data,
255                           int                            dataSize,
256                           AtomLocality                   atomLocality);
257
258         /*! \brief Performs the copy of data from device to host buffer.
259          *
260          * \param[in,out]  h_data        Host-side buffer.
261          * \param[in,out]  d_data        Device-side buffer.
262          * \param[in]      dataSize      Device-side data allocation size.
263          * \param[in]      atomLocality  If all, local or non-local ranges should be copied.
264          */
265         void copyFromDevice(gmx::ArrayRef<gmx::RVec>  h_data,
266                             DeviceBuffer<float>       d_data,
267                             int                       dataSize,
268                             AtomLocality              atomLocality);
269 };
270
271 }      // namespace gmx
272
273 #endif // GMX_MDTYPES_STATE_PROPAGATOR_DATA_GPU_IMPL_H