2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
39 * \author Szilard Pall <pall.szilard@gmail.com>
48 // TODO We would like to move this down, but the way NbnxmGpu
49 // is currently declared means this has to be before gpu_types.h
50 #include "nbnxm_cuda_types.h"
52 // TODO Remove this comment when the above order issue is resolved
53 #include "gromacs/gpu_utils/cudautils.cuh"
54 #include "gromacs/gpu_utils/device_context.h"
55 #include "gromacs/gpu_utils/device_stream_manager.h"
56 #include "gromacs/gpu_utils/gpu_utils.h"
57 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
58 #include "gromacs/gpu_utils/pmalloc_cuda.h"
59 #include "gromacs/hardware/device_information.h"
60 #include "gromacs/hardware/device_management.h"
61 #include "gromacs/math/vectypes.h"
62 #include "gromacs/mdlib/force_flags.h"
63 #include "gromacs/mdtypes/interaction_const.h"
64 #include "gromacs/mdtypes/md_enums.h"
65 #include "gromacs/nbnxm/atomdata.h"
66 #include "gromacs/nbnxm/gpu_data_mgmt.h"
67 #include "gromacs/nbnxm/gridset.h"
68 #include "gromacs/nbnxm/nbnxm.h"
69 #include "gromacs/nbnxm/nbnxm_gpu.h"
70 #include "gromacs/nbnxm/nbnxm_gpu_data_mgmt.h"
71 #include "gromacs/nbnxm/pairlistsets.h"
72 #include "gromacs/pbcutil/ishift.h"
73 #include "gromacs/timing/gpu_timing.h"
74 #include "gromacs/utility/basedefinitions.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/real.h"
78 #include "gromacs/utility/smalloc.h"
80 #include "nbnxm_cuda.h"
85 /* This is a heuristically determined parameter for the Kepler
86 * and Maxwell architectures for the minimum size of ci lists by multiplying
87 * this constant with the # of multiprocessors on the current device.
88 * Since the maximum number of blocks per multiprocessor is 16, the ideal
89 * count for small systems is 32 or 48 blocks per multiprocessor. Because
90 * there is a bit of fluctuations in the generated block counts, we use
91 * a target of 44 instead of the ideal value of 48.
93 static unsigned int gpu_min_ci_balanced_factor = 44;
96 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
98 /*! Initializes the atomdata structure first time, it only gets filled at
100 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes, const DeviceContext& deviceContext)
103 allocateDeviceBuffer(&ad->shift_vec, SHIFTS, deviceContext);
104 ad->bShiftVecUploaded = false;
106 allocateDeviceBuffer(&ad->fshift, SHIFTS, deviceContext);
107 allocateDeviceBuffer(&ad->e_lj, 1, deviceContext);
108 allocateDeviceBuffer(&ad->e_el, 1, deviceContext);
110 /* initialize to nullptr poiters to data that is not allocated here and will
111 need reallocation in nbnxn_cuda_init_atomdata */
115 /* size -1 indicates that the respective array hasn't been initialized yet */
120 /*! Initializes the nonbonded parameter data structure. */
121 static void init_nbparam(NBParamGpu* nbp,
122 const interaction_const_t* ic,
123 const PairlistParams& listParams,
124 const nbnxn_atomdata_t::Params& nbatParams,
125 const DeviceContext& deviceContext)
129 ntypes = nbatParams.numTypes;
131 set_cutoff_parameters(nbp, ic, listParams);
133 /* The kernel code supports LJ combination rules (geometric and LB) for
134 * all kernel types, but we only generate useful combination rule kernels.
135 * We currently only use LJ combination rule (geometric and LB) kernels
136 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
137 * with PME and 20% with RF, the other kernels speed up about half as much.
138 * For LJ force-switch the geometric rule would give 7% speed-up, but this
139 * combination is rarely used. LJ force-switch with LB rule is more common,
140 * but gives only 1% speed-up.
142 if (ic->vdwtype == evdwCUT)
144 switch (ic->vdw_modifier)
147 case eintmodPOTSHIFT:
148 switch (nbatParams.comb_rule)
150 case ljcrNONE: nbp->vdwtype = evdwTypeCUT; break;
151 case ljcrGEOM: nbp->vdwtype = evdwTypeCUTCOMBGEOM; break;
152 case ljcrLB: nbp->vdwtype = evdwTypeCUTCOMBLB; break;
155 "The requested LJ combination rule is not implemented in the CUDA "
156 "GPU accelerated kernels!");
159 case eintmodFORCESWITCH: nbp->vdwtype = evdwTypeFSWITCH; break;
160 case eintmodPOTSWITCH: nbp->vdwtype = evdwTypePSWITCH; break;
163 "The requested VdW interaction modifier is not implemented in the CUDA GPU "
164 "accelerated kernels!");
167 else if (ic->vdwtype == evdwPME)
169 if (ic->ljpme_comb_rule == ljcrGEOM)
171 assert(nbatParams.comb_rule == ljcrGEOM);
172 nbp->vdwtype = evdwTypeEWALDGEOM;
176 assert(nbatParams.comb_rule == ljcrLB);
177 nbp->vdwtype = evdwTypeEWALDLB;
183 "The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
186 if (ic->eeltype == eelCUT)
188 nbp->eeltype = eelTypeCUT;
190 else if (EEL_RF(ic->eeltype))
192 nbp->eeltype = eelTypeRF;
194 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
196 nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
200 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
202 "The requested electrostatics type is not implemented in the CUDA GPU accelerated "
206 /* generate table for PME */
207 nbp->coulomb_tab = nullptr;
208 if (nbp->eeltype == eelTypeEWALD_TAB || nbp->eeltype == eelTypeEWALD_TAB_TWIN)
210 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
211 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, deviceContext);
214 /* set up LJ parameter lookup table */
215 if (!useLjCombRule(nbp->vdwtype))
217 initParamLookupTable(&nbp->nbfp, &nbp->nbfp_texobj, nbatParams.nbfp.data(),
218 2 * ntypes * ntypes, deviceContext);
221 /* set up LJ-PME parameter lookup table */
222 if (ic->vdwtype == evdwPME)
224 initParamLookupTable(&nbp->nbfp_comb, &nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(),
225 2 * ntypes, deviceContext);
229 /*! Initializes simulation constant data. */
230 static void cuda_init_const(NbnxmGpu* nb,
231 const interaction_const_t* ic,
232 const PairlistParams& listParams,
233 const nbnxn_atomdata_t::Params& nbatParams)
235 init_atomdata_first(nb->atdat, nbatParams.numTypes, *nb->deviceContext_);
236 init_nbparam(nb->nbparam, ic, listParams, nbatParams, *nb->deviceContext_);
238 /* clear energy and shift force outputs */
239 nbnxn_cuda_clear_e_fshift(nb);
242 NbnxmGpu* gpu_init(const gmx::DeviceStreamManager& deviceStreamManager,
243 const interaction_const_t* ic,
244 const PairlistParams& listParams,
245 const nbnxn_atomdata_t* nbat,
246 bool bLocalAndNonlocal)
250 auto nb = new NbnxmGpu();
251 nb->deviceContext_ = &deviceStreamManager.context();
253 snew(nb->nbparam, 1);
254 snew(nb->plist[InteractionLocality::Local], 1);
255 if (bLocalAndNonlocal)
257 snew(nb->plist[InteractionLocality::NonLocal], 1);
260 nb->bUseTwoStreams = bLocalAndNonlocal;
262 nb->timers = new cu_timers_t();
263 snew(nb->timings, 1);
266 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
267 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
268 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
270 init_plist(nb->plist[InteractionLocality::Local]);
272 /* local/non-local GPU streams */
273 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedLocal),
274 "Local non-bonded stream should be initialized to use GPU for non-bonded.");
275 nb->deviceStreams[InteractionLocality::Local] =
276 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedLocal);
277 if (nb->bUseTwoStreams)
279 init_plist(nb->plist[InteractionLocality::NonLocal]);
281 /* Note that the device we're running on does not have to support
282 * priorities, because we are querying the priority range which in this
283 * case will be a single value.
285 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedNonLocal),
286 "Non-local non-bonded stream should be initialized to use GPU for "
287 "non-bonded with domain decomposition.");
288 nb->deviceStreams[InteractionLocality::NonLocal] =
289 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedNonLocal);
293 /* init events for sychronization (timing disabled for performance reasons!) */
294 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
295 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
296 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
297 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
299 nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
301 /* WARNING: CUDA timings are incorrect with multiple streams.
302 * This is the main reason why they are disabled by default.
304 // TODO: Consider turning on by default when we can detect nr of streams.
305 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
309 init_timings(nb->timings);
312 /* set the kernel type for the current GPU */
313 /* pick L1 cache configuration */
314 cuda_set_cacheconfig();
316 cuda_init_const(nb, ic, listParams, nbat->params());
318 nb->atomIndicesSize = 0;
319 nb->atomIndicesSize_alloc = 0;
321 nb->ncxy_na_alloc = 0;
323 nb->ncxy_ind_alloc = 0;
329 fprintf(debug, "Initialized CUDA data structures.\n");
335 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
337 cu_atomdata_t* adat = nb->atdat;
338 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
340 /* only if we have a dynamic box */
341 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
343 static_assert(sizeof(adat->shift_vec[0]) == sizeof(nbatom->shift_vec[0]),
344 "Sizes of host- and device-side shift vectors should be the same.");
345 copyToDeviceBuffer(&adat->shift_vec, reinterpret_cast<const float3*>(nbatom->shift_vec.data()),
346 0, SHIFTS, localStream, GpuApiCallBehavior::Async, nullptr);
347 adat->bShiftVecUploaded = true;
351 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
352 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
354 cu_atomdata_t* adat = nb->atdat;
355 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
356 clearDeviceBufferAsync(&adat->f, 0, natoms_clear, localStream);
359 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
360 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
362 cu_atomdata_t* adat = nb->atdat;
363 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
365 clearDeviceBufferAsync(&adat->fshift, 0, SHIFTS, localStream);
366 clearDeviceBufferAsync(&adat->e_lj, 0, 1, localStream);
367 clearDeviceBufferAsync(&adat->e_el, 0, 1, localStream);
370 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
372 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
373 /* clear shift force array and energies if the outputs were
374 used in the current step */
377 nbnxn_cuda_clear_e_fshift(nb);
381 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
385 bool bDoTime = nb->bDoTime;
386 cu_timers_t* timers = nb->timers;
387 cu_atomdata_t* d_atdat = nb->atdat;
388 const DeviceContext& deviceContext = *nb->deviceContext_;
389 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
391 natoms = nbat->numAtoms();
396 /* time async copy */
397 timers->atdat.openTimingRegion(localStream);
400 /* need to reallocate if we have to copy more atoms than the amount of space
401 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
402 if (natoms > d_atdat->nalloc)
404 nalloc = over_alloc_small(natoms);
406 /* free up first if the arrays have already been initialized */
407 if (d_atdat->nalloc != -1)
409 freeDeviceBuffer(&d_atdat->f);
410 freeDeviceBuffer(&d_atdat->xq);
411 freeDeviceBuffer(&d_atdat->atom_types);
412 freeDeviceBuffer(&d_atdat->lj_comb);
415 allocateDeviceBuffer(&d_atdat->f, nalloc, deviceContext);
416 allocateDeviceBuffer(&d_atdat->xq, nalloc, deviceContext);
417 if (useLjCombRule(nb->nbparam->vdwtype))
419 allocateDeviceBuffer(&d_atdat->lj_comb, nalloc, deviceContext);
423 allocateDeviceBuffer(&d_atdat->atom_types, nalloc, deviceContext);
426 d_atdat->nalloc = nalloc;
430 d_atdat->natoms = natoms;
431 d_atdat->natoms_local = nbat->natoms_local;
433 /* need to clear GPU f output if realloc happened */
436 nbnxn_cuda_clear_f(nb, nalloc);
439 if (useLjCombRule(nb->nbparam->vdwtype))
441 static_assert(sizeof(d_atdat->lj_comb[0]) == sizeof(float2),
442 "Size of the LJ parameters element should be equal to the size of float2.");
443 copyToDeviceBuffer(&d_atdat->lj_comb,
444 reinterpret_cast<const float2*>(nbat->params().lj_comb.data()), 0,
445 natoms, localStream, GpuApiCallBehavior::Async, nullptr);
449 static_assert(sizeof(d_atdat->atom_types[0]) == sizeof(nbat->params().type[0]),
450 "Sizes of host- and device-side atom types should be the same.");
451 copyToDeviceBuffer(&d_atdat->atom_types, nbat->params().type.data(), 0, natoms, localStream,
452 GpuApiCallBehavior::Async, nullptr);
457 timers->atdat.closeTimingRegion(localStream);
461 void gpu_free(NbnxmGpu* nb)
464 cu_atomdata_t* atdat;
473 nbparam = nb->nbparam;
475 if ((!nbparam->coulomb_tab)
476 && (nbparam->eeltype == eelTypeEWALD_TAB || nbparam->eeltype == eelTypeEWALD_TAB_TWIN))
478 destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
481 stat = cudaEventDestroy(nb->nonlocal_done);
482 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
483 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
484 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
488 if (!useLjCombRule(nb->nbparam->vdwtype))
490 destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
493 if (nbparam->vdwtype == evdwTypeEWALDGEOM || nbparam->vdwtype == evdwTypeEWALDLB)
495 destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
498 freeDeviceBuffer(&atdat->shift_vec);
499 freeDeviceBuffer(&atdat->fshift);
501 freeDeviceBuffer(&atdat->e_lj);
502 freeDeviceBuffer(&atdat->e_el);
504 freeDeviceBuffer(&atdat->f);
505 freeDeviceBuffer(&atdat->xq);
506 freeDeviceBuffer(&atdat->atom_types);
507 freeDeviceBuffer(&atdat->lj_comb);
510 auto* plist = nb->plist[InteractionLocality::Local];
511 freeDeviceBuffer(&plist->sci);
512 freeDeviceBuffer(&plist->cj4);
513 freeDeviceBuffer(&plist->imask);
514 freeDeviceBuffer(&plist->excl);
516 if (nb->bUseTwoStreams)
518 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
519 freeDeviceBuffer(&plist_nl->sci);
520 freeDeviceBuffer(&plist_nl->cj4);
521 freeDeviceBuffer(&plist_nl->imask);
522 freeDeviceBuffer(&plist_nl->excl);
527 pfree(nb->nbst.e_lj);
528 nb->nbst.e_lj = nullptr;
530 pfree(nb->nbst.e_el);
531 nb->nbst.e_el = nullptr;
533 pfree(nb->nbst.fshift);
534 nb->nbst.fshift = nullptr;
543 fprintf(debug, "Cleaned up CUDA data structures.\n");
547 int gpu_min_ci_balanced(NbnxmGpu* nb)
549 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
553 void* gpu_get_xq(NbnxmGpu* nb)
557 return static_cast<void*>(nb->atdat->xq);
560 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
564 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
567 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
571 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fshift);
574 /* Initialization for X buffer operations on GPU. */
575 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
576 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
578 const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
579 bool bDoTime = gpu_nbv->bDoTime;
580 const int maxNumColumns = gridSet.numColumnsMax();
582 reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
583 &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, *gpu_nbv->deviceContext_);
584 reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
585 &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, *gpu_nbv->deviceContext_);
587 for (unsigned int g = 0; g < gridSet.grids().size(); g++)
590 const Nbnxm::Grid& grid = gridSet.grids()[g];
592 const int numColumns = grid.numColumns();
593 const int* atomIndices = gridSet.atomIndices().data();
594 const int atomIndicesSize = gridSet.atomIndices().size();
595 const int* cxy_na = grid.cxy_na().data();
596 const int* cxy_ind = grid.cxy_ind().data();
598 reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
599 &gpu_nbv->atomIndicesSize_alloc, *gpu_nbv->deviceContext_);
601 if (atomIndicesSize > 0)
606 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
609 copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, deviceStream,
610 GpuApiCallBehavior::Async, nullptr);
614 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
622 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
625 int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
626 copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, deviceStream,
627 GpuApiCallBehavior::Async, nullptr);
631 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
636 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
639 destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
640 copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, deviceStream,
641 GpuApiCallBehavior::Async, nullptr);
645 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
650 // The above data is transferred on the local stream but is a
651 // dependency of the nonlocal stream (specifically the nonlocal X
652 // buf ops kernel). We therefore set a dependency to ensure
653 // that the nonlocal stream waits on the local stream here.
654 // This call records an event in the local stream:
655 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
656 // ...and this call instructs the nonlocal stream to wait on that event:
657 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
662 /* Initialization for F buffer operations on GPU. */
663 void nbnxn_gpu_init_add_nbat_f_to_f(const int* cell,
666 GpuEventSynchronizer* const localReductionDone)
669 const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
671 GMX_ASSERT(localReductionDone, "localReductionDone should be a valid pointer");
672 gpu_nbv->localFReductionDone = localReductionDone;
674 if (natoms_total > 0)
676 reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc,
677 *gpu_nbv->deviceContext_);
678 copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, deviceStream,
679 GpuApiCallBehavior::Async, nullptr);