b1d6774a26eb329362fc11cee409be22644c11a8
[alexxy/gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_cuda_data_mgmt.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
9  *
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.
14  *
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.
19  *
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.
24  *
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.
32  *
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.
35  */
36 /*! \file
37  *  \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
38  *
39  *  \author Szilard Pall <pall.szilard@gmail.com>
40  */
41 #include "gmxpre.h"
42
43 #include <assert.h>
44 #include <stdarg.h>
45 #include <stdio.h>
46 #include <stdlib.h>
47
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"
51
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"
79
80 #include "nbnxm_cuda.h"
81
82 namespace Nbnxm
83 {
84
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.
92  */
93 static unsigned int gpu_min_ci_balanced_factor = 44;
94
95 /* Fw. decl. */
96 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
97
98 /*! Initializes the atomdata structure first time, it only gets filled at
99     pair-search. */
100 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes, const DeviceContext& deviceContext)
101 {
102     ad->ntypes = ntypes;
103     allocateDeviceBuffer(&ad->shift_vec, SHIFTS, deviceContext);
104     ad->bShiftVecUploaded = false;
105
106     allocateDeviceBuffer(&ad->fshift, SHIFTS, deviceContext);
107     allocateDeviceBuffer(&ad->e_lj, 1, deviceContext);
108     allocateDeviceBuffer(&ad->e_el, 1, deviceContext);
109
110     /* initialize to nullptr poiters to data that is not allocated here and will
111        need reallocation in nbnxn_cuda_init_atomdata */
112     ad->xq = nullptr;
113     ad->f  = nullptr;
114
115     /* size -1 indicates that the respective array hasn't been initialized yet */
116     ad->natoms = -1;
117     ad->nalloc = -1;
118 }
119
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)
126 {
127     int ntypes;
128
129     ntypes = nbatParams.numTypes;
130
131     set_cutoff_parameters(nbp, ic, listParams);
132
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.
141      */
142     if (ic->vdwtype == evdwCUT)
143     {
144         switch (ic->vdw_modifier)
145         {
146             case eintmodNONE:
147             case eintmodPOTSHIFT:
148                 switch (nbatParams.comb_rule)
149                 {
150                     case ljcrNONE: nbp->vdwtype = evdwTypeCUT; break;
151                     case ljcrGEOM: nbp->vdwtype = evdwTypeCUTCOMBGEOM; break;
152                     case ljcrLB: nbp->vdwtype = evdwTypeCUTCOMBLB; break;
153                     default:
154                         gmx_incons(
155                                 "The requested LJ combination rule is not implemented in the CUDA "
156                                 "GPU accelerated kernels!");
157                 }
158                 break;
159             case eintmodFORCESWITCH: nbp->vdwtype = evdwTypeFSWITCH; break;
160             case eintmodPOTSWITCH: nbp->vdwtype = evdwTypePSWITCH; break;
161             default:
162                 gmx_incons(
163                         "The requested VdW interaction modifier is not implemented in the CUDA GPU "
164                         "accelerated kernels!");
165         }
166     }
167     else if (ic->vdwtype == evdwPME)
168     {
169         if (ic->ljpme_comb_rule == ljcrGEOM)
170         {
171             assert(nbatParams.comb_rule == ljcrGEOM);
172             nbp->vdwtype = evdwTypeEWALDGEOM;
173         }
174         else
175         {
176             assert(nbatParams.comb_rule == ljcrLB);
177             nbp->vdwtype = evdwTypeEWALDLB;
178         }
179     }
180     else
181     {
182         gmx_incons(
183                 "The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
184     }
185
186     if (ic->eeltype == eelCUT)
187     {
188         nbp->eeltype = eelTypeCUT;
189     }
190     else if (EEL_RF(ic->eeltype))
191     {
192         nbp->eeltype = eelTypeRF;
193     }
194     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
195     {
196         nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
197     }
198     else
199     {
200         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
201         gmx_incons(
202                 "The requested electrostatics type is not implemented in the CUDA GPU accelerated "
203                 "kernels!");
204     }
205
206     /* generate table for PME */
207     nbp->coulomb_tab = nullptr;
208     if (nbp->eeltype == eelTypeEWALD_TAB || nbp->eeltype == eelTypeEWALD_TAB_TWIN)
209     {
210         GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
211         init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, deviceContext);
212     }
213
214     /* set up LJ parameter lookup table */
215     if (!useLjCombRule(nbp->vdwtype))
216     {
217         initParamLookupTable(&nbp->nbfp, &nbp->nbfp_texobj, nbatParams.nbfp.data(),
218                              2 * ntypes * ntypes, deviceContext);
219     }
220
221     /* set up LJ-PME parameter lookup table */
222     if (ic->vdwtype == evdwPME)
223     {
224         initParamLookupTable(&nbp->nbfp_comb, &nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(),
225                              2 * ntypes, deviceContext);
226     }
227 }
228
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)
234 {
235     init_atomdata_first(nb->atdat, nbatParams.numTypes, *nb->deviceContext_);
236     init_nbparam(nb->nbparam, ic, listParams, nbatParams, *nb->deviceContext_);
237
238     /* clear energy and shift force outputs */
239     nbnxn_cuda_clear_e_fshift(nb);
240 }
241
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)
247 {
248     cudaError_t stat;
249
250     auto nb            = new NbnxmGpu();
251     nb->deviceContext_ = &deviceStreamManager.context();
252     snew(nb->atdat, 1);
253     snew(nb->nbparam, 1);
254     snew(nb->plist[InteractionLocality::Local], 1);
255     if (bLocalAndNonlocal)
256     {
257         snew(nb->plist[InteractionLocality::NonLocal], 1);
258     }
259
260     nb->bUseTwoStreams = bLocalAndNonlocal;
261
262     nb->timers = new cu_timers_t();
263     snew(nb->timings, 1);
264
265     /* init nbst */
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));
269
270     init_plist(nb->plist[InteractionLocality::Local]);
271
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)
278     {
279         init_plist(nb->plist[InteractionLocality::NonLocal]);
280
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.
284          */
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);
290         ;
291     }
292
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");
298
299     nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
300
301     /* WARNING: CUDA timings are incorrect with multiple streams.
302      *          This is the main reason why they are disabled by default.
303      */
304     // TODO: Consider turning on by default when we can detect nr of streams.
305     nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
306
307     if (nb->bDoTime)
308     {
309         init_timings(nb->timings);
310     }
311
312     /* set the kernel type for the current GPU */
313     /* pick L1 cache configuration */
314     cuda_set_cacheconfig();
315
316     cuda_init_const(nb, ic, listParams, nbat->params());
317
318     nb->atomIndicesSize       = 0;
319     nb->atomIndicesSize_alloc = 0;
320     nb->ncxy_na               = 0;
321     nb->ncxy_na_alloc         = 0;
322     nb->ncxy_ind              = 0;
323     nb->ncxy_ind_alloc        = 0;
324
325     if (debug)
326     {
327         fprintf(debug, "Initialized CUDA data structures.\n");
328     }
329
330     return nb;
331 }
332
333 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
334 {
335     cu_atomdata_t*      adat        = nb->atdat;
336     const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
337
338     /* only if we have a dynamic box */
339     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
340     {
341         static_assert(sizeof(adat->shift_vec[0]) == sizeof(nbatom->shift_vec[0]),
342                       "Sizes of host- and device-side shift vectors should be the same.");
343         copyToDeviceBuffer(&adat->shift_vec, reinterpret_cast<const float3*>(nbatom->shift_vec.data()),
344                            0, SHIFTS, localStream, GpuApiCallBehavior::Async, nullptr);
345         adat->bShiftVecUploaded = true;
346     }
347 }
348
349 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
350 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
351 {
352     cu_atomdata_t*      adat        = nb->atdat;
353     const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
354     clearDeviceBufferAsync(&adat->f, 0, natoms_clear, localStream);
355 }
356
357 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
358 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
359 {
360     cu_atomdata_t*      adat        = nb->atdat;
361     const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
362
363     clearDeviceBufferAsync(&adat->fshift, 0, SHIFTS, localStream);
364     clearDeviceBufferAsync(&adat->e_lj, 0, 1, localStream);
365     clearDeviceBufferAsync(&adat->e_el, 0, 1, localStream);
366 }
367
368 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
369 {
370     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
371     /* clear shift force array and energies if the outputs were
372        used in the current step */
373     if (computeVirial)
374     {
375         nbnxn_cuda_clear_e_fshift(nb);
376     }
377 }
378
379 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
380 {
381     int                  nalloc, natoms;
382     bool                 realloced;
383     bool                 bDoTime       = nb->bDoTime;
384     cu_timers_t*         timers        = nb->timers;
385     cu_atomdata_t*       d_atdat       = nb->atdat;
386     const DeviceContext& deviceContext = *nb->deviceContext_;
387     const DeviceStream&  localStream   = *nb->deviceStreams[InteractionLocality::Local];
388
389     natoms    = nbat->numAtoms();
390     realloced = false;
391
392     if (bDoTime)
393     {
394         /* time async copy */
395         timers->atdat.openTimingRegion(localStream);
396     }
397
398     /* need to reallocate if we have to copy more atoms than the amount of space
399        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
400     if (natoms > d_atdat->nalloc)
401     {
402         nalloc = over_alloc_small(natoms);
403
404         /* free up first if the arrays have already been initialized */
405         if (d_atdat->nalloc != -1)
406         {
407             freeDeviceBuffer(&d_atdat->f);
408             freeDeviceBuffer(&d_atdat->xq);
409             freeDeviceBuffer(&d_atdat->atom_types);
410             freeDeviceBuffer(&d_atdat->lj_comb);
411         }
412
413         allocateDeviceBuffer(&d_atdat->f, nalloc, deviceContext);
414         allocateDeviceBuffer(&d_atdat->xq, nalloc, deviceContext);
415         if (useLjCombRule(nb->nbparam->vdwtype))
416         {
417             allocateDeviceBuffer(&d_atdat->lj_comb, nalloc, deviceContext);
418         }
419         else
420         {
421             allocateDeviceBuffer(&d_atdat->atom_types, nalloc, deviceContext);
422         }
423
424         d_atdat->nalloc = nalloc;
425         realloced       = true;
426     }
427
428     d_atdat->natoms       = natoms;
429     d_atdat->natoms_local = nbat->natoms_local;
430
431     /* need to clear GPU f output if realloc happened */
432     if (realloced)
433     {
434         nbnxn_cuda_clear_f(nb, nalloc);
435     }
436
437     if (useLjCombRule(nb->nbparam->vdwtype))
438     {
439         static_assert(sizeof(d_atdat->lj_comb[0]) == sizeof(float2),
440                       "Size of the LJ parameters element should be equal to the size of float2.");
441         copyToDeviceBuffer(&d_atdat->lj_comb,
442                            reinterpret_cast<const float2*>(nbat->params().lj_comb.data()), 0,
443                            natoms, localStream, GpuApiCallBehavior::Async, nullptr);
444     }
445     else
446     {
447         static_assert(sizeof(d_atdat->atom_types[0]) == sizeof(nbat->params().type[0]),
448                       "Sizes of host- and device-side atom types should be the same.");
449         copyToDeviceBuffer(&d_atdat->atom_types, nbat->params().type.data(), 0, natoms, localStream,
450                            GpuApiCallBehavior::Async, nullptr);
451     }
452
453     if (bDoTime)
454     {
455         timers->atdat.closeTimingRegion(localStream);
456     }
457 }
458
459 void gpu_free(NbnxmGpu* nb)
460 {
461     cudaError_t    stat;
462     cu_atomdata_t* atdat;
463     NBParamGpu*    nbparam;
464
465     if (nb == nullptr)
466     {
467         return;
468     }
469
470     atdat   = nb->atdat;
471     nbparam = nb->nbparam;
472
473     if ((!nbparam->coulomb_tab)
474         && (nbparam->eeltype == eelTypeEWALD_TAB || nbparam->eeltype == eelTypeEWALD_TAB_TWIN))
475     {
476         destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
477     }
478
479     stat = cudaEventDestroy(nb->nonlocal_done);
480     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
481     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
482     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
483
484     delete nb->timers;
485
486     if (!useLjCombRule(nb->nbparam->vdwtype))
487     {
488         destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
489     }
490
491     if (nbparam->vdwtype == evdwTypeEWALDGEOM || nbparam->vdwtype == evdwTypeEWALDLB)
492     {
493         destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
494     }
495
496     freeDeviceBuffer(&atdat->shift_vec);
497     freeDeviceBuffer(&atdat->fshift);
498
499     freeDeviceBuffer(&atdat->e_lj);
500     freeDeviceBuffer(&atdat->e_el);
501
502     freeDeviceBuffer(&atdat->f);
503     freeDeviceBuffer(&atdat->xq);
504     freeDeviceBuffer(&atdat->atom_types);
505     freeDeviceBuffer(&atdat->lj_comb);
506
507     /* Free plist */
508     auto* plist = nb->plist[InteractionLocality::Local];
509     freeDeviceBuffer(&plist->sci);
510     freeDeviceBuffer(&plist->cj4);
511     freeDeviceBuffer(&plist->imask);
512     freeDeviceBuffer(&plist->excl);
513     sfree(plist);
514     if (nb->bUseTwoStreams)
515     {
516         auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
517         freeDeviceBuffer(&plist_nl->sci);
518         freeDeviceBuffer(&plist_nl->cj4);
519         freeDeviceBuffer(&plist_nl->imask);
520         freeDeviceBuffer(&plist_nl->excl);
521         sfree(plist_nl);
522     }
523
524     /* Free nbst */
525     pfree(nb->nbst.e_lj);
526     nb->nbst.e_lj = nullptr;
527
528     pfree(nb->nbst.e_el);
529     nb->nbst.e_el = nullptr;
530
531     pfree(nb->nbst.fshift);
532     nb->nbst.fshift = nullptr;
533
534     sfree(atdat);
535     sfree(nbparam);
536     sfree(nb->timings);
537     delete nb;
538
539     if (debug)
540     {
541         fprintf(debug, "Cleaned up CUDA data structures.\n");
542     }
543 }
544
545 int gpu_min_ci_balanced(NbnxmGpu* nb)
546 {
547     return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
548                          : 0;
549 }
550
551 void* gpu_get_xq(NbnxmGpu* nb)
552 {
553     assert(nb);
554
555     return static_cast<void*>(nb->atdat->xq);
556 }
557
558 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
559 {
560     assert(nb);
561
562     return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
563 }
564
565 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
566 {
567     assert(nb);
568
569     return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fshift);
570 }
571
572 /* Initialization for X buffer operations on GPU. */
573 /* TODO  Remove explicit pinning from host arrays from here and manage in a more natural way*/
574 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
575 {
576     const DeviceStream& deviceStream  = *gpu_nbv->deviceStreams[InteractionLocality::Local];
577     bool                bDoTime       = gpu_nbv->bDoTime;
578     const int           maxNumColumns = gridSet.numColumnsMax();
579
580     reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
581                            &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, *gpu_nbv->deviceContext_);
582     reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
583                            &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, *gpu_nbv->deviceContext_);
584
585     for (unsigned int g = 0; g < gridSet.grids().size(); g++)
586     {
587
588         const Nbnxm::Grid& grid = gridSet.grids()[g];
589
590         const int  numColumns      = grid.numColumns();
591         const int* atomIndices     = gridSet.atomIndices().data();
592         const int  atomIndicesSize = gridSet.atomIndices().size();
593         const int* cxy_na          = grid.cxy_na().data();
594         const int* cxy_ind         = grid.cxy_ind().data();
595
596         reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
597                                &gpu_nbv->atomIndicesSize_alloc, *gpu_nbv->deviceContext_);
598
599         if (atomIndicesSize > 0)
600         {
601
602             if (bDoTime)
603             {
604                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
605             }
606
607             copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, deviceStream,
608                                GpuApiCallBehavior::Async, nullptr);
609
610             if (bDoTime)
611             {
612                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
613             }
614         }
615
616         if (numColumns > 0)
617         {
618             if (bDoTime)
619             {
620                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
621             }
622
623             int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
624             copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, deviceStream,
625                                GpuApiCallBehavior::Async, nullptr);
626
627             if (bDoTime)
628             {
629                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
630             }
631
632             if (bDoTime)
633             {
634                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
635             }
636
637             destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
638             copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, deviceStream,
639                                GpuApiCallBehavior::Async, nullptr);
640
641             if (bDoTime)
642             {
643                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
644             }
645         }
646     }
647
648     // The above data is transferred on the local stream but is a
649     // dependency of the nonlocal stream (specifically the nonlocal X
650     // buf ops kernel).  We therefore set a dependency to ensure
651     // that the nonlocal stream waits on the local stream here.
652     // This call records an event in the local stream:
653     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
654     // ...and this call instructs the nonlocal stream to wait on that event:
655     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
656
657     return;
658 }
659
660 } // namespace Nbnxm