Store DeviceInformation inside the DeviceContext class
[alexxy/gromacs.git] / src / gromacs / mdtypes / state_propagator_data_gpu_impl_gpu.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  *
37  * \brief Definitions of interfaces for GPU state data propagator object.
38  *
39  * \author Artem Zhmurov <zhmurov@gmail.com>
40  *
41  * \ingroup module_mdtypes
42  */
43 #include "gmxpre.h"
44
45 #include "config.h"
46
47 #if GMX_GPU != GMX_GPU_NONE
48
49 #    if GMX_GPU == GMX_GPU_CUDA
50 #        include "gromacs/gpu_utils/cudautils.cuh"
51 #    endif
52 #    include "gromacs/gpu_utils/devicebuffer.h"
53 #    include "gromacs/gpu_utils/gputraits.h"
54 #    if GMX_GPU == GMX_GPU_OPENCL
55 #        include "gromacs/gpu_utils/oclutils.h"
56 #    endif
57 #    include "gromacs/math/vectypes.h"
58 #    include "gromacs/mdtypes/state_propagator_data_gpu.h"
59 #    include "gromacs/timing/wallcycle.h"
60 #    include "gromacs/utility/classhelpers.h"
61
62 #    include "state_propagator_data_gpu_impl.h"
63
64
65 namespace gmx
66 {
67
68 StatePropagatorDataGpu::Impl::Impl(const DeviceStream*  pmeStream,
69                                    const DeviceStream*  localStream,
70                                    const DeviceStream*  nonLocalStream,
71                                    const DeviceContext& deviceContext,
72                                    GpuApiCallBehavior   transferKind,
73                                    int                  allocationBlockSizeDivisor,
74                                    gmx_wallcycle*       wcycle) :
75     deviceContext_(deviceContext),
76     transferKind_(transferKind),
77     allocationBlockSizeDivisor_(allocationBlockSizeDivisor),
78     wcycle_(wcycle)
79 {
80     static_assert(GMX_GPU != GMX_GPU_NONE,
81                   "This object should only be constructed on the GPU code-paths.");
82
83     // TODO: Refactor when the StreamManager is introduced.
84     if (GMX_GPU == GMX_GPU_OPENCL)
85     {
86         GMX_ASSERT(pmeStream != nullptr, "GPU PME stream should be set in OpenCL builds.");
87
88         // The update stream is set to the PME stream in OpenCL, since PME stream is the only stream created in the PME context.
89         pmeStream_    = pmeStream;
90         updateStream_ = pmeStream;
91         GMX_UNUSED_VALUE(localStream);
92         GMX_UNUSED_VALUE(nonLocalStream);
93     }
94
95     if (GMX_GPU == GMX_GPU_CUDA)
96     {
97         if (pmeStream != nullptr)
98         {
99             pmeStream_ = pmeStream;
100         }
101         if (localStream != nullptr)
102         {
103             localStream_ = localStream;
104         }
105         if (nonLocalStream != nullptr)
106         {
107             nonLocalStream_ = nonLocalStream;
108         }
109
110         // TODO: The update stream should be created only when it is needed.
111 #    if (GMX_GPU == GMX_GPU_CUDA)
112         // In CUDA we only need priority to create stream.
113         // (note that this will be moved from here in the follow-up patch)
114         updateStreamOwn_.init(deviceContext, DeviceStreamPriority::Normal, false);
115         updateStream_ = &updateStreamOwn_;
116 #    endif
117     }
118
119     // Map the atom locality to the stream that will be used for coordinates,
120     // velocities and forces transfers. Same streams are used for H2D and D2H copies.
121     // Note, that nullptr stream is used here to indicate that the copy is not supported.
122     xCopyStreams_[AtomLocality::Local]    = updateStream_;
123     xCopyStreams_[AtomLocality::NonLocal] = nonLocalStream_;
124     xCopyStreams_[AtomLocality::All]      = updateStream_;
125
126     vCopyStreams_[AtomLocality::Local]    = updateStream_;
127     vCopyStreams_[AtomLocality::NonLocal] = nullptr;
128     vCopyStreams_[AtomLocality::All]      = updateStream_;
129
130     fCopyStreams_[AtomLocality::Local]    = localStream_;
131     fCopyStreams_[AtomLocality::NonLocal] = nonLocalStream_;
132     fCopyStreams_[AtomLocality::All]      = updateStream_;
133 }
134
135 StatePropagatorDataGpu::Impl::Impl(const DeviceStream*  pmeStream,
136                                    const DeviceContext& deviceContext,
137                                    GpuApiCallBehavior   transferKind,
138                                    int                  allocationBlockSizeDivisor,
139                                    gmx_wallcycle*       wcycle) :
140     deviceContext_(deviceContext),
141     transferKind_(transferKind),
142     allocationBlockSizeDivisor_(allocationBlockSizeDivisor),
143     wcycle_(wcycle)
144 {
145     static_assert(GMX_GPU != GMX_GPU_NONE,
146                   "This object should only be constructed on the GPU code-paths.");
147
148     GMX_ASSERT(pmeStream != nullptr, "GPU PME stream should be set.");
149     pmeStream_      = pmeStream;
150     localStream_    = pmeStream; // For clearing the force buffer
151     nonLocalStream_ = nullptr;
152     updateStream_   = nullptr;
153
154
155     // Only local/all coordinates are allowed to be copied in PME-only rank/ PME tests.
156     // This it temporary measure to make it safe to use this class in those cases.
157     xCopyStreams_[AtomLocality::Local]    = pmeStream_;
158     xCopyStreams_[AtomLocality::NonLocal] = nullptr;
159     xCopyStreams_[AtomLocality::All]      = pmeStream_;
160
161     vCopyStreams_[AtomLocality::Local]    = nullptr;
162     vCopyStreams_[AtomLocality::NonLocal] = nullptr;
163     vCopyStreams_[AtomLocality::All]      = nullptr;
164
165     fCopyStreams_[AtomLocality::Local]    = nullptr;
166     fCopyStreams_[AtomLocality::NonLocal] = nullptr;
167     fCopyStreams_[AtomLocality::All]      = nullptr;
168 }
169
170 StatePropagatorDataGpu::Impl::~Impl() {}
171
172 void StatePropagatorDataGpu::Impl::reinit(int numAtomsLocal, int numAtomsAll)
173 {
174     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
175     wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
176
177     numAtomsLocal_ = numAtomsLocal;
178     numAtomsAll_   = numAtomsAll;
179
180     int numAtomsPadded;
181     if (allocationBlockSizeDivisor_ > 0)
182     {
183         numAtomsPadded = ((numAtomsAll_ + allocationBlockSizeDivisor_ - 1) / allocationBlockSizeDivisor_)
184                          * allocationBlockSizeDivisor_;
185     }
186     else
187     {
188         numAtomsPadded = numAtomsAll_;
189     }
190
191     reallocateDeviceBuffer(&d_x_, numAtomsPadded, &d_xSize_, &d_xCapacity_, deviceContext_);
192
193     const size_t paddingAllocationSize = numAtomsPadded - numAtomsAll_;
194     if (paddingAllocationSize > 0)
195     {
196         // The PME stream is used here because the padding region of d_x_ is only in the PME task.
197         clearDeviceBufferAsync(&d_x_, numAtomsAll_, paddingAllocationSize, *pmeStream_);
198     }
199
200     reallocateDeviceBuffer(&d_v_, numAtomsAll_, &d_vSize_, &d_vCapacity_, deviceContext_);
201     const int d_fOldCapacity = d_fCapacity_;
202     reallocateDeviceBuffer(&d_f_, numAtomsAll_, &d_fSize_, &d_fCapacity_, deviceContext_);
203     // Clearing of the forces can be done in local stream since the nonlocal stream cannot reach
204     // the force accumulation stage before syncing with the local stream. Only done in CUDA,
205     // since the force buffer ops are not implemented in OpenCL.
206     if (GMX_GPU == GMX_GPU_CUDA && d_fCapacity_ != d_fOldCapacity)
207     {
208         clearDeviceBufferAsync(&d_f_, 0, d_fCapacity_, *localStream_);
209     }
210
211     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
212     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
213 }
214
215 std::tuple<int, int> StatePropagatorDataGpu::Impl::getAtomRangesFromAtomLocality(AtomLocality atomLocality)
216 {
217     int atomsStartAt   = 0;
218     int numAtomsToCopy = 0;
219     switch (atomLocality)
220     {
221         case AtomLocality::All:
222             atomsStartAt   = 0;
223             numAtomsToCopy = numAtomsAll_;
224             break;
225         case AtomLocality::Local:
226             atomsStartAt   = 0;
227             numAtomsToCopy = numAtomsLocal_;
228             break;
229         case AtomLocality::NonLocal:
230             atomsStartAt   = numAtomsLocal_;
231             numAtomsToCopy = numAtomsAll_ - numAtomsLocal_;
232             break;
233         default:
234             GMX_RELEASE_ASSERT(false,
235                                "Wrong range of atoms requested in GPU state data manager. Should "
236                                "be All, Local or NonLocal.");
237     }
238     GMX_ASSERT(atomsStartAt >= 0,
239                "The first elemtnt to copy has negative index. Probably, the GPU propagator state "
240                "was not initialized.");
241     GMX_ASSERT(numAtomsToCopy >= 0,
242                "Number of atoms to copy is negative. Probably, the GPU propagator state was not "
243                "initialized.");
244     return std::make_tuple(atomsStartAt, numAtomsToCopy);
245 }
246
247 void StatePropagatorDataGpu::Impl::copyToDevice(DeviceBuffer<RVec>                   d_data,
248                                                 const gmx::ArrayRef<const gmx::RVec> h_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 to device buffer before it was allocated.");
258
259     GMX_ASSERT(deviceStream.stream() != nullptr,
260                "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         copyToDeviceBuffer(&d_data, reinterpret_cast<const RVec*>(&h_data.data()[atomsStartAt]),
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 void StatePropagatorDataGpu::Impl::copyFromDevice(gmx::ArrayRef<gmx::RVec> h_data,
283                                                   DeviceBuffer<RVec>       d_data,
284                                                   int                      dataSize,
285                                                   AtomLocality             atomLocality,
286                                                   const DeviceStream&      deviceStream)
287 {
288     GMX_UNUSED_VALUE(dataSize);
289
290     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
291
292     GMX_ASSERT(dataSize >= 0, "Trying to copy from device buffer before it was allocated.");
293
294     GMX_ASSERT(deviceStream.stream() != nullptr,
295                "No stream is valid for copying with given atom locality.");
296     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
297     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
298
299     int atomsStartAt, numAtomsToCopy;
300     std::tie(atomsStartAt, numAtomsToCopy) = getAtomRangesFromAtomLocality(atomLocality);
301
302     if (numAtomsToCopy != 0)
303     {
304         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= dataSize,
305                    "The device allocation is smaller than requested copy range.");
306         GMX_ASSERT(atomsStartAt + numAtomsToCopy <= h_data.ssize(),
307                    "The host buffer is smaller than the requested copy range.");
308
309         copyFromDeviceBuffer(reinterpret_cast<RVec*>(&h_data.data()[atomsStartAt]), &d_data,
310                              atomsStartAt, numAtomsToCopy, deviceStream, transferKind_, nullptr);
311     }
312
313     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
314     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
315 }
316
317 DeviceBuffer<RVec> StatePropagatorDataGpu::Impl::getCoordinates()
318 {
319     return d_x_;
320 }
321
322 void StatePropagatorDataGpu::Impl::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec> h_x,
323                                                         AtomLocality atomLocality)
324 {
325     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
326     const DeviceStream* deviceStream = xCopyStreams_[atomLocality];
327     GMX_ASSERT(deviceStream != nullptr,
328                "No stream is valid for copying positions with given atom locality.");
329
330     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
331     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
332
333     copyToDevice(d_x_, h_x, d_xSize_, atomLocality, *deviceStream);
334
335     // markEvent is skipped in OpenCL as:
336     //   - it's not needed, copy is done in the same stream as the only consumer task (PME)
337     //   - we don't consume the events in OpenCL which is not allowed by GpuEventSynchronizer (would leak memory).
338     // TODO: remove this by adding an event-mark free flavor of this function
339     if (GMX_GPU == GMX_GPU_CUDA)
340     {
341         xReadyOnDevice_[atomLocality].markEvent(*deviceStream);
342     }
343
344     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
345     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
346 }
347
348 GpuEventSynchronizer*
349 StatePropagatorDataGpu::Impl::getCoordinatesReadyOnDeviceEvent(AtomLocality atomLocality,
350                                                                const SimulationWorkload& simulationWork,
351                                                                const StepWorkload&       stepWork)
352 {
353     // The provider of the coordinates may be different for local atoms. If the update is offloaded
354     // and this is not a neighbor search step, then the consumer needs to wait for the update
355     // to complete. Otherwise, the coordinates are copied from the host and we need to wait for
356     // the copy event. Non-local coordinates are always provided by the H2D copy.
357     //
358     // TODO: This should be reconsidered to support the halo exchange.
359     //
360     // In OpenCL no events are used as coordinate sync is not necessary
361     if (GMX_GPU == GMX_GPU_OPENCL)
362     {
363         return nullptr;
364     }
365     if (atomLocality == AtomLocality::Local && simulationWork.useGpuUpdate && !stepWork.doNeighborSearch)
366     {
367         return &xUpdatedOnDevice_;
368     }
369     else
370     {
371         return &xReadyOnDevice_[atomLocality];
372     }
373 }
374
375 void StatePropagatorDataGpu::Impl::waitCoordinatesCopiedToDevice(AtomLocality atomLocality)
376 {
377     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
378     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
379     xReadyOnDevice_[atomLocality].waitForEvent();
380     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
381 }
382
383 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::xUpdatedOnDevice()
384 {
385     return &xUpdatedOnDevice_;
386 }
387
388 void StatePropagatorDataGpu::Impl::copyCoordinatesFromGpu(gmx::ArrayRef<gmx::RVec> h_x, AtomLocality atomLocality)
389 {
390     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
391     const DeviceStream* deviceStream = xCopyStreams_[atomLocality];
392     GMX_ASSERT(deviceStream != nullptr,
393                "No stream is valid for copying positions with given atom locality.");
394
395     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
396     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
397
398     copyFromDevice(h_x, d_x_, d_xSize_, atomLocality, *deviceStream);
399     // Note: unlike copyCoordinatesToGpu this is not used in OpenCL, and the conditional is not needed.
400     xReadyOnHost_[atomLocality].markEvent(*deviceStream);
401
402     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
403     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
404 }
405
406 void StatePropagatorDataGpu::Impl::waitCoordinatesReadyOnHost(AtomLocality atomLocality)
407 {
408     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
409     xReadyOnHost_[atomLocality].waitForEvent();
410     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
411 }
412
413
414 DeviceBuffer<RVec> StatePropagatorDataGpu::Impl::getVelocities()
415 {
416     return d_v_;
417 }
418
419 void StatePropagatorDataGpu::Impl::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec> h_v,
420                                                        AtomLocality atomLocality)
421 {
422     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
423     const DeviceStream* deviceStream = vCopyStreams_[atomLocality];
424     GMX_ASSERT(deviceStream != nullptr,
425                "No stream is valid for copying velocities with given atom locality.");
426
427     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
428     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
429
430     copyToDevice(d_v_, h_v, d_vSize_, atomLocality, *deviceStream);
431     vReadyOnDevice_[atomLocality].markEvent(*deviceStream);
432
433     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
434     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
435 }
436
437 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality)
438 {
439     return &vReadyOnDevice_[atomLocality];
440 }
441
442
443 void StatePropagatorDataGpu::Impl::copyVelocitiesFromGpu(gmx::ArrayRef<gmx::RVec> h_v, AtomLocality atomLocality)
444 {
445     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
446     const DeviceStream* deviceStream = vCopyStreams_[atomLocality];
447     GMX_ASSERT(deviceStream != nullptr,
448                "No stream is valid for copying velocities with given atom locality.");
449
450     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
451     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
452
453     copyFromDevice(h_v, d_v_, d_vSize_, atomLocality, *deviceStream);
454     vReadyOnHost_[atomLocality].markEvent(*deviceStream);
455
456     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
457     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
458 }
459
460 void StatePropagatorDataGpu::Impl::waitVelocitiesReadyOnHost(AtomLocality atomLocality)
461 {
462     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
463     vReadyOnHost_[atomLocality].waitForEvent();
464     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
465 }
466
467
468 DeviceBuffer<RVec> StatePropagatorDataGpu::Impl::getForces()
469 {
470     return d_f_;
471 }
472
473 void StatePropagatorDataGpu::Impl::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec> h_f,
474                                                    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     copyToDevice(d_f_, h_f, d_fSize_, atomLocality, *deviceStream);
485     fReadyOnDevice_[atomLocality].markEvent(*deviceStream);
486
487     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
488     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
489 }
490
491 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::getForcesReadyOnDeviceEvent(AtomLocality atomLocality,
492                                                                                 bool useGpuFBufferOps)
493 {
494     if ((atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal) && useGpuFBufferOps)
495     {
496         return &fReducedOnDevice_;
497     }
498     else
499     {
500         return &fReadyOnDevice_[atomLocality];
501     }
502 }
503
504 GpuEventSynchronizer* StatePropagatorDataGpu::Impl::fReducedOnDevice()
505 {
506     return &fReducedOnDevice_;
507 }
508
509 void StatePropagatorDataGpu::Impl::copyForcesFromGpu(gmx::ArrayRef<gmx::RVec> h_f, AtomLocality atomLocality)
510 {
511     GMX_ASSERT(atomLocality < AtomLocality::Count, "Wrong atom locality.");
512     const DeviceStream* deviceStream = fCopyStreams_[atomLocality];
513     GMX_ASSERT(deviceStream != nullptr,
514                "No stream is valid for copying forces with given atom locality.");
515
516     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
517     wallcycle_sub_start(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
518
519     copyFromDevice(h_f, d_f_, d_fSize_, atomLocality, *deviceStream);
520     fReadyOnHost_[atomLocality].markEvent(*deviceStream);
521
522     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_STATE_PROPAGATOR_DATA);
523     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
524 }
525
526 void StatePropagatorDataGpu::Impl::waitForcesReadyOnHost(AtomLocality atomLocality)
527 {
528     wallcycle_start(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
529     fReadyOnHost_[atomLocality].waitForEvent();
530     wallcycle_stop(wcycle_, ewcWAIT_GPU_STATE_PROPAGATOR_DATA);
531 }
532
533 const DeviceStream* StatePropagatorDataGpu::Impl::getUpdateStream()
534 {
535     return updateStream_;
536 }
537
538 int StatePropagatorDataGpu::Impl::numAtomsLocal()
539 {
540     return numAtomsLocal_;
541 }
542
543 int StatePropagatorDataGpu::Impl::numAtomsAll()
544 {
545     return numAtomsAll_;
546 }
547
548
549 StatePropagatorDataGpu::StatePropagatorDataGpu(const DeviceStream*  pmeStream,
550                                                const DeviceStream*  localStream,
551                                                const DeviceStream*  nonLocalStream,
552                                                const DeviceContext& deviceContext,
553                                                GpuApiCallBehavior   transferKind,
554                                                int                  allocationBlockSizeDivisor,
555                                                gmx_wallcycle*       wcycle) :
556     impl_(new Impl(pmeStream, localStream, nonLocalStream, deviceContext, transferKind, allocationBlockSizeDivisor, wcycle))
557 {
558 }
559
560 StatePropagatorDataGpu::StatePropagatorDataGpu(const DeviceStream*  pmeStream,
561                                                const DeviceContext& deviceContext,
562                                                GpuApiCallBehavior   transferKind,
563                                                int                  allocationBlockSizeDivisor,
564                                                gmx_wallcycle*       wcycle) :
565     impl_(new Impl(pmeStream, deviceContext, transferKind, allocationBlockSizeDivisor, wcycle))
566 {
567 }
568
569 StatePropagatorDataGpu::StatePropagatorDataGpu(StatePropagatorDataGpu&& /* other */) noexcept = default;
570
571 StatePropagatorDataGpu& StatePropagatorDataGpu::operator=(StatePropagatorDataGpu&& /* other */) noexcept = default;
572
573 StatePropagatorDataGpu::~StatePropagatorDataGpu() = default;
574
575
576 void StatePropagatorDataGpu::reinit(int numAtomsLocal, int numAtomsAll)
577 {
578     return impl_->reinit(numAtomsLocal, numAtomsAll);
579 }
580
581 std::tuple<int, int> StatePropagatorDataGpu::getAtomRangesFromAtomLocality(AtomLocality atomLocality)
582 {
583     return impl_->getAtomRangesFromAtomLocality(atomLocality);
584 }
585
586
587 DeviceBuffer<RVec> StatePropagatorDataGpu::getCoordinates()
588 {
589     return impl_->getCoordinates();
590 }
591
592 void StatePropagatorDataGpu::copyCoordinatesToGpu(const gmx::ArrayRef<const gmx::RVec> h_x,
593                                                   AtomLocality                         atomLocality)
594 {
595     return impl_->copyCoordinatesToGpu(h_x, atomLocality);
596 }
597
598 GpuEventSynchronizer*
599 StatePropagatorDataGpu::getCoordinatesReadyOnDeviceEvent(AtomLocality              atomLocality,
600                                                          const SimulationWorkload& simulationWork,
601                                                          const StepWorkload&       stepWork)
602 {
603     return impl_->getCoordinatesReadyOnDeviceEvent(atomLocality, simulationWork, stepWork);
604 }
605
606 void StatePropagatorDataGpu::waitCoordinatesCopiedToDevice(AtomLocality atomLocality)
607 {
608     return impl_->waitCoordinatesCopiedToDevice(atomLocality);
609 }
610
611 GpuEventSynchronizer* StatePropagatorDataGpu::xUpdatedOnDevice()
612 {
613     return impl_->xUpdatedOnDevice();
614 }
615
616 void StatePropagatorDataGpu::copyCoordinatesFromGpu(gmx::ArrayRef<RVec> h_x, AtomLocality atomLocality)
617 {
618     return impl_->copyCoordinatesFromGpu(h_x, atomLocality);
619 }
620
621 void StatePropagatorDataGpu::waitCoordinatesReadyOnHost(AtomLocality atomLocality)
622 {
623     return impl_->waitCoordinatesReadyOnHost(atomLocality);
624 }
625
626
627 DeviceBuffer<RVec> StatePropagatorDataGpu::getVelocities()
628 {
629     return impl_->getVelocities();
630 }
631
632 void StatePropagatorDataGpu::copyVelocitiesToGpu(const gmx::ArrayRef<const gmx::RVec> h_v,
633                                                  AtomLocality                         atomLocality)
634 {
635     return impl_->copyVelocitiesToGpu(h_v, atomLocality);
636 }
637
638 GpuEventSynchronizer* StatePropagatorDataGpu::getVelocitiesReadyOnDeviceEvent(AtomLocality atomLocality)
639 {
640     return impl_->getVelocitiesReadyOnDeviceEvent(atomLocality);
641 }
642
643 void StatePropagatorDataGpu::copyVelocitiesFromGpu(gmx::ArrayRef<RVec> h_v, AtomLocality atomLocality)
644 {
645     return impl_->copyVelocitiesFromGpu(h_v, atomLocality);
646 }
647
648 void StatePropagatorDataGpu::waitVelocitiesReadyOnHost(AtomLocality atomLocality)
649 {
650     return impl_->waitVelocitiesReadyOnHost(atomLocality);
651 }
652
653
654 DeviceBuffer<RVec> StatePropagatorDataGpu::getForces()
655 {
656     return impl_->getForces();
657 }
658
659 void StatePropagatorDataGpu::copyForcesToGpu(const gmx::ArrayRef<const gmx::RVec> h_f, AtomLocality atomLocality)
660 {
661     return impl_->copyForcesToGpu(h_f, atomLocality);
662 }
663
664 GpuEventSynchronizer* StatePropagatorDataGpu::getForcesReadyOnDeviceEvent(AtomLocality atomLocality,
665                                                                           bool useGpuFBufferOps)
666 {
667     return impl_->getForcesReadyOnDeviceEvent(atomLocality, useGpuFBufferOps);
668 }
669
670 GpuEventSynchronizer* StatePropagatorDataGpu::fReducedOnDevice()
671 {
672     return impl_->fReducedOnDevice();
673 }
674
675 void StatePropagatorDataGpu::copyForcesFromGpu(gmx::ArrayRef<RVec> h_f, AtomLocality atomLocality)
676 {
677     return impl_->copyForcesFromGpu(h_f, atomLocality);
678 }
679
680 void StatePropagatorDataGpu::waitForcesReadyOnHost(AtomLocality atomLocality)
681 {
682     return impl_->waitForcesReadyOnHost(atomLocality);
683 }
684
685
686 const DeviceStream* StatePropagatorDataGpu::getUpdateStream()
687 {
688     return impl_->getUpdateStream();
689 }
690
691 int StatePropagatorDataGpu::numAtomsLocal()
692 {
693     return impl_->numAtomsLocal();
694 }
695
696 int StatePropagatorDataGpu::numAtomsAll()
697 {
698     return impl_->numAtomsAll();
699 }
700
701 } // namespace gmx
702
703 #endif // GMX_GPU == GMX_GPU_NONE