Pass the GPU streams to StatePropagatorDataGpu constructor
[alexxy/gromacs.git] / src / gromacs / mdtypes / state_propagator_data_gpu_impl_gpu.cpp
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 Definitions of interfaces for GPU state data propagator object.
38  *
39  * \author Artem Zhmurov <zhmurov@gmail.com>
40  *
41  * \ingroup module_mdtypes
42  */
43 #include "gmxpre.h"
44
45 #include "config.h"
46
47 #if GMX_GPU != GMX_GPU_NONE
48
49 #if GMX_GPU == GMX_GPU_CUDA
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #endif
52 #include "gromacs/gpu_utils/devicebuffer.h"
53 #if GMX_GPU == GMX_GPU_OPENCL
54 #include "gromacs/gpu_utils/oclutils.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
60 #include "state_propagator_data_gpu_impl.h"
61
62 namespace gmx
63 {
64
65 StatePropagatorDataGpu::Impl::Impl(const void            *pmeStream,
66                                    const void            *localStream,
67                                    const void            *nonLocalStream,
68                                    const void            *deviceContext,
69                                    GpuApiCallBehavior     transferKind,
70                                    int                    paddingSize) :
71     transferKind_(transferKind),
72     paddingSize_(paddingSize)
73 {
74
75     GMX_RELEASE_ASSERT(getenv("GMX_USE_GPU_BUFFER_OPS") == nullptr, "GPU buffer ops are not supported in this build.");
76
77     if (pmeStream != nullptr)
78     {
79         pmeStream_ = *static_cast<const CommandStream*>(pmeStream);
80     }
81     if (localStream != nullptr)
82     {
83         localStream_ = *static_cast<const CommandStream*>(localStream);
84     }
85     if (nonLocalStream != nullptr)
86     {
87         nonLocalStream_ = *static_cast<const CommandStream*>(nonLocalStream);
88     }
89 // The OpenCL build will never use the updateStream
90 // TODO: The update stream should be created only when it is needed.
91 #if GMX_GPU == GMX_GPU_CUDA
92     cudaError_t stat;
93     stat = cudaStreamCreate(&updateStream_);
94     CU_RET_ERR(stat, "CUDA stream creation failed in StatePropagatorDataGpu");
95 #endif
96     if (deviceContext != nullptr)
97     {
98         deviceContext_ = *static_cast<const DeviceContext*>(deviceContext);
99     }
100
101 }
102
103 StatePropagatorDataGpu::Impl::~Impl()
104 {
105 }
106
107 void StatePropagatorDataGpu::Impl::reinit(int numAtomsLocal, int numAtomsAll)
108 {
109 #if GMX_GPU == GMX_GPU_OPENCL
110     GMX_ASSERT(deviceContext_ != nullptr, "GPU context should be set in OpenCL builds.");
111 #endif
112     numAtomsLocal_ = numAtomsLocal;
113     numAtomsAll_   = numAtomsAll;
114
115     int numAtomsPadded;
116     if (paddingSize_ > 0)
117     {
118         numAtomsPadded = ((numAtomsAll_ + paddingSize_ - 1 ) / paddingSize_ )*paddingSize_;
119     }
120     else
121     {
122         numAtomsPadded = numAtomsAll_;
123     }
124
125     reallocateDeviceBuffer(&d_x_, DIM*numAtomsPadded, &d_xSize_, &d_xCapacity_, deviceContext_);
126
127     const size_t paddingAllocationSize = numAtomsPadded - numAtomsAll_;
128     if (paddingAllocationSize > 0)
129     {
130         // The PME stream is used here because:
131         // 1. The padding clearing is only needed by PME.
132         // 2. It is the stream that is created in the PME OpenCL context.
133         clearDeviceBufferAsync(&d_x_, DIM*numAtomsAll_, DIM*paddingAllocationSize, pmeStream_);
134     }
135
136     reallocateDeviceBuffer(&d_v_, DIM*numAtomsAll_, &d_vSize_, &d_vCapacity_, deviceContext_);
137     reallocateDeviceBuffer(&d_f_, DIM*numAtomsAll_, &d_fSize_, &d_fCapacity_, deviceContext_);
138
139 }
140
141 std::tuple<int, int> StatePropagatorDataGpu::Impl::getAtomRangesFromAtomLocality(AtomLocality  atomLocality)
142 {
143     int atomsStartAt   = 0;
144     int numAtomsToCopy = 0;
145     switch (atomLocality)
146     {
147         case AtomLocality::All:
148             atomsStartAt    = 0;
149             numAtomsToCopy  = numAtomsAll_;
150             break;
151         case AtomLocality::Local:
152             atomsStartAt    = 0;
153             numAtomsToCopy  = numAtomsLocal_;
154             break;
155         case AtomLocality::NonLocal:
156             atomsStartAt    = numAtomsLocal_;
157             numAtomsToCopy  = numAtomsAll_ - numAtomsLocal_;
158             break;
159         default:
160             GMX_RELEASE_ASSERT(false, "Wrong range of atoms requested in GPU state data manager. Should be All, Local or NonLocal.");
161     }
162     GMX_ASSERT(atomsStartAt   >= 0, "The first elemtnt to copy has negative index. Probably, the GPU propagator state was not initialized.");
163     GMX_ASSERT(numAtomsToCopy >= 0, "Number of atoms to copy is negative. Probably, the GPU propagator state was not initialized.");
164     return std::make_tuple(atomsStartAt, numAtomsToCopy);
165 }
166
167 void StatePropagatorDataGpu::Impl::copyToDevice(DeviceBuffer<float>                   d_data,
168                                                 const gmx::ArrayRef<const gmx::RVec>  h_data,
169                                                 int                                   dataSize,
170                                                 AtomLocality                          atomLocality,
171                                                 CommandStream                         commandStream)
172 {
173
174 #if GMX_GPU == GMX_GPU_OPENCL
175     GMX_ASSERT(deviceContext_ != nullptr, "GPU context should be set in OpenCL builds.");
176     // The PME stream is used for OpenCL builds, because it is the context that it associated with the
177     // PME task which requires the coordinates managed here in OpenCL.
178     // TODO: This will have to be changed when the OpenCL implementation will be extended.
179     commandStream = pmeStream_;
180 #endif
181
182     GMX_UNUSED_VALUE(dataSize);
183
184     GMX_ASSERT(dataSize >= 0, "Trying to copy to device buffer before it was allocated.");
185
186     int atomsStartAt, numAtomsToCopy;
187     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
188
189     int elementsStartAt   = atomsStartAt*DIM;
190     int numElementsToCopy = numAtomsToCopy*DIM;
191
192     if (numAtomsToCopy != 0)
193     {
194         GMX_ASSERT(elementsStartAt + numElementsToCopy <= dataSize, "The device allocation is smaller than requested copy range.");
195         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(), "The host buffer is smaller than the requested copy range.");
196
197         copyToDeviceBuffer(&d_data, reinterpret_cast<const float *>(&h_data.data()[atomsStartAt]),
198                            elementsStartAt, numElementsToCopy,
199                            commandStream, transferKind_, nullptr);
200     }
201 }
202
203 void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec>  h_data,
204                                                   DeviceBuffer<float>       d_data,
205                                                   int                       dataSize,
206                                                   AtomLocality              atomLocality,
207                                                   CommandStream             commandStream)
208 {
209
210 #if GMX_GPU == GMX_GPU_OPENCL
211     GMX_ASSERT(deviceContext_ != nullptr, "GPU context should be set in OpenCL builds.");
212     // The PME stream is used for OpenCL builds, because it is the context that it associated with the
213     // PME task which requires the coordinates managed here in OpenCL.
214     // TODO: This will have to be changed when the OpenCL implementation will be extended.
215     commandStream = pmeStream_;
216 #endif
217
218     GMX_UNUSED_VALUE(dataSize);
219
220     GMX_ASSERT(dataSize >= 0, "Trying to copy from device buffer before it was allocated.");
221
222     int atomsStartAt, numAtomsToCopy;
223     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
224
225     int elementsStartAt   = atomsStartAt*DIM;
226     int numElementsToCopy = numAtomsToCopy*DIM;
227
228     if (numAtomsToCopy != 0)
229     {
230         GMX_ASSERT(elementsStartAt + numElementsToCopy <= dataSize, "The device allocation is smaller than requested copy range.");
231         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(), "The host buffer is smaller than the requested copy range.");
232
233         copyFromDeviceBuffer(reinterpret_cast<float*>(&h_data.data()[atomsStartAt]), &d_data,
234                              elementsStartAt, numElementsToCopy,
235                              commandStream, transferKind_, nullptr);
236     }
237 }
238
239 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getCoordinates()
240 {
241     return d_x_;
242 }
243
244 void StatePropagatorDataGpu::Impl::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_x,
245                                                         AtomLocality                          atomLocality)
246 {
247     // TODO: Use the correct stream
248     copyToDevice(d_x_, h_x, d_xSize_, atomLocality, nullptr);
249 }
250
251 void StatePropagatorDataGpu::Impl::copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec>  h_x,
252                                                           AtomLocality              atomLocality)
253 {
254     // TODO: Use the correct stream
255     copyFromDevice(h_x, d_x_, d_xSize_, atomLocality, nullptr);
256 }
257
258
259 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getVelocities()
260 {
261     return d_v_;
262 }
263
264 void StatePropagatorDataGpu::Impl::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_v,
265                                                        AtomLocality                          atomLocality)
266 {
267     // TODO: Use the correct stream
268     copyToDevice(d_v_, h_v, d_vSize_, atomLocality, nullptr);
269 }
270
271 void StatePropagatorDataGpu::Impl::copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec>  h_v,
272                                                          AtomLocality              atomLocality)
273 {
274     // TODO: Use the correct stream
275     copyFromDevice(h_v, d_v_, d_vSize_, atomLocality, nullptr);
276 }
277
278
279 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getForces()
280 {
281     return d_f_;
282 }
283
284 void StatePropagatorDataGpu::Impl::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_f,
285                                                    AtomLocality                          atomLocality)
286 {
287     // TODO: Use the correct stream
288     copyToDevice(d_f_, h_f, d_fSize_, atomLocality, nullptr);
289 }
290
291 void StatePropagatorDataGpu::Impl::copyForcesFromGpu(gmx::ArrayRef<gmx::RVec>  h_f,
292                                                      AtomLocality              atomLocality)
293 {
294     // TODO: Use the correct stream
295     copyFromDevice(h_f, d_f_, d_fSize_, atomLocality, nullptr);
296 }
297
298 void* StatePropagatorDataGpu::Impl::getUpdateStream()
299 {
300     return updateStream_;
301 }
302
303 int StatePropagatorDataGpu::Impl::numAtomsLocal()
304 {
305     return numAtomsLocal_;
306 }
307
308 int StatePropagatorDataGpu::Impl::numAtomsAll()
309 {
310     return numAtomsAll_;
311 }
312
313
314
315 StatePropagatorDataGpu::StatePropagatorDataGpu(const void        *pmeStream,
316                                                const void        *localStream,
317                                                const void        *nonLocalStream,
318                                                const void        *deviceContext,
319                                                GpuApiCallBehavior transferKind,
320                                                int                paddingSize)
321     : impl_(new Impl(pmeStream,
322                      localStream,
323                      nonLocalStream,
324                      deviceContext,
325                      transferKind,
326                      paddingSize))
327 {
328 }
329
330 StatePropagatorDataGpu::StatePropagatorDataGpu(StatePropagatorDataGpu && /* other */) noexcept = default;
331
332 StatePropagatorDataGpu &StatePropagatorDataGpu::operator=(StatePropagatorDataGpu && /* other */) noexcept = default;
333
334 StatePropagatorDataGpu::~StatePropagatorDataGpu() = default;
335
336
337 void StatePropagatorDataGpu::reinit(int numAtomsLocal, int numAtomsAll)
338 {
339     return impl_->reinit(numAtomsLocal, numAtomsAll);
340 }
341
342 std::tuple<int, int> StatePropagatorDataGpu::getAtomRangesFromAtomLocality(AtomLocality  atomLocality)
343 {
344     return impl_->getAtomRangesFromAtomLocality(atomLocality);
345 }
346
347
348 DeviceBuffer<float> StatePropagatorDataGpu::getCoordinates()
349 {
350     return impl_->getCoordinates();
351 }
352
353 void StatePropagatorDataGpu::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_x,
354                                                   AtomLocality                          atomLocality)
355 {
356     return impl_->copyCoordinatesToGpu(h_x, atomLocality);
357 }
358
359 void StatePropagatorDataGpu::copyCoordinatesFromGpu(gmx::ArrayRef<RVec>  h_x,
360                                                     AtomLocality         atomLocality)
361 {
362     return impl_->copyCoordinatesFromGpu(h_x, atomLocality);
363 }
364
365
366 DeviceBuffer<float> StatePropagatorDataGpu::getVelocities()
367 {
368     return impl_->getVelocities();
369 }
370
371 void StatePropagatorDataGpu::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_v,
372                                                  AtomLocality                          atomLocality)
373 {
374     return impl_->copyVelocitiesToGpu(h_v, atomLocality);
375 }
376
377 void StatePropagatorDataGpu::copyVelocitiesFromGpu(gmx::ArrayRef<RVec>  h_v,
378                                                    AtomLocality         atomLocality)
379 {
380     return impl_->copyVelocitiesFromGpu(h_v, atomLocality);
381 }
382
383
384 DeviceBuffer<float> StatePropagatorDataGpu::getForces()
385 {
386     return impl_->getForces();
387 }
388
389 void StatePropagatorDataGpu::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_f,
390                                              AtomLocality                          atomLocality)
391 {
392     return impl_->copyForcesToGpu(h_f, atomLocality);
393 }
394
395 void StatePropagatorDataGpu::copyForcesFromGpu(gmx::ArrayRef<RVec>  h_f,
396                                                AtomLocality         atomLocality)
397 {
398     return impl_->copyForcesFromGpu(h_f, atomLocality);
399 }
400
401 void* StatePropagatorDataGpu::getUpdateStream()
402 {
403     return impl_->getUpdateStream();
404 }
405
406 int StatePropagatorDataGpu::numAtomsLocal()
407 {
408     return impl_->numAtomsLocal();
409 }
410
411 int StatePropagatorDataGpu::numAtomsAll()
412 {
413     return impl_->numAtomsAll();
414 }
415
416 }      // namespace gmx
417
418 #endif // GMX_GPU == GMX_GPU_NONE