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