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/gpu_hw_info.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/pairlistsets.h"
69 #include "gromacs/pbcutil/ishift.h"
70 #include "gromacs/timing/gpu_timing.h"
71 #include "gromacs/utility/basedefinitions.h"
72 #include "gromacs/utility/cstringutil.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/real.h"
75 #include "gromacs/utility/smalloc.h"
77 #include "nbnxm_cuda.h"
82 /* This is a heuristically determined parameter for the Kepler
83 * and Maxwell architectures for the minimum size of ci lists by multiplying
84 * this constant with the # of multiprocessors on the current device.
85 * Since the maximum number of blocks per multiprocessor is 16, the ideal
86 * count for small systems is 32 or 48 blocks per multiprocessor. Because
87 * there is a bit of fluctuations in the generated block counts, we use
88 * a target of 44 instead of the ideal value of 48.
90 static unsigned int gpu_min_ci_balanced_factor = 44;
93 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
96 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t* nbparam);
98 /*! \brief Return whether combination rules are used.
100 * \param[in] pointer to nonbonded paramter struct
101 * \return true if combination rules are used in this run, false otherwise
103 static inline bool useLjCombRule(const cu_nbparam_t* nbparam)
105 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM || nbparam->vdwtype == evdwCuCUTCOMBLB);
108 /*! \brief Initialized the Ewald Coulomb correction GPU table.
110 Tabulates the Ewald Coulomb force and initializes the size/scale
111 and the table GPU array. If called with an already allocated table,
112 it just re-uploads the table.
114 static void init_ewald_coulomb_force_table(const EwaldCorrectionTables& tables,
116 const DeviceContext& deviceContext)
118 if (nbp->coulomb_tab != nullptr)
120 nbnxn_cuda_free_nbparam_table(nbp);
123 nbp->coulomb_tab_scale = tables.scale;
124 initParamLookupTable(&nbp->coulomb_tab, &nbp->coulomb_tab_texobj, tables.tableF.data(),
125 tables.tableF.size(), deviceContext);
129 /*! Initializes the atomdata structure first time, it only gets filled at
131 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes, const DeviceContext& deviceContext)
134 allocateDeviceBuffer(&ad->shift_vec, SHIFTS, deviceContext);
135 ad->bShiftVecUploaded = false;
137 allocateDeviceBuffer(&ad->fshift, SHIFTS, deviceContext);
138 allocateDeviceBuffer(&ad->e_lj, 1, deviceContext);
139 allocateDeviceBuffer(&ad->e_el, 1, deviceContext);
141 /* initialize to nullptr poiters to data that is not allocated here and will
142 need reallocation in nbnxn_cuda_init_atomdata */
146 /* size -1 indicates that the respective array hasn't been initialized yet */
151 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
152 earlier GPUs, single or twin cut-off. */
153 static int pick_ewald_kernel_type(const interaction_const_t& ic)
155 bool bTwinCut = (ic.rcoulomb != ic.rvdw);
156 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
159 /* Benchmarking/development environment variables to force the use of
160 analytical or tabulated Ewald kernel. */
161 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
162 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
164 if (bForceAnalyticalEwald && bForceTabulatedEwald)
167 "Both analytical and tabulated Ewald CUDA non-bonded kernels "
168 "requested through environment variables.");
171 /* By default use analytical Ewald. */
172 bUseAnalyticalEwald = true;
173 if (bForceAnalyticalEwald)
177 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
180 else if (bForceTabulatedEwald)
182 bUseAnalyticalEwald = false;
186 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
190 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
191 forces it (use it for debugging/benchmarking only). */
192 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
194 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
198 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
204 /*! Copies all parameters related to the cut-off from ic to nbp */
205 static void set_cutoff_parameters(cu_nbparam_t* nbp, const interaction_const_t* ic, const PairlistParams& listParams)
207 nbp->ewald_beta = ic->ewaldcoeff_q;
208 nbp->sh_ewald = ic->sh_ewald;
209 nbp->epsfac = ic->epsfac;
210 nbp->two_k_rf = 2.0 * ic->k_rf;
211 nbp->c_rf = ic->c_rf;
212 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
213 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
214 nbp->rlistOuter_sq = listParams.rlistOuter * listParams.rlistOuter;
215 nbp->rlistInner_sq = listParams.rlistInner * listParams.rlistInner;
216 nbp->useDynamicPruning = listParams.useDynamicPruning;
218 nbp->sh_lj_ewald = ic->sh_lj_ewald;
219 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
221 nbp->rvdw_switch = ic->rvdw_switch;
222 nbp->dispersion_shift = ic->dispersion_shift;
223 nbp->repulsion_shift = ic->repulsion_shift;
224 nbp->vdw_switch = ic->vdw_switch;
227 /*! Initializes the nonbonded parameter data structure. */
228 static void init_nbparam(cu_nbparam_t* nbp,
229 const interaction_const_t* ic,
230 const PairlistParams& listParams,
231 const nbnxn_atomdata_t::Params& nbatParams,
232 const DeviceContext& deviceContext)
236 ntypes = nbatParams.numTypes;
238 set_cutoff_parameters(nbp, ic, listParams);
240 /* The kernel code supports LJ combination rules (geometric and LB) for
241 * all kernel types, but we only generate useful combination rule kernels.
242 * We currently only use LJ combination rule (geometric and LB) kernels
243 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
244 * with PME and 20% with RF, the other kernels speed up about half as much.
245 * For LJ force-switch the geometric rule would give 7% speed-up, but this
246 * combination is rarely used. LJ force-switch with LB rule is more common,
247 * but gives only 1% speed-up.
249 if (ic->vdwtype == evdwCUT)
251 switch (ic->vdw_modifier)
254 case eintmodPOTSHIFT:
255 switch (nbatParams.comb_rule)
257 case ljcrNONE: nbp->vdwtype = evdwCuCUT; break;
258 case ljcrGEOM: nbp->vdwtype = evdwCuCUTCOMBGEOM; break;
259 case ljcrLB: nbp->vdwtype = evdwCuCUTCOMBLB; break;
262 "The requested LJ combination rule is not implemented in the CUDA "
263 "GPU accelerated kernels!");
266 case eintmodFORCESWITCH: nbp->vdwtype = evdwCuFSWITCH; break;
267 case eintmodPOTSWITCH: nbp->vdwtype = evdwCuPSWITCH; break;
270 "The requested VdW interaction modifier is not implemented in the CUDA GPU "
271 "accelerated kernels!");
274 else if (ic->vdwtype == evdwPME)
276 if (ic->ljpme_comb_rule == ljcrGEOM)
278 assert(nbatParams.comb_rule == ljcrGEOM);
279 nbp->vdwtype = evdwCuEWALDGEOM;
283 assert(nbatParams.comb_rule == ljcrLB);
284 nbp->vdwtype = evdwCuEWALDLB;
290 "The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
293 if (ic->eeltype == eelCUT)
295 nbp->eeltype = eelCuCUT;
297 else if (EEL_RF(ic->eeltype))
299 nbp->eeltype = eelCuRF;
301 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
303 nbp->eeltype = pick_ewald_kernel_type(*ic);
307 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
309 "The requested electrostatics type is not implemented in the CUDA GPU accelerated "
313 /* generate table for PME */
314 nbp->coulomb_tab = nullptr;
315 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
317 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
318 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, deviceContext);
321 /* set up LJ parameter lookup table */
322 if (!useLjCombRule(nbp))
324 initParamLookupTable(&nbp->nbfp, &nbp->nbfp_texobj, nbatParams.nbfp.data(),
325 2 * ntypes * ntypes, deviceContext);
328 /* set up LJ-PME parameter lookup table */
329 if (ic->vdwtype == evdwPME)
331 initParamLookupTable(&nbp->nbfp_comb, &nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(),
332 2 * ntypes, deviceContext);
336 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
337 * electrostatic type switching to twin cut-off (or back) if needed. */
338 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
340 if (!nbv || !nbv->useGpu())
344 NbnxmGpu* nb = nbv->gpu_nbv;
345 cu_nbparam_t* nbp = nbv->gpu_nbv->nbparam;
347 set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
349 nbp->eeltype = pick_ewald_kernel_type(*ic);
351 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
352 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, *nb->deviceContext_);
355 /*! Initializes the pair list data structure. */
356 static void init_plist(cu_plist_t* pl)
358 /* initialize to nullptr pointers to data that is not allocated here and will
359 need reallocation in nbnxn_gpu_init_pairlist */
365 /* size -1 indicates that the respective array hasn't been initialized yet */
372 pl->imask_nalloc = -1;
374 pl->excl_nalloc = -1;
375 pl->haveFreshList = false;
378 /*! Initializes the timings data structure. */
379 static void init_timings(gmx_wallclock_gpu_nbnxn_t* t)
388 for (i = 0; i < 2; i++)
390 for (j = 0; j < 2; j++)
392 t->ktime[i][j].t = 0.0;
393 t->ktime[i][j].c = 0;
397 t->pruneTime.t = 0.0;
398 t->dynamicPruneTime.c = 0;
399 t->dynamicPruneTime.t = 0.0;
402 /*! Initializes simulation constant data. */
403 static void cuda_init_const(NbnxmGpu* nb,
404 const interaction_const_t* ic,
405 const PairlistParams& listParams,
406 const nbnxn_atomdata_t::Params& nbatParams)
408 init_atomdata_first(nb->atdat, nbatParams.numTypes, *nb->deviceContext_);
409 init_nbparam(nb->nbparam, ic, listParams, nbatParams, *nb->deviceContext_);
411 /* clear energy and shift force outputs */
412 nbnxn_cuda_clear_e_fshift(nb);
415 NbnxmGpu* gpu_init(const gmx::DeviceStreamManager& deviceStreamManager,
416 const interaction_const_t* ic,
417 const PairlistParams& listParams,
418 const nbnxn_atomdata_t* nbat,
419 bool bLocalAndNonlocal)
423 auto nb = new NbnxmGpu();
424 nb->deviceContext_ = &deviceStreamManager.context();
426 snew(nb->nbparam, 1);
427 snew(nb->plist[InteractionLocality::Local], 1);
428 if (bLocalAndNonlocal)
430 snew(nb->plist[InteractionLocality::NonLocal], 1);
433 nb->bUseTwoStreams = bLocalAndNonlocal;
435 nb->timers = new cu_timers_t();
436 snew(nb->timings, 1);
439 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
440 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
441 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
443 init_plist(nb->plist[InteractionLocality::Local]);
445 /* local/non-local GPU streams */
446 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedLocal),
447 "Local non-bonded stream should be initialized to use GPU for non-bonded.");
448 nb->deviceStreams[InteractionLocality::Local] =
449 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedLocal);
450 if (nb->bUseTwoStreams)
452 init_plist(nb->plist[InteractionLocality::NonLocal]);
454 /* Note that the device we're running on does not have to support
455 * priorities, because we are querying the priority range which in this
456 * case will be a single value.
458 GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedNonLocal),
459 "Non-local non-bonded stream should be initialized to use GPU for "
460 "non-bonded with domain decomposition.");
461 nb->deviceStreams[InteractionLocality::NonLocal] =
462 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedNonLocal);
466 /* init events for sychronization (timing disabled for performance reasons!) */
467 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
468 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
469 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
470 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
472 nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
474 /* WARNING: CUDA timings are incorrect with multiple streams.
475 * This is the main reason why they are disabled by default.
477 // TODO: Consider turning on by default when we can detect nr of streams.
478 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
482 init_timings(nb->timings);
485 /* set the kernel type for the current GPU */
486 /* pick L1 cache configuration */
487 cuda_set_cacheconfig();
489 cuda_init_const(nb, ic, listParams, nbat->params());
491 nb->atomIndicesSize = 0;
492 nb->atomIndicesSize_alloc = 0;
494 nb->ncxy_na_alloc = 0;
496 nb->ncxy_ind_alloc = 0;
502 fprintf(debug, "Initialized CUDA data structures.\n");
508 void gpu_init_pairlist(NbnxmGpu* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
511 bool bDoTime = (nb->bDoTime && !h_plist->sci.empty());
512 const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
513 cu_plist_t* d_plist = nb->plist[iloc];
515 if (d_plist->na_c < 0)
517 d_plist->na_c = h_plist->na_ci;
521 if (d_plist->na_c != h_plist->na_ci)
523 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
524 d_plist->na_c, h_plist->na_ci);
529 gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
533 iTimers.pl_h2d.openTimingRegion(deviceStream);
534 iTimers.didPairlistH2D = true;
537 const DeviceContext& deviceContext = *nb->deviceContext_;
539 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc,
541 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), deviceStream,
542 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
544 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc,
546 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), deviceStream,
547 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
549 reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
550 &d_plist->nimask, &d_plist->imask_nalloc, deviceContext);
552 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
553 &d_plist->excl_nalloc, deviceContext);
554 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), deviceStream,
555 GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
559 iTimers.pl_h2d.closeTimingRegion(deviceStream);
562 /* the next use of thist list we be the first one, so we need to prune */
563 d_plist->haveFreshList = true;
566 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
568 cu_atomdata_t* adat = nb->atdat;
569 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
571 /* only if we have a dynamic box */
572 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
574 static_assert(sizeof(adat->shift_vec[0]) == sizeof(nbatom->shift_vec[0]),
575 "Sizes of host- and device-side shift vectors should be the same.");
576 copyToDeviceBuffer(&adat->shift_vec, reinterpret_cast<const float3*>(nbatom->shift_vec.data()),
577 0, SHIFTS, localStream, GpuApiCallBehavior::Async, nullptr);
578 adat->bShiftVecUploaded = true;
582 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
583 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
585 cu_atomdata_t* adat = nb->atdat;
586 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
587 clearDeviceBufferAsync(&adat->f, 0, natoms_clear, localStream);
590 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
591 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
593 cu_atomdata_t* adat = nb->atdat;
594 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
596 clearDeviceBufferAsync(&adat->fshift, 0, SHIFTS, localStream);
597 clearDeviceBufferAsync(&adat->e_lj, 0, 1, localStream);
598 clearDeviceBufferAsync(&adat->e_el, 0, 1, localStream);
601 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
603 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
604 /* clear shift force array and energies if the outputs were
605 used in the current step */
608 nbnxn_cuda_clear_e_fshift(nb);
612 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
616 bool bDoTime = nb->bDoTime;
617 cu_timers_t* timers = nb->timers;
618 cu_atomdata_t* d_atdat = nb->atdat;
619 const DeviceContext& deviceContext = *nb->deviceContext_;
620 const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
622 natoms = nbat->numAtoms();
627 /* time async copy */
628 timers->atdat.openTimingRegion(localStream);
631 /* need to reallocate if we have to copy more atoms than the amount of space
632 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
633 if (natoms > d_atdat->nalloc)
635 nalloc = over_alloc_small(natoms);
637 /* free up first if the arrays have already been initialized */
638 if (d_atdat->nalloc != -1)
640 freeDeviceBuffer(&d_atdat->f);
641 freeDeviceBuffer(&d_atdat->xq);
642 freeDeviceBuffer(&d_atdat->atom_types);
643 freeDeviceBuffer(&d_atdat->lj_comb);
646 allocateDeviceBuffer(&d_atdat->f, nalloc, deviceContext);
647 allocateDeviceBuffer(&d_atdat->xq, nalloc, deviceContext);
648 if (useLjCombRule(nb->nbparam))
650 allocateDeviceBuffer(&d_atdat->lj_comb, nalloc, deviceContext);
654 allocateDeviceBuffer(&d_atdat->atom_types, nalloc, deviceContext);
657 d_atdat->nalloc = nalloc;
661 d_atdat->natoms = natoms;
662 d_atdat->natoms_local = nbat->natoms_local;
664 /* need to clear GPU f output if realloc happened */
667 nbnxn_cuda_clear_f(nb, nalloc);
670 if (useLjCombRule(nb->nbparam))
672 static_assert(sizeof(d_atdat->lj_comb[0]) == sizeof(float2),
673 "Size of the LJ parameters element should be equal to the size of float2.");
674 copyToDeviceBuffer(&d_atdat->lj_comb,
675 reinterpret_cast<const float2*>(nbat->params().lj_comb.data()), 0,
676 natoms, localStream, GpuApiCallBehavior::Async, nullptr);
680 static_assert(sizeof(d_atdat->atom_types[0]) == sizeof(nbat->params().type[0]),
681 "Sizes of host- and device-side atom types should be the same.");
682 copyToDeviceBuffer(&d_atdat->atom_types, nbat->params().type.data(), 0, natoms, localStream,
683 GpuApiCallBehavior::Async, nullptr);
688 timers->atdat.closeTimingRegion(localStream);
692 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t* nbparam)
694 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
696 destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
700 void gpu_free(NbnxmGpu* nb)
703 cu_atomdata_t* atdat;
704 cu_nbparam_t* nbparam;
712 nbparam = nb->nbparam;
714 nbnxn_cuda_free_nbparam_table(nbparam);
716 stat = cudaEventDestroy(nb->nonlocal_done);
717 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
718 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
719 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
723 if (!useLjCombRule(nb->nbparam))
725 destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
728 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
730 destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
733 freeDeviceBuffer(&atdat->shift_vec);
734 freeDeviceBuffer(&atdat->fshift);
736 freeDeviceBuffer(&atdat->e_lj);
737 freeDeviceBuffer(&atdat->e_el);
739 freeDeviceBuffer(&atdat->f);
740 freeDeviceBuffer(&atdat->xq);
741 freeDeviceBuffer(&atdat->atom_types);
742 freeDeviceBuffer(&atdat->lj_comb);
745 auto* plist = nb->plist[InteractionLocality::Local];
746 freeDeviceBuffer(&plist->sci);
747 freeDeviceBuffer(&plist->cj4);
748 freeDeviceBuffer(&plist->imask);
749 freeDeviceBuffer(&plist->excl);
751 if (nb->bUseTwoStreams)
753 auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
754 freeDeviceBuffer(&plist_nl->sci);
755 freeDeviceBuffer(&plist_nl->cj4);
756 freeDeviceBuffer(&plist_nl->imask);
757 freeDeviceBuffer(&plist_nl->excl);
762 pfree(nb->nbst.e_lj);
763 nb->nbst.e_lj = nullptr;
765 pfree(nb->nbst.e_el);
766 nb->nbst.e_el = nullptr;
768 pfree(nb->nbst.fshift);
769 nb->nbst.fshift = nullptr;
778 fprintf(debug, "Cleaned up CUDA data structures.\n");
782 //! This function is documented in the header file
783 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(NbnxmGpu* nb)
785 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
788 void gpu_reset_timings(nonbonded_verlet_t* nbv)
790 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
792 init_timings(nbv->gpu_nbv->timings);
796 int gpu_min_ci_balanced(NbnxmGpu* nb)
798 return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
802 gmx_bool gpu_is_kernel_ewald_analytical(const NbnxmGpu* nb)
804 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) || (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
807 void* gpu_get_xq(NbnxmGpu* nb)
811 return static_cast<void*>(nb->atdat->xq);
814 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
818 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
821 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
825 return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fshift);
828 /* Initialization for X buffer operations on GPU. */
829 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
830 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
832 const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
833 bool bDoTime = gpu_nbv->bDoTime;
834 const int maxNumColumns = gridSet.numColumnsMax();
836 reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
837 &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, *gpu_nbv->deviceContext_);
838 reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
839 &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, *gpu_nbv->deviceContext_);
841 for (unsigned int g = 0; g < gridSet.grids().size(); g++)
844 const Nbnxm::Grid& grid = gridSet.grids()[g];
846 const int numColumns = grid.numColumns();
847 const int* atomIndices = gridSet.atomIndices().data();
848 const int atomIndicesSize = gridSet.atomIndices().size();
849 const int* cxy_na = grid.cxy_na().data();
850 const int* cxy_ind = grid.cxy_ind().data();
852 reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
853 &gpu_nbv->atomIndicesSize_alloc, *gpu_nbv->deviceContext_);
855 if (atomIndicesSize > 0)
860 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
863 copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, deviceStream,
864 GpuApiCallBehavior::Async, nullptr);
868 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
876 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
879 int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
880 copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, deviceStream,
881 GpuApiCallBehavior::Async, nullptr);
885 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
890 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
893 destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
894 copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, deviceStream,
895 GpuApiCallBehavior::Async, nullptr);
899 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
904 // The above data is transferred on the local stream but is a
905 // dependency of the nonlocal stream (specifically the nonlocal X
906 // buf ops kernel). We therefore set a dependency to ensure
907 // that the nonlocal stream waits on the local stream here.
908 // This call records an event in the local stream:
909 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
910 // ...and this call instructs the nonlocal stream to wait on that event:
911 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
916 /* Initialization for F buffer operations on GPU. */
917 void nbnxn_gpu_init_add_nbat_f_to_f(const int* cell,
920 GpuEventSynchronizer* const localReductionDone)
923 const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
925 GMX_ASSERT(localReductionDone, "localReductionDone should be a valid pointer");
926 gpu_nbv->localFReductionDone = localReductionDone;
928 if (natoms_total > 0)
930 reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc,
931 *gpu_nbv->deviceContext_);
932 copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, deviceStream,
933 GpuApiCallBehavior::Async, nullptr);