Add Flake8 linting to gmxapi tests.
[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,
75                   "This object should only be constructed on the GPU code-paths.");
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,
116     // velocities and forces transfers. Same streams are used for H2D and D2H copies.
117     // Note, that nullptr stream is used here to indicate that the copy is not supported.
118     xCopyStreams_[AtomLocality::Local]    = updateStream_;
119     xCopyStreams_[AtomLocality::NonLocal] = nonLocalStream_;
120     xCopyStreams_[AtomLocality::All]      = updateStream_;
121
122     vCopyStreams_[AtomLocality::Local]    = updateStream_;
123     vCopyStreams_[AtomLocality::NonLocal] = nullptr;
124     vCopyStreams_[AtomLocality::All]      = updateStream_;
125
126     fCopyStreams_[AtomLocality::Local]    = localStream_;
127     fCopyStreams_[AtomLocality::NonLocal] = nonLocalStream_;
128     fCopyStreams_[AtomLocality::All]      = updateStream_;
129 }
130
131 StatePropagatorDataGpu::Impl::Impl(const void*        pmeStream,
132                                    const void*        deviceContext,
133                                    GpuApiCallBehavior transferKind,
134                                    int                paddingSize) :
135     transferKind_(transferKind),
136     paddingSize_(paddingSize)
137 {
138     static_assert(GMX_GPU != GMX_GPU_NONE,
139                   "This object should only be constructed on the GPU code-paths.");
140
141     if (GMX_GPU == GMX_GPU_OPENCL)
142     {
143         GMX_ASSERT(deviceContext != nullptr, "GPU context should be set in OpenCL builds.");
144         deviceContext_ = *static_cast<const DeviceContext*>(deviceContext);
145     }
146
147     GMX_ASSERT(pmeStream != nullptr, "GPU PME stream should be set.");
148     pmeStream_ = *static_cast<const CommandStream*>(pmeStream);
149
150     localStream_    = nullptr;
151     nonLocalStream_ = nullptr;
152     updateStream_   = nullptr;
153
154
155     // Only local/all coordinates are allowed to be copied in PME-only rank/ PME tests.
156     // This it temporary measure to make it safe to use this class in those cases.
157     xCopyStreams_[AtomLocality::Local]    = pmeStream_;
158     xCopyStreams_[AtomLocality::NonLocal] = nullptr;
159     xCopyStreams_[AtomLocality::All]      = pmeStream_;
160
161     vCopyStreams_[AtomLocality::Local]    = nullptr;
162     vCopyStreams_[AtomLocality::NonLocal] = nullptr;
163     vCopyStreams_[AtomLocality::All]      = nullptr;
164
165     fCopyStreams_[AtomLocality::Local]    = nullptr;
166     fCopyStreams_[AtomLocality::NonLocal] = nullptr;
167     fCopyStreams_[AtomLocality::All]      = nullptr;
168 }
169
170 StatePropagatorDataGpu::Impl::~Impl() {}
171
172 void StatePropagatorDataGpu::Impl::reinit(int numAtomsLocal, int numAtomsAll)
173 {
174     numAtomsLocal_ = numAtomsLocal;
175     numAtomsAll_   = numAtomsAll;
176
177     int numAtomsPadded;
178     if (paddingSize_ > 0)
179     {
180         numAtomsPadded = ((numAtomsAll_ + paddingSize_ - 1) / paddingSize_) * paddingSize_;
181     }
182     else
183     {
184         numAtomsPadded = numAtomsAll_;
185     }
186
187     reallocateDeviceBuffer(&d_x_, DIM * numAtomsPadded, &d_xSize_, &d_xCapacity_, deviceContext_);
188
189     const size_t paddingAllocationSize = numAtomsPadded - numAtomsAll_;
190     if (paddingAllocationSize > 0)
191     {
192         // The PME stream is used here because the padding region of d_x_ is only in the PME task.
193         clearDeviceBufferAsync(&d_x_, DIM * numAtomsAll_, DIM * paddingAllocationSize, pmeStream_);
194     }
195
196     reallocateDeviceBuffer(&d_v_, DIM * numAtomsAll_, &d_vSize_, &d_vCapacity_, deviceContext_);
197     const int d_fOldCapacity = d_fCapacity_;
198     reallocateDeviceBuffer(&d_f_, DIM * numAtomsAll_, &d_fSize_, &d_fCapacity_, deviceContext_);
199     // Clearing of the forces can be done in local stream since the nonlocal stream cannot reach
200     // the force accumulation stage before syncing with the local stream. Only done in CUDA,
201     // since the force buffer ops are not implemented in OpenCL.
202     if (GMX_GPU == GMX_GPU_CUDA && d_fCapacity_ != d_fOldCapacity)
203     {
204         clearDeviceBufferAsync(&d_f_, 0, d_fCapacity_, localStream_);
205     }
206 }
207
208 std::tuple<int, int> StatePropagatorDataGpu::Impl::getAtomRangesFromAtomLocality(AtomLocality atomLocality)
209 {
210     int atomsStartAt   = 0;
211     int numAtomsToCopy = 0;
212     switch (atomLocality)
213     {
214         case AtomLocality::All:
215             atomsStartAt   = 0;
216             numAtomsToCopy = numAtomsAll_;
217             break;
218         case AtomLocality::Local:
219             atomsStartAt   = 0;
220             numAtomsToCopy = numAtomsLocal_;
221             break;
222         case AtomLocality::NonLocal:
223             atomsStartAt   = numAtomsLocal_;
224             numAtomsToCopy = numAtomsAll_ - numAtomsLocal_;
225             break;
226         default:
227             GMX_RELEASE_ASSERT(false,
228                                "Wrong range of atoms requested in GPU state data manager. Should "
229                                "be All, Local or NonLocal.");
230     }
231     GMX_ASSERT(atomsStartAt >= 0,
232                "The first elemtnt to copy has negative index. Probably, the GPU propagator state "
233                "was not initialized.");
234     GMX_ASSERT(numAtomsToCopy >= 0,
235                "Number of atoms to copy is negative. Probably, the GPU propagator state was not "
236                "initialized.");
237     return std::make_tuple(atomsStartAt, numAtomsToCopy);
238 }
239
240 void StatePropagatorDataGpu::Impl::copyToDevice(DeviceBuffer<float>                  d_data,
241                                                 const gmx::ArrayRef<const gmx::RVec> h_data,
242                                                 int                                  dataSize,
243                                                 AtomLocality                         atomLocality,
244                                                 CommandStream                        commandStream)
245 {
246
247     GMX_UNUSED_VALUE(dataSize);
248
249     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
250
251     GMX_ASSERT(dataSize >= 0, "Trying to copy to device buffer before it was allocated.");
252
253     GMX_ASSERT(commandStream != nullptr,
254                "No stream is valid for copying with given atom locality.");
255
256     int atomsStartAt, numAtomsToCopy;
257     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
258
259     int elementsStartAt   = atomsStartAt * DIM;
260     int numElementsToCopy = numAtomsToCopy * DIM;
261
262     if (numAtomsToCopy != 0)
263     {
264         GMX_ASSERT(elementsStartAt + numElementsToCopy <= dataSize,
265                    "The device allocation is smaller than requested copy range.");
266         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(),
267                    "The host buffer is smaller than the requested copy range.");
268
269         copyToDeviceBuffer(&d_data, reinterpret_cast<const float*>(&h_data.data()[atomsStartAt]),
270                            elementsStartAt, numElementsToCopy, commandStream, transferKind_, nullptr);
271     }
272 }
273
274 void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
275                                                   DeviceBuffer<float>      d_data,
276                                                   int                      dataSize,
277                                                   AtomLocality             atomLocality,
278                                                   CommandStream            commandStream)
279 {
280
281     GMX_UNUSED_VALUE(dataSize);
282
283     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
284
285     GMX_ASSERT(dataSize >= 0, "Trying to copy from device buffer before it was allocated.");
286
287     GMX_ASSERT(commandStream != nullptr,
288                "No stream is valid for copying with given atom locality.");
289
290     int atomsStartAt, numAtomsToCopy;
291     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
292
293     int elementsStartAt   = atomsStartAt * DIM;
294     int numElementsToCopy = numAtomsToCopy * DIM;
295
296     if (numAtomsToCopy != 0)
297     {
298         GMX_ASSERT(elementsStartAt + numElementsToCopy <= dataSize,
299                    "The device allocation is smaller than requested copy range.");
300         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(),
301                    "The host buffer is smaller than the requested copy range.");
302
303         copyFromDeviceBuffer(reinterpret_cast<float*>(&h_data.data()[atomsStartAt]), &d_data,
304                              elementsStartAt, numElementsToCopy, commandStream, transferKind_, nullptr);
305     }
306 }
307
308 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getCoordinates()
309 {
310     return d_x_;
311 }
312
313 void StatePropagatorDataGpu::Impl::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec> h_x,
314                                                         AtomLocality atomLocality)
315 {
316     copyToDevice(d_x_, h_x, d_xSize_, atomLocality, xCopyStreams_[atomLocality]);
317
318     // markEvent is skipped in OpenCL as:
319     //   - it's not needed, copy is done in the same stream as the only consumer task (PME)
320     //   - we don't consume the events in OpenCL which is not allowed by GpuEventSynchronizer (would leak memory).
321     // TODO: remove this by adding an event-mark free flavor of this function
322     if (GMX_GPU == GMX_GPU_CUDA)
323     {
324         xReadyOnDevice_[atomLocality].markEvent(xCopyStreams_[atomLocality]);
325     }
326 }
327
328 GpuEventSynchronizer*
329 StatePropagatorDataGpu::Impl::getCoordinatesReadyOnDeviceEvent(AtomLocality atomLocality,
330                                                                const SimulationWorkload& simulationWork,
331                                                                const StepWorkload&       stepWork)
332 {
333     // The provider of the coordinates may be different for local atoms. If the update is offloaded
334     // and this is not a neighbor search step, then the consumer needs to wait for the update
335     // to complete. Otherwise, the coordinates are copied from the host and we need to wait for
336     // the copy event. Non-local coordinates are always provided by the H2D copy.
337     //
338     // TODO: This should be reconsidered to support the halo exchange.
339     //
340     // In OpenCL no events are used as coordinate sync is not necessary
341     if (GMX_GPU == GMX_GPU_OPENCL)
342     {
343         return nullptr;
344     }
345     if (atomLocality == AtomLocality::Local && simulationWork.useGpuUpdate && !stepWork.doNeighborSearch)
346     {
347         return &xUpdatedOnDevice_;
348     }
349     else
350     {
351         return &xReadyOnDevice_[atomLocality];
352     }
353 }
354
355 void StatePropagatorDataGpu::Impl::waitCoordinatesCopiedToDevice(AtomLocality atomLocality)
356 {
357     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
358     xReadyOnDevice_[atomLocality].waitForEvent();
359 }
360
361 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::xUpdatedOnDevice()
362 {
363     return &xUpdatedOnDevice_;
364 }
365
366 void StatePropagatorDataGpu::Impl::copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality)
367 {
368     copyFromDevice(h_x, d_x_, d_xSize_, atomLocality, xCopyStreams_[atomLocality]);
369     // Note: unlike copyCoordinatesToGpu this is not used in OpenCL, and the conditional is not needed.
370     xReadyOnHost_[atomLocality].markEvent(xCopyStreams_[atomLocality]);
371 }
372
373 void StatePropagatorDataGpu::Impl::waitCoordinatesReadyOnHost(AtomLocality atomLocality)
374 {
375     xReadyOnHost_[atomLocality].waitForEvent();
376 }
377
378
379 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getVelocities()
380 {
381     return d_v_;
382 }
383
384 void StatePropagatorDataGpu::Impl::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec> h_v,
385                                                        AtomLocality atomLocality)
386 {
387     copyToDevice(d_v_, h_v, d_vSize_, atomLocality, vCopyStreams_[atomLocality]);
388     vReadyOnDevice_[atomLocality].markEvent(vCopyStreams_[atomLocality]);
389 }
390
391 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality)
392 {
393     return &vReadyOnDevice_[atomLocality];
394 }
395
396
397 void StatePropagatorDataGpu::Impl::copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality)
398 {
399     copyFromDevice(h_v, d_v_, d_vSize_, atomLocality, vCopyStreams_[atomLocality]);
400     vReadyOnHost_[atomLocality].markEvent(vCopyStreams_[atomLocality]);
401 }
402
403 void StatePropagatorDataGpu::Impl::waitVelocitiesReadyOnHost(AtomLocality atomLocality)
404 {
405     vReadyOnHost_[atomLocality].waitForEvent();
406 }
407
408
409 DeviceBuffer<float> StatePropagatorDataGpu::Impl::getForces()
410 {
411     return d_f_;
412 }
413
414 void StatePropagatorDataGpu::Impl::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec> h_f,
415                                                    AtomLocality atomLocality)
416 {
417     copyToDevice(d_f_, h_f, d_fSize_, atomLocality, fCopyStreams_[atomLocality]);
418     fReadyOnDevice_[atomLocality].markEvent(fCopyStreams_[atomLocality]);
419 }
420
421 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getForcesReadyOnDeviceEvent(AtomLocality atomLocality,
422                                                                                 bool useGpuFBufferOps)
423 {
424     if ((atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal) && useGpuFBufferOps)
425     {
426         return &fReducedOnDevice_;
427     }
428     else
429     {
430         return &fReadyOnDevice_[atomLocality];
431     }
432 }
433
434 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::fReducedOnDevice()
435 {
436     return &fReducedOnDevice_;
437 }
438
439 void StatePropagatorDataGpu::Impl::copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality)
440 {
441     copyFromDevice(h_f, d_f_, d_fSize_, atomLocality, fCopyStreams_[atomLocality]);
442     fReadyOnHost_[atomLocality].markEvent(fCopyStreams_[atomLocality]);
443 }
444
445 void StatePropagatorDataGpu::Impl::waitForcesReadyOnHost(AtomLocality atomLocality)
446 {
447     fReadyOnHost_[atomLocality].waitForEvent();
448 }
449
450 void* StatePropagatorDataGpu::Impl::getUpdateStream()
451 {
452     return &updateStream_;
453 }
454
455 int StatePropagatorDataGpu::Impl::numAtomsLocal()
456 {
457     return numAtomsLocal_;
458 }
459
460 int StatePropagatorDataGpu::Impl::numAtomsAll()
461 {
462     return numAtomsAll_;
463 }
464
465
466 StatePropagatorDataGpu::StatePropagatorDataGpu(const void*        pmeStream,
467                                                const void*        localStream,
468                                                const void*        nonLocalStream,
469                                                const void*        deviceContext,
470                                                GpuApiCallBehavior transferKind,
471                                                int                paddingSize) :
472     impl_(new Impl(pmeStream, localStream, nonLocalStream, deviceContext, transferKind, paddingSize))
473 {
474 }
475
476 StatePropagatorDataGpu::StatePropagatorDataGpu(const void*        pmeStream,
477                                                const void*        deviceContext,
478                                                GpuApiCallBehavior transferKind,
479                                                int                paddingSize) :
480     impl_(new Impl(pmeStream, deviceContext, transferKind, paddingSize))
481 {
482 }
483
484 StatePropagatorDataGpu::StatePropagatorDataGpu(StatePropagatorDataGpu&& /* other */) noexcept = default;
485
486 StatePropagatorDataGpu& StatePropagatorDataGpu::operator=(StatePropagatorDataGpu&& /* other */) noexcept = default;
487
488 StatePropagatorDataGpu::~StatePropagatorDataGpu() = default;
489
490
491 void StatePropagatorDataGpu::reinit(int numAtomsLocal, int numAtomsAll)
492 {
493     return impl_->reinit(numAtomsLocal, numAtomsAll);
494 }
495
496 std::tuple<int, int> StatePropagatorDataGpu::getAtomRangesFromAtomLocality(AtomLocality atomLocality)
497 {
498     return impl_->getAtomRangesFromAtomLocality(atomLocality);
499 }
500
501
502 DeviceBuffer<float> StatePropagatorDataGpu::getCoordinates()
503 {
504     return impl_->getCoordinates();
505 }
506
507 void StatePropagatorDataGpu::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec> h_x,
508                                                   AtomLocality                         atomLocality)
509 {
510     return impl_->copyCoordinatesToGpu(h_x, atomLocality);
511 }
512
513 GpuEventSynchronizer*
514 StatePropagatorDataGpu::getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
515                                                          const SimulationWorkload& simulationWork,
516                                                          const StepWorkload&       stepWork)
517 {
518     return impl_->getCoordinatesReadyOnDeviceEvent(atomLocality, simulationWork, stepWork);
519 }
520
521 void StatePropagatorDataGpu::waitCoordinatesCopiedToDevice(AtomLocality atomLocality)
522 {
523     return impl_->waitCoordinatesCopiedToDevice(atomLocality);
524 }
525
526 GpuEventSynchronizer* StatePropagatorDataGpu::xUpdatedOnDevice()
527 {
528     return impl_->xUpdatedOnDevice();
529 }
530
531 void StatePropagatorDataGpu::copyCoordinatesFromGpu(gmx::ArrayRef<RVec> h_x, AtomLocality atomLocality)
532 {
533     return impl_->copyCoordinatesFromGpu(h_x, atomLocality);
534 }
535
536 void StatePropagatorDataGpu::waitCoordinatesReadyOnHost(AtomLocality atomLocality)
537 {
538     return impl_->waitCoordinatesReadyOnHost(atomLocality);
539 }
540
541
542 DeviceBuffer<float> StatePropagatorDataGpu::getVelocities()
543 {
544     return impl_->getVelocities();
545 }
546
547 void StatePropagatorDataGpu::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec> h_v,
548                                                  AtomLocality                         atomLocality)
549 {
550     return impl_->copyVelocitiesToGpu(h_v, atomLocality);
551 }
552
553 GpuEventSynchronizer* StatePropagatorDataGpu::getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality)
554 {
555     return impl_->getVelocitiesReadyOnDeviceEvent(atomLocality);
556 }
557
558 void StatePropagatorDataGpu::copyVelocitiesFromGpu(gmx::ArrayRef<RVec> h_v, AtomLocality atomLocality)
559 {
560     return impl_->copyVelocitiesFromGpu(h_v, atomLocality);
561 }
562
563 void StatePropagatorDataGpu::waitVelocitiesReadyOnHost(AtomLocality atomLocality)
564 {
565     return impl_->waitVelocitiesReadyOnHost(atomLocality);
566 }
567
568
569 DeviceBuffer<float> StatePropagatorDataGpu::getForces()
570 {
571     return impl_->getForces();
572 }
573
574 void StatePropagatorDataGpu::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality)
575 {
576     return impl_->copyForcesToGpu(h_f, atomLocality);
577 }
578
579 GpuEventSynchronizer* StatePropagatorDataGpu::getForcesReadyOnDeviceEvent(AtomLocality atomLocality,
580                                                                           bool useGpuFBufferOps)
581 {
582     return impl_->getForcesReadyOnDeviceEvent(atomLocality, useGpuFBufferOps);
583 }
584
585 GpuEventSynchronizer* StatePropagatorDataGpu::fReducedOnDevice()
586 {
587     return impl_->fReducedOnDevice();
588 }
589
590 void StatePropagatorDataGpu::copyForcesFromGpu(gmx::ArrayRef<RVec> h_f, AtomLocality atomLocality)
591 {
592     return impl_->copyForcesFromGpu(h_f, atomLocality);
593 }
594
595 void StatePropagatorDataGpu::waitForcesReadyOnHost(AtomLocality atomLocality)
596 {
597     return impl_->waitForcesReadyOnHost(atomLocality);
598 }
599
600
601 void* StatePropagatorDataGpu::getUpdateStream()
602 {
603     return impl_->getUpdateStream();
604 }
605
606 int StatePropagatorDataGpu::numAtomsLocal()
607 {
608     return impl_->numAtomsLocal();
609 }
610
611 int StatePropagatorDataGpu::numAtomsAll()
612 {
613     return impl_->numAtomsAll();
614 }
615
616 } // namespace gmx
617
618 #endif // GMX_GPU == GMX_GPU_NONE