Add management for coordinates copy events into StatePropagatorDataGpu
[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     static_assert(GMX_GPU != GMX_GPU_NONE, "This object should only be constructed on the GPU code-paths.");
75     GMX_RELEASE_ASSERT(getenv("GMX_USE_GPU_BUFFER_OPS") == nullptr, "GPU buffer ops are not supported in this build.");
76
77     // TODO: Refactor when the StreamManager is introduced.
78     if (GMX_GPU == GMX_GPU_OPENCL)
79     {
80         GMX_ASSERT(deviceContext != nullptr, "GPU context should be set in OpenCL builds.");
81         GMX_ASSERT(pmeStream != nullptr, "GPU PME stream should be set in OpenCL builds.");
82
83         // The update stream is set to the PME stream in OpenCL, since PME stream is the only stream created in the PME context.
84         pmeStream_      = *static_cast<const CommandStream*>(pmeStream);
85         updateStream_   = *static_cast<const CommandStream*>(pmeStream);
86         deviceContext_  = *static_cast<const DeviceContext*>(deviceContext);
87         GMX_UNUSED_VALUE(localStream);
88         GMX_UNUSED_VALUE(nonLocalStream);
89     }
90
91     if (GMX_GPU == GMX_GPU_CUDA)
92     {
93         if (pmeStream != nullptr)
94         {
95             pmeStream_ = *static_cast<const CommandStream*>(pmeStream);
96         }
97         if (localStream != nullptr)
98         {
99             localStream_ = *static_cast<const CommandStream*>(localStream);
100         }
101         if (nonLocalStream != nullptr)
102         {
103             nonLocalStream_ = *static_cast<const CommandStream*>(nonLocalStream);
104         }
105
106         // TODO: The update stream should be created only when it is needed.
107 #if (GMX_GPU == GMX_GPU_CUDA)
108         cudaError_t stat;
109         stat = cudaStreamCreate(&updateStream_);
110         CU_RET_ERR(stat, "CUDA stream creation failed in StatePropagatorDataGpu");
111 #endif
112         GMX_UNUSED_VALUE(deviceContext);
113     }
114
115     // Map the atom locality to the stream that will be used for coordinates transfer.
116     // Same streams are used for H2D and D2H copies
117     xCopyStreams_[AtomLocality::Local]    = updateStream_;
118     xCopyStreams_[AtomLocality::NonLocal] = nonLocalStream_;
119     xCopyStreams_[AtomLocality::All]      = updateStream_;
120 }
121
122 StatePropagatorDataGpu::Impl::~Impl()
123 {
124 }
125
126 void StatePropagatorDataGpu::Impl::reinit(int numAtomsLocal, int numAtomsAll)
127 {
128     numAtomsLocal_ = numAtomsLocal;
129     numAtomsAll_   = numAtomsAll;
130
131     int numAtomsPadded;
132     if (paddingSize_ > 0)
133     {
134         numAtomsPadded = ((numAtomsAll_ + paddingSize_ - 1 ) / paddingSize_ )*paddingSize_;
135     }
136     else
137     {
138         numAtomsPadded = numAtomsAll_;
139     }
140
141     reallocateDeviceBuffer(&d_x_, DIM*numAtomsPadded, &d_xSize_, &d_xCapacity_, deviceContext_);
142
143     const size_t paddingAllocationSize = numAtomsPadded - numAtomsAll_;
144     if (paddingAllocationSize > 0)
145     {
146         // The PME stream is used here because the padding region of d_x_ is only in the PME task.
147         clearDeviceBufferAsync(&d_x_, DIM*numAtomsAll_, DIM*paddingAllocationSize, pmeStream_);
148     }
149
150     reallocateDeviceBuffer(&d_v_, DIM*numAtomsAll_, &d_vSize_, &d_vCapacity_, deviceContext_);
151     reallocateDeviceBuffer(&d_f_, DIM*numAtomsAll_, &d_fSize_, &d_fCapacity_, deviceContext_);
152
153 }
154
155 std::tuple<int, int> StatePropagatorDataGpu::Impl::getAtomRangesFromAtomLocality(AtomLocality  atomLocality)
156 {
157     int atomsStartAt   = 0;
158     int numAtomsToCopy = 0;
159     switch (atomLocality)
160     {
161         case AtomLocality::All:
162             atomsStartAt    = 0;
163             numAtomsToCopy  = numAtomsAll_;
164             break;
165         case AtomLocality::Local:
166             atomsStartAt    = 0;
167             numAtomsToCopy  = numAtomsLocal_;
168             break;
169         case AtomLocality::NonLocal:
170             atomsStartAt    = numAtomsLocal_;
171             numAtomsToCopy  = numAtomsAll_ - numAtomsLocal_;
172             break;
173         default:
174             GMX_RELEASE_ASSERT(false, "Wrong range of atoms requested in GPU state data manager. Should be All, Local or NonLocal.");
175     }
176     GMX_ASSERT(atomsStartAt   >= 0, "The first elemtnt to copy has negative index. Probably, the GPU propagator state was not initialized.");
177     GMX_ASSERT(numAtomsToCopy >= 0, "Number of atoms to copy is negative. Probably, the GPU propagator state was not initialized.");
178     return std::make_tuple(atomsStartAt, numAtomsToCopy);
179 }
180
181 void StatePropagatorDataGpu::Impl::copyToDevice(DeviceBuffer<float>                   d_data,
182                                                 const gmx::ArrayRef<const gmx::RVec>  h_data,
183                                                 int                                   dataSize,
184                                                 AtomLocality                          atomLocality,
185                                                 CommandStream                         commandStream)
186 {
187
188     GMX_UNUSED_VALUE(dataSize);
189
190     GMX_ASSERT(dataSize >= 0, "Trying to copy to device buffer before it was allocated.");
191
192     int atomsStartAt, numAtomsToCopy;
193     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
194
195     int elementsStartAt   = atomsStartAt*DIM;
196     int numElementsToCopy = numAtomsToCopy*DIM;
197
198     if (numAtomsToCopy != 0)
199     {
200         GMX_ASSERT(elementsStartAt + numElementsToCopy <= dataSize, "The device allocation is smaller than requested copy range.");
201         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(), "The host buffer is smaller than the requested copy range.");
202
203         copyToDeviceBuffer(&d_data, reinterpret_cast<const float *>(&h_data.data()[atomsStartAt]),
204                            elementsStartAt, numElementsToCopy,
205                            commandStream, transferKind_, nullptr);
206     }
207 }
208
209 void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec>  h_data,
210                                                   DeviceBuffer<float>       d_data,
211                                                   int                       dataSize,
212                                                   AtomLocality              atomLocality,
213                                                   CommandStream             commandStream)
214 {
215
216     GMX_UNUSED_VALUE(dataSize);
217
218     GMX_ASSERT(dataSize >= 0, "Trying to copy from device buffer before it was allocated.");
219
220     int atomsStartAt, numAtomsToCopy;
221     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
222
223     int elementsStartAt   = atomsStartAt*DIM;
224     int numElementsToCopy = numAtomsToCopy*DIM;
225
226     if (numAtomsToCopy != 0)
227     {
228         GMX_ASSERT(elementsStartAt + numElementsToCopy <= dataSize, "The device allocation is smaller than requested copy range.");
229         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(), "The host buffer is smaller than the requested copy range.");
230
231         copyFromDeviceBuffer(reinterpret_cast<float*>(&h_data.data()[atomsStartAt]), &d_data,
232                              elementsStartAt, numElementsToCopy,
233                              commandStream, transferKind_, nullptr);
234     }
235 }
236
237 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getCoordinates()
238 {
239     return d_x_;
240 }
241
242 void StatePropagatorDataGpu::Impl::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_x,
243                                                         AtomLocality                          atomLocality)
244 {
245     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
246     CommandStream commandStream = xCopyStreams_[atomLocality];
247     GMX_ASSERT(commandStream != nullptr, "No stream is valid for copying positions with given atom locality.");
248
249     copyToDevice(d_x_, h_x, d_xSize_, atomLocality, commandStream);
250
251     // markEvent is skipped in OpenCL as:
252     //   - it's not needed, copy is done in the same stream as the only consumer task (PME)
253     //   - we don't consume the events in OpenCL which is not allowed by GpuEventSynchronizer (would leak memory).
254     // TODO: remove this by adding an event-mark free flavor of this function
255     if (GMX_GPU == GMX_GPU_CUDA)
256     {
257         xReadyOnDevice_[atomLocality].markEvent(commandStream);
258         // TODO: Remove When event-based synchronization is introduced
259         gpuStreamSynchronize(commandStream);
260     }
261 }
262
263 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getCoordinatesReadyOnDeviceEvent(AtomLocality  atomLocality)
264 {
265     return &xReadyOnDevice_[atomLocality];
266 }
267
268 void StatePropagatorDataGpu::Impl::copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec>  h_x,
269                                                           AtomLocality              atomLocality)
270 {
271     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
272     CommandStream commandStream = xCopyStreams_[atomLocality];
273     GMX_ASSERT(commandStream != nullptr, "No stream is valid for copying positions with given atom locality.");
274
275     copyFromDevice(h_x, d_x_, d_xSize_, atomLocality, commandStream);
276     // TODO: Remove When event-based synchronization is introduced
277     gpuStreamSynchronize(commandStream);
278     // Note: unlike copyCoordinatesToGpu this is not used in OpenCL, and the conditional is not needed.
279     xReadyOnHost_[atomLocality].markEvent(commandStream);
280 }
281
282 void StatePropagatorDataGpu::Impl::waitCoordinatesReadyOnHost(AtomLocality  atomLocality)
283 {
284     xReadyOnHost_[atomLocality].waitForEvent();
285 }
286
287
288 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getVelocities()
289 {
290     return d_v_;
291 }
292
293 void StatePropagatorDataGpu::Impl::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_v,
294                                                        AtomLocality                          atomLocality)
295 {
296     // TODO: Use the correct stream
297     copyToDevice(d_v_, h_v, d_vSize_, atomLocality, nullptr);
298 }
299
300 void StatePropagatorDataGpu::Impl::copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec>  h_v,
301                                                          AtomLocality              atomLocality)
302 {
303     // TODO: Use the correct stream
304     copyFromDevice(h_v, d_v_, d_vSize_, atomLocality, nullptr);
305 }
306
307
308 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getForces()
309 {
310     return d_f_;
311 }
312
313 void StatePropagatorDataGpu::Impl::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_f,
314                                                    AtomLocality                          atomLocality)
315 {
316     // TODO: Use the correct stream
317     copyToDevice(d_f_, h_f, d_fSize_, atomLocality, nullptr);
318 }
319
320 void StatePropagatorDataGpu::Impl::copyForcesFromGpu(gmx::ArrayRef<gmx::RVec>  h_f,
321                                                      AtomLocality              atomLocality)
322 {
323     // TODO: Use the correct stream
324     copyFromDevice(h_f, d_f_, d_fSize_, atomLocality, nullptr);
325 }
326
327 void* StatePropagatorDataGpu::Impl::getUpdateStream()
328 {
329     return updateStream_;
330 }
331
332 int StatePropagatorDataGpu::Impl::numAtomsLocal()
333 {
334     return numAtomsLocal_;
335 }
336
337 int StatePropagatorDataGpu::Impl::numAtomsAll()
338 {
339     return numAtomsAll_;
340 }
341
342
343
344 StatePropagatorDataGpu::StatePropagatorDataGpu(const void        *pmeStream,
345                                                const void        *localStream,
346                                                const void        *nonLocalStream,
347                                                const void        *deviceContext,
348                                                GpuApiCallBehavior transferKind,
349                                                int                paddingSize)
350     : impl_(new Impl(pmeStream,
351                      localStream,
352                      nonLocalStream,
353                      deviceContext,
354                      transferKind,
355                      paddingSize))
356 {
357 }
358
359 StatePropagatorDataGpu::StatePropagatorDataGpu(StatePropagatorDataGpu && /* other */) noexcept = default;
360
361 StatePropagatorDataGpu &StatePropagatorDataGpu::operator=(StatePropagatorDataGpu && /* other */) noexcept = default;
362
363 StatePropagatorDataGpu::~StatePropagatorDataGpu() = default;
364
365
366 void StatePropagatorDataGpu::reinit(int numAtomsLocal, int numAtomsAll)
367 {
368     return impl_->reinit(numAtomsLocal, numAtomsAll);
369 }
370
371 std::tuple<int, int> StatePropagatorDataGpu::getAtomRangesFromAtomLocality(AtomLocality  atomLocality)
372 {
373     return impl_->getAtomRangesFromAtomLocality(atomLocality);
374 }
375
376
377 DeviceBuffer<float> StatePropagatorDataGpu::getCoordinates()
378 {
379     return impl_->getCoordinates();
380 }
381
382 void StatePropagatorDataGpu::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_x,
383                                                   AtomLocality                          atomLocality)
384 {
385     return impl_->copyCoordinatesToGpu(h_x, atomLocality);
386 }
387
388 GpuEventSynchronizer* StatePropagatorDataGpu::getCoordinatesReadyOnDeviceEvent(AtomLocality  atomLocality)
389 {
390     return impl_->getCoordinatesReadyOnDeviceEvent(atomLocality);
391 }
392
393 void StatePropagatorDataGpu::copyCoordinatesFromGpu(gmx::ArrayRef<RVec>  h_x,
394                                                     AtomLocality         atomLocality)
395 {
396     return impl_->copyCoordinatesFromGpu(h_x, atomLocality);
397 }
398
399 void StatePropagatorDataGpu::waitCoordinatesReadyOnHost(AtomLocality  atomLocality)
400 {
401     return impl_->waitCoordinatesReadyOnHost(atomLocality);
402 }
403
404
405 DeviceBuffer<float> StatePropagatorDataGpu::getVelocities()
406 {
407     return impl_->getVelocities();
408 }
409
410 void StatePropagatorDataGpu::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_v,
411                                                  AtomLocality                          atomLocality)
412 {
413     return impl_->copyVelocitiesToGpu(h_v, atomLocality);
414 }
415
416 void StatePropagatorDataGpu::copyVelocitiesFromGpu(gmx::ArrayRef<RVec>  h_v,
417                                                    AtomLocality         atomLocality)
418 {
419     return impl_->copyVelocitiesFromGpu(h_v, atomLocality);
420 }
421
422
423 DeviceBuffer<float> StatePropagatorDataGpu::getForces()
424 {
425     return impl_->getForces();
426 }
427
428 void StatePropagatorDataGpu::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec>  h_f,
429                                              AtomLocality                          atomLocality)
430 {
431     return impl_->copyForcesToGpu(h_f, atomLocality);
432 }
433
434 void StatePropagatorDataGpu::copyForcesFromGpu(gmx::ArrayRef<RVec>  h_f,
435                                                AtomLocality         atomLocality)
436 {
437     return impl_->copyForcesFromGpu(h_f, atomLocality);
438 }
439
440 void* StatePropagatorDataGpu::getUpdateStream()
441 {
442     return impl_->getUpdateStream();
443 }
444
445 int StatePropagatorDataGpu::numAtomsLocal()
446 {
447     return impl_->numAtomsLocal();
448 }
449
450 int StatePropagatorDataGpu::numAtomsAll()
451 {
452     return impl_->numAtomsAll();
453 }
454
455 }      // namespace gmx
456
457 #endif // GMX_GPU == GMX_GPU_NONE