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