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_stream_manager.h"
55 #include "gromacs/gpu_utils/gpu_utils.h"
56 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
57 #include "gromacs/gpu_utils/pmalloc_cuda.h"
58 #include "gromacs/hardware/device_information.h"
59 #include "gromacs/math/vectypes.h"
60 #include "gromacs/mdlib/force_flags.h"
61 #include "gromacs/mdtypes/interaction_const.h"
62 #include "gromacs/mdtypes/md_enums.h"
63 #include "gromacs/nbnxm/atomdata.h"
64 #include "gromacs/nbnxm/gpu_data_mgmt.h"
65 #include "gromacs/nbnxm/gridset.h"
66 #include "gromacs/nbnxm/nbnxm.h"
67 #include "gromacs/nbnxm/nbnxm_gpu.h"
68 #include "gromacs/nbnxm/nbnxm_gpu_data_mgmt.h"
69 #include "gromacs/nbnxm/pairlistsets.h"
70 #include "gromacs/pbcutil/ishift.h"
71 #include "gromacs/timing/gpu_timing.h"
72 #include "gromacs/utility/basedefinitions.h"
73 #include "gromacs/utility/cstringutil.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/real.h"
76 #include "gromacs/utility/smalloc.h"
78 #include "nbnxm_cuda.h"
83 /* This is a heuristically determined parameter for the Kepler
84 * and Maxwell architectures for the minimum size of ci lists by multiplying
85 * this constant with the # of multiprocessors on the current device.
86 * Since the maximum number of blocks per multiprocessor is 16, the ideal
87 * count for small systems is 32 or 48 blocks per multiprocessor. Because
88 * there is a bit of fluctuations in the generated block counts, we use
89 * a target of 44 instead of the ideal value of 48.
91 static unsigned int gpu_min_ci_balanced_factor = 44;
94 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
96 /*! Initializes the atomdata structure first time, it only gets filled at
98 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes, const DeviceContext& deviceContext)
101 allocateDeviceBuffer(&ad->shift_vec, SHIFTS, deviceContext);
102 ad->bShiftVecUploaded = false;
104 allocateDeviceBuffer(&ad->fshift, SHIFTS, deviceContext);
105 allocateDeviceBuffer(&ad->e_lj, 1, deviceContext);
106 allocateDeviceBuffer(&ad->e_el, 1, deviceContext);
108 /* initialize to nullptr poiters to data that is not allocated here and will
109 need reallocation in nbnxn_cuda_init_atomdata */
113 /* size -1 indicates that the respective array hasn't been initialized yet */
118 /*! Initializes the nonbonded parameter data structure. */
119 static void init_nbparam(NBParamGpu* nbp,
120 const interaction_const_t* ic,
121 const PairlistParams& listParams,
122 const nbnxn_atomdata_t::Params& nbatParams,
123 const DeviceContext& deviceContext)
127 ntypes = nbatParams.numTypes;
129 set_cutoff_parameters(nbp, ic, listParams);
131 /* The kernel code supports LJ combination rules (geometric and LB) for
132 * all kernel types, but we only generate useful combination rule kernels.
133 * We currently only use LJ combination rule (geometric and LB) kernels
134 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
135 * with PME and 20% with RF, the other kernels speed up about half as much.
136 * For LJ force-switch the geometric rule would give 7% speed-up, but this
137 * combination is rarely used. LJ force-switch with LB rule is more common,
138 * but gives only 1% speed-up.
140 if (ic->vdwtype == evdwCUT)
142 switch (ic->vdw_modifier)
145 case eintmodPOTSHIFT:
146 switch (nbatParams.comb_rule)
148 case ljcrNONE: nbp->vdwtype = evdwTypeCUT; break;
149 case ljcrGEOM: nbp->vdwtype = evdwTypeCUTCOMBGEOM; break;
150 case ljcrLB: nbp->vdwtype = evdwTypeCUTCOMBLB; break;
153 "The requested LJ combination rule is not implemented in the CUDA "
154 "GPU accelerated kernels!");
157 case eintmodFORCESWITCH: nbp->vdwtype = evdwTypeFSWITCH; break;
158 case eintmodPOTSWITCH: nbp->vdwtype = evdwTypePSWITCH; break;
161 "The requested VdW interaction modifier is not implemented in the CUDA GPU "
162 "accelerated kernels!");
165 else if (ic->vdwtype == evdwPME)
167 if (ic->ljpme_comb_rule == ljcrGEOM)
169 assert(nbatParams.comb_rule == ljcrGEOM);
170 nbp->vdwtype = evdwTypeEWALDGEOM;
174 assert(nbatParams.comb_rule == ljcrLB);
175 nbp->vdwtype = evdwTypeEWALDLB;
181 "The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
184 if (ic->eeltype == eelCUT)
186 nbp->eeltype = eelTypeCUT;
188 else if (EEL_RF(ic->eeltype))
190 nbp->eeltype = eelTypeRF;
192 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
194 nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
198 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
200 "The requested electrostatics type is not implemented in the CUDA GPU accelerated "
204 /* generate table for PME */
205 nbp->coulomb_tab = nullptr;
206 if (nbp->eeltype == eelTypeEWALD_TAB || nbp->eeltype == eelTypeEWALD_TAB_TWIN)
208 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
209 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, deviceContext);
212 /* set up LJ parameter lookup table */
213 if (!useLjCombRule(nbp->vdwtype))
215 initParamLookupTable(&nbp->nbfp, &nbp->nbfp_texobj, nbatParams.nbfp.data(),
216 2 * ntypes * ntypes, deviceContext);
219 /* set up LJ-PME parameter lookup table */
220 if (ic->vdwtype == evdwPME)
222 initParamLookupTable(&nbp->nbfp_comb, &nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(),
223 2 * ntypes, deviceContext);
227 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
228 * electrostatic type switching to twin cut-off (or back) if needed. */
229 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
231 if (!nbv || !nbv->useGpu())
235 NbnxmGpu* nb = nbv->gpu_nbv;
236 NBParamGpu* nbp = nbv->gpu_nbv->nbparam;
238 set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
240 nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
242 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
243 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, *nb->deviceContext_);
246 /*! Initializes simulation constant data. */
247 static void cuda_init_const(NbnxmGpu* nb,
248 const interaction_const_t* ic,
249 const PairlistParams& listParams,
250 const nbnxn_atomdata_t::Params& nbatParams)
252 init_atomdata_first(nb->atdat, nbatParams.numTypes, *nb->deviceContext_);
253 init_nbparam(nb->nbparam, ic, listParams, nbatParams, *nb->deviceContext_);
255 /* clear energy and shift force outputs */
256 nbnxn_cuda_clear_e_fshift(nb);
259 NbnxmGpu* gpu_init(const gmx::DeviceStreamManager& deviceStreamManager,
260 const interaction_const_t* ic,
261 const PairlistParams& listParams,
262 const nbnxn_atomdata_t* nbat,
263 bool bLocalAndNonlocal)
267 auto nb = new NbnxmGpu();
268 nb->deviceContext_ = &deviceStreamManager.context();
270 snew(nb->nbparam, 1);
271 snew(nb->plist[InteractionLocality::Local], 1);
272 if (bLocalAndNonlocal)
274 snew(nb->plist[InteractionLocality::NonLocal], 1);
277 nb->bUseTwoStreams = bLocalAndNonlocal;
279 nb->timers = new cu_timers_t();
280 snew(nb->timings, 1);
283 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
284 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
285 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
287 init_plist(nb->plist[InteractionLocality::Local]);
289 /* local/non-local GPU streams */
290 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedLocal),
291 "Local non-bonded stream should be initialized to use GPU for non-bonded.");
292 nb->deviceStreams[InteractionLocality::Local] =
293 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedLocal);
294 if (nb->bUseTwoStreams)
296 init_plist(nb->plist[InteractionLocality::NonLocal]);
298 /* Note that the device we're running on does not have to support
299 * priorities, because we are querying the priority range which in this
300 * case will be a single value.
302 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedNonLocal),
303 "Non-local non-bonded stream should be initialized to use GPU for "
304 "non-bonded with domain decomposition.");
305 nb->deviceStreams[InteractionLocality::NonLocal] =
306 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedNonLocal);
310 /* init events for sychronization (timing disabled for performance reasons!) */
311 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
312 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
313 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
314 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
316 nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
318 /* WARNING: CUDA timings are incorrect with multiple streams.
319 * This is the main reason why they are disabled by default.
321 // TODO: Consider turning on by default when we can detect nr of streams.
322 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
326 init_timings(nb->timings);
329 /* set the kernel type for the current GPU */
330 /* pick L1 cache configuration */
331 cuda_set_cacheconfig();
333 cuda_init_const(nb, ic, listParams, nbat->params());
335 nb->atomIndicesSize = 0;
336 nb->atomIndicesSize_alloc = 0;
338 nb->ncxy_na_alloc = 0;
340 nb->ncxy_ind_alloc = 0;
346 fprintf(debug, "Initialized CUDA data structures.\n");
352 void gpu_init_pairlist(NbnxmGpu* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
355 bool bDoTime = (nb->bDoTime && !h_plist->sci.empty());
356 const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
357 gpu_plist* d_plist = nb->plist[iloc];
359 if (d_plist->na_c < 0)
361 d_plist->na_c = h_plist->na_ci;
365 if (d_plist->na_c != h_plist->na_ci)
367 sprintf(sbuf, "In init_plist: the #atoms per cell has changed (from %d to %d)",
368 d_plist->na_c, h_plist->na_ci);
373 gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
377 iTimers.pl_h2d.openTimingRegion(deviceStream);
378 iTimers.didPairlistH2D = true;
381 const DeviceContext& deviceContext = *nb->deviceContext_;
383 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc,
385 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), deviceStream,
386 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
388 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc,
390 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), deviceStream,
391 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
393 reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
394 &d_plist->nimask, &d_plist->imask_nalloc, deviceContext);
396 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
397 &d_plist->excl_nalloc, deviceContext);
398 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), deviceStream,
399 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
403 iTimers.pl_h2d.closeTimingRegion(deviceStream);
406 /* the next use of thist list we be the first one, so we need to prune */
407 d_plist->haveFreshList = true;
410 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
412 cu_atomdata_t* adat = nb->atdat;
413 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
415 /* only if we have a dynamic box */
416 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
418 static_assert(sizeof(adat->shift_vec[0]) == sizeof(nbatom->shift_vec[0]),
419 "Sizes of host- and device-side shift vectors should be the same.");
420 copyToDeviceBuffer(&adat->shift_vec, reinterpret_cast<const float3*>(nbatom->shift_vec.data()),
421 0, SHIFTS, localStream, GpuApiCallBehavior::Async, nullptr);
422 adat->bShiftVecUploaded = true;
426 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
427 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
429 cu_atomdata_t* adat = nb->atdat;
430 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
431 clearDeviceBufferAsync(&adat->f, 0, natoms_clear, localStream);
434 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
435 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
437 cu_atomdata_t* adat = nb->atdat;
438 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
440 clearDeviceBufferAsync(&adat->fshift, 0, SHIFTS, localStream);
441 clearDeviceBufferAsync(&adat->e_lj, 0, 1, localStream);
442 clearDeviceBufferAsync(&adat->e_el, 0, 1, localStream);
445 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
447 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
448 /* clear shift force array and energies if the outputs were
449 used in the current step */
452 nbnxn_cuda_clear_e_fshift(nb);
456 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
460 bool bDoTime = nb->bDoTime;
461 cu_timers_t* timers = nb->timers;
462 cu_atomdata_t* d_atdat = nb->atdat;
463 const DeviceContext& deviceContext = *nb->deviceContext_;
464 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
466 natoms = nbat->numAtoms();
471 /* time async copy */
472 timers->atdat.openTimingRegion(localStream);
475 /* need to reallocate if we have to copy more atoms than the amount of space
476 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
477 if (natoms > d_atdat->nalloc)
479 nalloc = over_alloc_small(natoms);
481 /* free up first if the arrays have already been initialized */
482 if (d_atdat->nalloc != -1)
484 freeDeviceBuffer(&d_atdat->f);
485 freeDeviceBuffer(&d_atdat->xq);
486 freeDeviceBuffer(&d_atdat->atom_types);
487 freeDeviceBuffer(&d_atdat->lj_comb);
490 allocateDeviceBuffer(&d_atdat->f, nalloc, deviceContext);
491 allocateDeviceBuffer(&d_atdat->xq, nalloc, deviceContext);
492 if (useLjCombRule(nb->nbparam->vdwtype))
494 allocateDeviceBuffer(&d_atdat->lj_comb, nalloc, deviceContext);
498 allocateDeviceBuffer(&d_atdat->atom_types, nalloc, deviceContext);
501 d_atdat->nalloc = nalloc;
505 d_atdat->natoms = natoms;
506 d_atdat->natoms_local = nbat->natoms_local;
508 /* need to clear GPU f output if realloc happened */
511 nbnxn_cuda_clear_f(nb, nalloc);
514 if (useLjCombRule(nb->nbparam->vdwtype))
516 static_assert(sizeof(d_atdat->lj_comb[0]) == sizeof(float2),
517 "Size of the LJ parameters element should be equal to the size of float2.");
518 copyToDeviceBuffer(&d_atdat->lj_comb,
519 reinterpret_cast<const float2*>(nbat->params().lj_comb.data()), 0,
520 natoms, localStream, GpuApiCallBehavior::Async, nullptr);
524 static_assert(sizeof(d_atdat->atom_types[0]) == sizeof(nbat->params().type[0]),
525 "Sizes of host- and device-side atom types should be the same.");
526 copyToDeviceBuffer(&d_atdat->atom_types, nbat->params().type.data(), 0, natoms, localStream,
527 GpuApiCallBehavior::Async, nullptr);
532 timers->atdat.closeTimingRegion(localStream);
536 void gpu_free(NbnxmGpu* nb)
539 cu_atomdata_t* atdat;
548 nbparam = nb->nbparam;
550 if ((!nbparam->coulomb_tab)
551 && (nbparam->eeltype == eelTypeEWALD_TAB || nbparam->eeltype == eelTypeEWALD_TAB_TWIN))
553 destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
556 stat = cudaEventDestroy(nb->nonlocal_done);
557 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
558 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
559 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
563 if (!useLjCombRule(nb->nbparam->vdwtype))
565 destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
568 if (nbparam->vdwtype == evdwTypeEWALDGEOM || nbparam->vdwtype == evdwTypeEWALDLB)
570 destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
573 freeDeviceBuffer(&atdat->shift_vec);
574 freeDeviceBuffer(&atdat->fshift);
576 freeDeviceBuffer(&atdat->e_lj);
577 freeDeviceBuffer(&atdat->e_el);
579 freeDeviceBuffer(&atdat->f);
580 freeDeviceBuffer(&atdat->xq);
581 freeDeviceBuffer(&atdat->atom_types);
582 freeDeviceBuffer(&atdat->lj_comb);
585 auto* plist = nb->plist[InteractionLocality::Local];
586 freeDeviceBuffer(&plist->sci);
587 freeDeviceBuffer(&plist->cj4);
588 freeDeviceBuffer(&plist->imask);
589 freeDeviceBuffer(&plist->excl);
591 if (nb->bUseTwoStreams)
593 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
594 freeDeviceBuffer(&plist_nl->sci);
595 freeDeviceBuffer(&plist_nl->cj4);
596 freeDeviceBuffer(&plist_nl->imask);
597 freeDeviceBuffer(&plist_nl->excl);
602 pfree(nb->nbst.e_lj);
603 nb->nbst.e_lj = nullptr;
605 pfree(nb->nbst.e_el);
606 nb->nbst.e_el = nullptr;
608 pfree(nb->nbst.fshift);
609 nb->nbst.fshift = nullptr;
618 fprintf(debug, "Cleaned up CUDA data structures.\n");
622 //! This function is documented in the header file
623 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(NbnxmGpu* nb)
625 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
628 void gpu_reset_timings(nonbonded_verlet_t* nbv)
630 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
632 init_timings(nbv->gpu_nbv->timings);
636 int gpu_min_ci_balanced(NbnxmGpu* nb)
638 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
642 gmx_bool gpu_is_kernel_ewald_analytical(const NbnxmGpu* nb)
644 return ((nb->nbparam->eeltype == eelTypeEWALD_ANA) || (nb->nbparam->eeltype == eelTypeEWALD_ANA_TWIN));
647 void* gpu_get_xq(NbnxmGpu* nb)
651 return static_cast<void*>(nb->atdat->xq);
654 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
658 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
661 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
665 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fshift);
668 /* Initialization for X buffer operations on GPU. */
669 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
670 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
672 const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
673 bool bDoTime = gpu_nbv->bDoTime;
674 const int maxNumColumns = gridSet.numColumnsMax();
676 reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
677 &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, *gpu_nbv->deviceContext_);
678 reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
679 &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, *gpu_nbv->deviceContext_);
681 for (unsigned int g = 0; g < gridSet.grids().size(); g++)
684 const Nbnxm::Grid& grid = gridSet.grids()[g];
686 const int numColumns = grid.numColumns();
687 const int* atomIndices = gridSet.atomIndices().data();
688 const int atomIndicesSize = gridSet.atomIndices().size();
689 const int* cxy_na = grid.cxy_na().data();
690 const int* cxy_ind = grid.cxy_ind().data();
692 reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
693 &gpu_nbv->atomIndicesSize_alloc, *gpu_nbv->deviceContext_);
695 if (atomIndicesSize > 0)
700 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
703 copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, deviceStream,
704 GpuApiCallBehavior::Async, nullptr);
708 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
716 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
719 int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
720 copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, deviceStream,
721 GpuApiCallBehavior::Async, nullptr);
725 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
730 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
733 destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
734 copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, deviceStream,
735 GpuApiCallBehavior::Async, nullptr);
739 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
744 // The above data is transferred on the local stream but is a
745 // dependency of the nonlocal stream (specifically the nonlocal X
746 // buf ops kernel). We therefore set a dependency to ensure
747 // that the nonlocal stream waits on the local stream here.
748 // This call records an event in the local stream:
749 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
750 // ...and this call instructs the nonlocal stream to wait on that event:
751 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
756 /* Initialization for F buffer operations on GPU. */
757 void nbnxn_gpu_init_add_nbat_f_to_f(const int* cell,
760 GpuEventSynchronizer* const localReductionDone)
763 const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
765 GMX_ASSERT(localReductionDone, "localReductionDone should be a valid pointer");
766 gpu_nbv->localFReductionDone = localReductionDone;
768 if (natoms_total > 0)
770 reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc,
771 *gpu_nbv->deviceContext_);
772 copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, deviceStream,
773 GpuApiCallBehavior::Async, nullptr);