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