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