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