bf927f2da297084caae44ef98a499038fabbccc7
[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,2020, 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 #    include "gromacs/gpu_utils/device_stream_manager.h"
50 #    include "gromacs/gpu_utils/devicebuffer.h"
51 #    include "gromacs/gpu_utils/gputraits.h"
52 #    include "gromacs/math/vectypes.h"
53 #    include "gromacs/mdtypes/state_propagator_data_gpu.h"
54 #    include "gromacs/timing/wallcycle.h"
55 #    include "gromacs/utility/classhelpers.h"
56
57 #    include "state_propagator_data_gpu_impl.h"
58
59
60 namespace gmx
61 {
62
63 StatePropagatorDataGpu::Impl::Impl(const DeviceStreamManager& deviceStreamManager,
64                                    GpuApiCallBehavior         transferKind,
65                                    int                        allocationBlockSizeDivisor,
66                                    gmx_wallcycle*             wcycle) :
67     deviceContext_(deviceStreamManager.context()),
68     transferKind_(transferKind),
69     allocationBlockSizeDivisor_(allocationBlockSizeDivisor),
70     wcycle_(wcycle)
71 {
72     static_assert(
73             GMX_GPU != GMX_GPU_NONE,
74             "GPU state propagator data object should only be constructed on the GPU code-paths.");
75
76     // We need to keep local copies for re-initialization.
77     pmeStream_      = &deviceStreamManager.stream(DeviceStreamType::Pme);
78     localStream_    = &deviceStreamManager.stream(DeviceStreamType::NonBondedLocal);
79     nonLocalStream_ = &deviceStreamManager.stream(DeviceStreamType::NonBondedNonLocal);
80     // PME stream is used in OpenCL for H2D coordinate transfer
81     if (GMX_GPU == GMX_GPU_OPENCL)
82     {
83         updateStream_ = &deviceStreamManager.stream(DeviceStreamType::Pme);
84     }
85     else
86     {
87         updateStream_ = &deviceStreamManager.stream(DeviceStreamType::UpdateAndConstraints);
88     }
89
90     // Map the atom locality to the stream that will be used for coordinates,
91     // velocities and forces transfers. Same streams are used for H2D and D2H copies.
92     // Note, that nullptr stream is used here to indicate that the copy is not supported.
93     xCopyStreams_[AtomLocality::Local]    = updateStream_;
94     xCopyStreams_[AtomLocality::NonLocal] = nonLocalStream_;
95     xCopyStreams_[AtomLocality::All]      = updateStream_;
96
97     vCopyStreams_[AtomLocality::Local]    = updateStream_;
98     vCopyStreams_[AtomLocality::NonLocal] = nullptr;
99     vCopyStreams_[AtomLocality::All]      = updateStream_;
100
101     fCopyStreams_[AtomLocality::Local]    = localStream_;
102     fCopyStreams_[AtomLocality::NonLocal] = nonLocalStream_;
103     fCopyStreams_[AtomLocality::All]      = updateStream_;
104 }
105
106 StatePropagatorDataGpu::Impl::Impl(const DeviceStream*  pmeStream,
107                                    const DeviceContext& deviceContext,
108                                    GpuApiCallBehavior   transferKind,
109                                    int                  allocationBlockSizeDivisor,
110                                    gmx_wallcycle*       wcycle) :
111     deviceContext_(deviceContext),
112     transferKind_(transferKind),
113     allocationBlockSizeDivisor_(allocationBlockSizeDivisor),
114     wcycle_(wcycle)
115 {
116     static_assert(
117             GMX_GPU != GMX_GPU_NONE,
118             "GPU state propagator data object should only be constructed on the GPU code-paths.");
119
120     GMX_ASSERT(pmeStream->isValid(), "GPU PME stream should be valid.");
121     pmeStream_      = pmeStream;
122     localStream_    = pmeStream; // For clearing the force buffer
123     nonLocalStream_ = nullptr;
124     updateStream_   = nullptr;
125
126
127     // Only local/all coordinates are allowed to be copied in PME-only rank/ PME tests.
128     // This it temporary measure to make it safe to use this class in those cases.
129     xCopyStreams_[AtomLocality::Local]    = pmeStream_;
130     xCopyStreams_[AtomLocality::NonLocal] = nullptr;
131     xCopyStreams_[AtomLocality::All]      = pmeStream_;
132
133     vCopyStreams_[AtomLocality::Local]    = nullptr;
134     vCopyStreams_[AtomLocality::NonLocal] = nullptr;
135     vCopyStreams_[AtomLocality::All]      = nullptr;
136
137     fCopyStreams_[AtomLocality::Local]    = nullptr;
138     fCopyStreams_[AtomLocality::NonLocal] = nullptr;
139     fCopyStreams_[AtomLocality::All]      = nullptr;
140 }
141
142 StatePropagatorDataGpu::Impl::~Impl() {}
143
144 void StatePropagatorDataGpu::Impl::reinit(int numAtomsLocal, int numAtomsAll)
145 {
146     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
147     wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
148
149     numAtomsLocal_ = numAtomsLocal;
150     numAtomsAll_   = numAtomsAll;
151
152     int numAtomsPadded;
153     if (allocationBlockSizeDivisor_ > 0)
154     {
155         numAtomsPadded = ((numAtomsAll_ + allocationBlockSizeDivisor_ - 1) / allocationBlockSizeDivisor_)
156                          * allocationBlockSizeDivisor_;
157     }
158     else
159     {
160         numAtomsPadded = numAtomsAll_;
161     }
162
163     reallocateDeviceBuffer(&d_x_, numAtomsPadded, &d_xSize_, &d_xCapacity_, deviceContext_);
164
165     const size_t paddingAllocationSize = numAtomsPadded - numAtomsAll_;
166     if (paddingAllocationSize > 0)
167     {
168         // The PME stream is used here because the padding region of d_x_ is only in the PME task.
169         clearDeviceBufferAsync(&d_x_, numAtomsAll_, paddingAllocationSize, *pmeStream_);
170     }
171
172     reallocateDeviceBuffer(&d_v_, numAtomsAll_, &d_vSize_, &d_vCapacity_, deviceContext_);
173     const int d_fOldCapacity = d_fCapacity_;
174     reallocateDeviceBuffer(&d_f_, numAtomsAll_, &d_fSize_, &d_fCapacity_, deviceContext_);
175     // Clearing of the forces can be done in local stream since the nonlocal stream cannot reach
176     // the force accumulation stage before syncing with the local stream. Only done in CUDA,
177     // since the force buffer ops are not implemented in OpenCL.
178     if (GMX_GPU == GMX_GPU_CUDA && d_fCapacity_ != d_fOldCapacity)
179     {
180         clearDeviceBufferAsync(&d_f_, 0, d_fCapacity_, *localStream_);
181     }
182
183     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
184     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
185 }
186
187 std::tuple<int, int> StatePropagatorDataGpu::Impl::getAtomRangesFromAtomLocality(AtomLocality atomLocality)
188 {
189     int atomsStartAt   = 0;
190     int numAtomsToCopy = 0;
191     switch (atomLocality)
192     {
193         case AtomLocality::All:
194             atomsStartAt   = 0;
195             numAtomsToCopy = numAtomsAll_;
196             break;
197         case AtomLocality::Local:
198             atomsStartAt   = 0;
199             numAtomsToCopy = numAtomsLocal_;
200             break;
201         case AtomLocality::NonLocal:
202             atomsStartAt   = numAtomsLocal_;
203             numAtomsToCopy = numAtomsAll_ - numAtomsLocal_;
204             break;
205         default:
206             GMX_RELEASE_ASSERT(false,
207                                "Wrong range of atoms requested in GPU state data manager. Should "
208                                "be All, Local or NonLocal.");
209     }
210     GMX_ASSERT(atomsStartAt >= 0,
211                "The first elemtnt to copy has negative index. Probably, the GPU propagator state "
212                "was not initialized.");
213     GMX_ASSERT(numAtomsToCopy >= 0,
214                "Number of atoms to copy is negative. Probably, the GPU propagator state was not "
215                "initialized.");
216     return std::make_tuple(atomsStartAt, numAtomsToCopy);
217 }
218
219 void StatePropagatorDataGpu::Impl::copyToDevice(DeviceBuffer<RVec>                   d_data,
220                                                 const gmx::ArrayRef<const gmx::RVec> h_data,
221                                                 int                                  dataSize,
222                                                 AtomLocality                         atomLocality,
223                                                 const DeviceStream&                  deviceStream)
224 {
225     GMX_UNUSED_VALUE(dataSize);
226
227     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
228
229     GMX_ASSERT(dataSize >= 0, "Trying to copy to device buffer before it was allocated.");
230
231     GMX_ASSERT(deviceStream.isValid(), "No stream is valid for copying with given atom locality.");
232     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
233     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
234
235     int atomsStartAt, numAtomsToCopy;
236     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
237
238     if (numAtomsToCopy != 0)
239     {
240         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= dataSize,
241                    "The device allocation is smaller than requested copy range.");
242         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(),
243                    "The host buffer is smaller than the requested copy range.");
244
245         copyToDeviceBuffer(&d_data, reinterpret_cast<const RVec*>(&h_data.data()[atomsStartAt]),
246                            atomsStartAt, numAtomsToCopy, deviceStream, transferKind_, nullptr);
247     }
248
249     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
250     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
251 }
252
253 void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
254                                                   DeviceBuffer<RVec>       d_data,
255                                                   int                      dataSize,
256                                                   AtomLocality             atomLocality,
257                                                   const DeviceStream&      deviceStream)
258 {
259     GMX_UNUSED_VALUE(dataSize);
260
261     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
262
263     GMX_ASSERT(dataSize >= 0, "Trying to copy from device buffer before it was allocated.");
264
265     GMX_ASSERT(deviceStream.isValid(), "No stream is valid for copying with given atom locality.");
266     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
267     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
268
269     int atomsStartAt, numAtomsToCopy;
270     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
271
272     if (numAtomsToCopy != 0)
273     {
274         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= dataSize,
275                    "The device allocation is smaller than requested copy range.");
276         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(),
277                    "The host buffer is smaller than the requested copy range.");
278
279         copyFromDeviceBuffer(reinterpret_cast<RVec*>(&h_data.data()[atomsStartAt]), &d_data,
280                              atomsStartAt, numAtomsToCopy, deviceStream, transferKind_, nullptr);
281     }
282
283     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
284     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
285 }
286
287 DeviceBuffer<RVec> StatePropagatorDataGpu::Impl::getCoordinates()
288 {
289     return d_x_;
290 }
291
292 void StatePropagatorDataGpu::Impl::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec> h_x,
293                                                         AtomLocality atomLocality)
294 {
295     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
296     const DeviceStream* deviceStream = xCopyStreams_[atomLocality];
297     GMX_ASSERT(deviceStream != nullptr,
298                "No stream is valid for copying positions with given atom locality.");
299
300     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
301     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
302
303     copyToDevice(d_x_, h_x, d_xSize_, atomLocality, *deviceStream);
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(*deviceStream);
312     }
313
314     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
315     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
316 }
317
318 GpuEventSynchronizer*
319 StatePropagatorDataGpu::Impl::getCoordinatesReadyOnDeviceEvent(AtomLocality atomLocality,
320                                                                const SimulationWorkload& simulationWork,
321                                                                const StepWorkload&       stepWork)
322 {
323     // The provider of the coordinates may be different for local atoms. If the update is offloaded
324     // and this is not a neighbor search step, then the consumer needs to wait for the update
325     // to complete. Otherwise, the coordinates are copied from the host and we need to wait for
326     // the copy event. Non-local coordinates are always provided by the H2D copy.
327     //
328     // TODO: This should be reconsidered to support the halo exchange.
329     //
330     // In OpenCL no events are used as coordinate sync is not necessary
331     if (GMX_GPU == GMX_GPU_OPENCL)
332     {
333         return nullptr;
334     }
335     if (atomLocality == AtomLocality::Local && simulationWork.useGpuUpdate && !stepWork.doNeighborSearch)
336     {
337         return &xUpdatedOnDevice_;
338     }
339     else
340     {
341         return &xReadyOnDevice_[atomLocality];
342     }
343 }
344
345 void StatePropagatorDataGpu::Impl::waitCoordinatesCopiedToDevice(AtomLocality atomLocality)
346 {
347     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
348     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
349     xReadyOnDevice_[atomLocality].waitForEvent();
350     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
351 }
352
353 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::xUpdatedOnDevice()
354 {
355     return &xUpdatedOnDevice_;
356 }
357
358 void StatePropagatorDataGpu::Impl::copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality)
359 {
360     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
361     const DeviceStream* deviceStream = xCopyStreams_[atomLocality];
362     GMX_ASSERT(deviceStream != nullptr,
363                "No stream is valid for copying positions with given atom locality.");
364
365     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
366     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
367
368     copyFromDevice(h_x, d_x_, d_xSize_, atomLocality, *deviceStream);
369     // Note: unlike copyCoordinatesToGpu this is not used in OpenCL, and the conditional is not needed.
370     xReadyOnHost_[atomLocality].markEvent(*deviceStream);
371
372     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
373     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
374 }
375
376 void StatePropagatorDataGpu::Impl::waitCoordinatesReadyOnHost(AtomLocality atomLocality)
377 {
378     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
379     xReadyOnHost_[atomLocality].waitForEvent();
380     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
381 }
382
383
384 DeviceBuffer<RVec> StatePropagatorDataGpu::Impl::getVelocities()
385 {
386     return d_v_;
387 }
388
389 void StatePropagatorDataGpu::Impl::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec> h_v,
390                                                        AtomLocality atomLocality)
391 {
392     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
393     const DeviceStream* deviceStream = vCopyStreams_[atomLocality];
394     GMX_ASSERT(deviceStream != nullptr,
395                "No stream is valid for copying velocities with given atom locality.");
396
397     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
398     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
399
400     copyToDevice(d_v_, h_v, d_vSize_, atomLocality, *deviceStream);
401     vReadyOnDevice_[atomLocality].markEvent(*deviceStream);
402
403     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
404     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
405 }
406
407 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality)
408 {
409     return &vReadyOnDevice_[atomLocality];
410 }
411
412
413 void StatePropagatorDataGpu::Impl::copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality)
414 {
415     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
416     const DeviceStream* deviceStream = vCopyStreams_[atomLocality];
417     GMX_ASSERT(deviceStream != nullptr,
418                "No stream is valid for copying velocities with given atom locality.");
419
420     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
421     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
422
423     copyFromDevice(h_v, d_v_, d_vSize_, atomLocality, *deviceStream);
424     vReadyOnHost_[atomLocality].markEvent(*deviceStream);
425
426     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
427     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
428 }
429
430 void StatePropagatorDataGpu::Impl::waitVelocitiesReadyOnHost(AtomLocality atomLocality)
431 {
432     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
433     vReadyOnHost_[atomLocality].waitForEvent();
434     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
435 }
436
437
438 DeviceBuffer<RVec> StatePropagatorDataGpu::Impl::getForces()
439 {
440     return d_f_;
441 }
442
443 void StatePropagatorDataGpu::Impl::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec> h_f,
444                                                    AtomLocality atomLocality)
445 {
446     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
447     const DeviceStream* deviceStream = fCopyStreams_[atomLocality];
448     GMX_ASSERT(deviceStream != nullptr,
449                "No stream is valid for copying forces with given atom locality.");
450
451     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
452     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
453
454     copyToDevice(d_f_, h_f, d_fSize_, atomLocality, *deviceStream);
455     fReadyOnDevice_[atomLocality].markEvent(*deviceStream);
456
457     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
458     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
459 }
460
461 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getForcesReadyOnDeviceEvent(AtomLocality atomLocality,
462                                                                                 bool useGpuFBufferOps)
463 {
464     if ((atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal) && useGpuFBufferOps)
465     {
466         return &fReducedOnDevice_;
467     }
468     else
469     {
470         return &fReadyOnDevice_[atomLocality];
471     }
472 }
473
474 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::fReducedOnDevice()
475 {
476     return &fReducedOnDevice_;
477 }
478
479 void StatePropagatorDataGpu::Impl::copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality)
480 {
481     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
482     const DeviceStream* deviceStream = fCopyStreams_[atomLocality];
483     GMX_ASSERT(deviceStream != nullptr,
484                "No stream is valid for copying forces with given atom locality.");
485
486     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
487     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
488
489     copyFromDevice(h_f, d_f_, d_fSize_, atomLocality, *deviceStream);
490     fReadyOnHost_[atomLocality].markEvent(*deviceStream);
491
492     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
493     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
494 }
495
496 void StatePropagatorDataGpu::Impl::waitForcesReadyOnHost(AtomLocality atomLocality)
497 {
498     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
499     fReadyOnHost_[atomLocality].waitForEvent();
500     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
501 }
502
503 const DeviceStream* StatePropagatorDataGpu::Impl::getUpdateStream()
504 {
505     return updateStream_;
506 }
507
508 int StatePropagatorDataGpu::Impl::numAtomsLocal()
509 {
510     return numAtomsLocal_;
511 }
512
513 int StatePropagatorDataGpu::Impl::numAtomsAll()
514 {
515     return numAtomsAll_;
516 }
517
518
519 StatePropagatorDataGpu::StatePropagatorDataGpu(const DeviceStreamManager& deviceStreamManager,
520                                                GpuApiCallBehavior         transferKind,
521                                                int            allocationBlockSizeDivisor,
522                                                gmx_wallcycle* wcycle) :
523     impl_(new Impl(deviceStreamManager, transferKind, allocationBlockSizeDivisor, wcycle))
524 {
525 }
526
527 StatePropagatorDataGpu::StatePropagatorDataGpu(const DeviceStream*  pmeStream,
528                                                const DeviceContext& deviceContext,
529                                                GpuApiCallBehavior   transferKind,
530                                                int                  allocationBlockSizeDivisor,
531                                                gmx_wallcycle*       wcycle) :
532     impl_(new Impl(pmeStream, deviceContext, transferKind, allocationBlockSizeDivisor, wcycle))
533 {
534 }
535
536 StatePropagatorDataGpu::StatePropagatorDataGpu(StatePropagatorDataGpu&& /* other */) noexcept = default;
537
538 StatePropagatorDataGpu& StatePropagatorDataGpu::operator=(StatePropagatorDataGpu&& /* other */) noexcept = default;
539
540 StatePropagatorDataGpu::~StatePropagatorDataGpu() = default;
541
542
543 void StatePropagatorDataGpu::reinit(int numAtomsLocal, int numAtomsAll)
544 {
545     return impl_->reinit(numAtomsLocal, numAtomsAll);
546 }
547
548 std::tuple<int, int> StatePropagatorDataGpu::getAtomRangesFromAtomLocality(AtomLocality atomLocality)
549 {
550     return impl_->getAtomRangesFromAtomLocality(atomLocality);
551 }
552
553
554 DeviceBuffer<RVec> StatePropagatorDataGpu::getCoordinates()
555 {
556     return impl_->getCoordinates();
557 }
558
559 void StatePropagatorDataGpu::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec> h_x,
560                                                   AtomLocality                         atomLocality)
561 {
562     return impl_->copyCoordinatesToGpu(h_x, atomLocality);
563 }
564
565 GpuEventSynchronizer*
566 StatePropagatorDataGpu::getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
567                                                          const SimulationWorkload& simulationWork,
568                                                          const StepWorkload&       stepWork)
569 {
570     return impl_->getCoordinatesReadyOnDeviceEvent(atomLocality, simulationWork, stepWork);
571 }
572
573 void StatePropagatorDataGpu::waitCoordinatesCopiedToDevice(AtomLocality atomLocality)
574 {
575     return impl_->waitCoordinatesCopiedToDevice(atomLocality);
576 }
577
578 GpuEventSynchronizer* StatePropagatorDataGpu::xUpdatedOnDevice()
579 {
580     return impl_->xUpdatedOnDevice();
581 }
582
583 void StatePropagatorDataGpu::copyCoordinatesFromGpu(gmx::ArrayRef<RVec> h_x, AtomLocality atomLocality)
584 {
585     return impl_->copyCoordinatesFromGpu(h_x, atomLocality);
586 }
587
588 void StatePropagatorDataGpu::waitCoordinatesReadyOnHost(AtomLocality atomLocality)
589 {
590     return impl_->waitCoordinatesReadyOnHost(atomLocality);
591 }
592
593
594 DeviceBuffer<RVec> StatePropagatorDataGpu::getVelocities()
595 {
596     return impl_->getVelocities();
597 }
598
599 void StatePropagatorDataGpu::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec> h_v,
600                                                  AtomLocality                         atomLocality)
601 {
602     return impl_->copyVelocitiesToGpu(h_v, atomLocality);
603 }
604
605 GpuEventSynchronizer* StatePropagatorDataGpu::getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality)
606 {
607     return impl_->getVelocitiesReadyOnDeviceEvent(atomLocality);
608 }
609
610 void StatePropagatorDataGpu::copyVelocitiesFromGpu(gmx::ArrayRef<RVec> h_v, AtomLocality atomLocality)
611 {
612     return impl_->copyVelocitiesFromGpu(h_v, atomLocality);
613 }
614
615 void StatePropagatorDataGpu::waitVelocitiesReadyOnHost(AtomLocality atomLocality)
616 {
617     return impl_->waitVelocitiesReadyOnHost(atomLocality);
618 }
619
620
621 DeviceBuffer<RVec> StatePropagatorDataGpu::getForces()
622 {
623     return impl_->getForces();
624 }
625
626 void StatePropagatorDataGpu::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality)
627 {
628     return impl_->copyForcesToGpu(h_f, atomLocality);
629 }
630
631 GpuEventSynchronizer* StatePropagatorDataGpu::getForcesReadyOnDeviceEvent(AtomLocality atomLocality,
632                                                                           bool useGpuFBufferOps)
633 {
634     return impl_->getForcesReadyOnDeviceEvent(atomLocality, useGpuFBufferOps);
635 }
636
637 GpuEventSynchronizer* StatePropagatorDataGpu::fReducedOnDevice()
638 {
639     return impl_->fReducedOnDevice();
640 }
641
642 void StatePropagatorDataGpu::copyForcesFromGpu(gmx::ArrayRef<RVec> h_f, AtomLocality atomLocality)
643 {
644     return impl_->copyForcesFromGpu(h_f, atomLocality);
645 }
646
647 void StatePropagatorDataGpu::waitForcesReadyOnHost(AtomLocality atomLocality)
648 {
649     return impl_->waitForcesReadyOnHost(atomLocality);
650 }
651
652
653 const DeviceStream* StatePropagatorDataGpu::getUpdateStream()
654 {
655     return impl_->getUpdateStream();
656 }
657
658 int StatePropagatorDataGpu::numAtomsLocal()
659 {
660     return impl_->numAtomsLocal();
661 }
662
663 int StatePropagatorDataGpu::numAtomsAll()
664 {
665     return impl_->numAtomsAll();
666 }
667
668 } // namespace gmx
669
670 #endif // GMX_GPU == GMX_GPU_NONE