2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
38 * \author Szilard Pall <pall.szilard@gmail.com>
47 // TODO We would like to move this down, but the way gmx_nbnxn_gpu_t
48 // is currently declared means this has to be before gpu_types.h
49 #include "nbnxm_cuda_types.h"
51 // TODO Remove this comment when the above order issue is resolved
52 #include "gromacs/gpu_utils/cudautils.cuh"
53 #include "gromacs/gpu_utils/gpu_utils.h"
54 #include "gromacs/gpu_utils/pmalloc_cuda.h"
55 #include "gromacs/hardware/gpu_hw_info.h"
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/mdlib/force_flags.h"
58 #include "gromacs/mdtypes/interaction_const.h"
59 #include "gromacs/mdtypes/md_enums.h"
60 #include "gromacs/nbnxm/atomdata.h"
61 #include "gromacs/nbnxm/gpu_data_mgmt.h"
62 #include "gromacs/nbnxm/gridset.h"
63 #include "gromacs/nbnxm/nbnxm.h"
64 #include "gromacs/nbnxm/nbnxm_gpu.h"
65 #include "gromacs/nbnxm/pairlistsets.h"
66 #include "gromacs/pbcutil/ishift.h"
67 #include "gromacs/timing/gpu_timing.h"
68 #include "gromacs/utility/basedefinitions.h"
69 #include "gromacs/utility/cstringutil.h"
70 #include "gromacs/utility/fatalerror.h"
71 #include "gromacs/utility/real.h"
72 #include "gromacs/utility/smalloc.h"
74 #include "nbnxm_cuda.h"
79 /* This is a heuristically determined parameter for the Kepler
80 * and Maxwell architectures for the minimum size of ci lists by multiplying
81 * this constant with the # of multiprocessors on the current device.
82 * Since the maximum number of blocks per multiprocessor is 16, the ideal
83 * count for small systems is 32 or 48 blocks per multiprocessor. Because
84 * there is a bit of fluctuations in the generated block counts, we use
85 * a target of 44 instead of the ideal value of 48.
87 static unsigned int gpu_min_ci_balanced_factor = 44;
90 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
93 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam);
95 /*! \brief Return whether combination rules are used.
97 * \param[in] pointer to nonbonded paramter struct
98 * \return true if combination rules are used in this run, false otherwise
100 static inline bool useLjCombRule(const cu_nbparam_t *nbparam)
102 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
103 nbparam->vdwtype == evdwCuCUTCOMBLB);
106 /*! \brief Initialized the Ewald Coulomb correction GPU table.
108 Tabulates the Ewald Coulomb force and initializes the size/scale
109 and the table GPU array. If called with an already allocated table,
110 it just re-uploads the table.
112 static void init_ewald_coulomb_force_table(const EwaldCorrectionTables &tables,
115 if (nbp->coulomb_tab != nullptr)
117 nbnxn_cuda_free_nbparam_table(nbp);
120 nbp->coulomb_tab_scale = tables.scale;
121 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
122 tables.tableF.data(), tables.tableF.size());
126 /*! Initializes the atomdata structure first time, it only gets filled at
128 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
133 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
134 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
135 ad->bShiftVecUploaded = false;
137 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
138 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
140 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
141 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
142 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
143 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
145 /* initialize to nullptr poiters to data that is not allocated here and will
146 need reallocation in nbnxn_cuda_init_atomdata */
150 /* size -1 indicates that the respective array hasn't been initialized yet */
155 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
156 earlier GPUs, single or twin cut-off. */
157 static int pick_ewald_kernel_type(bool bTwinCut)
159 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
162 /* Benchmarking/development environment variables to force the use of
163 analytical or tabulated Ewald kernel. */
164 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
165 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
167 if (bForceAnalyticalEwald && bForceTabulatedEwald)
169 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
170 "requested through environment variables.");
173 /* By default use analytical Ewald. */
174 bUseAnalyticalEwald = true;
175 if (bForceAnalyticalEwald)
179 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
182 else if (bForceTabulatedEwald)
184 bUseAnalyticalEwald = false;
188 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
192 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
193 forces it (use it for debugging/benchmarking only). */
194 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
196 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
200 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
206 /*! Copies all parameters related to the cut-off from ic to nbp */
207 static void set_cutoff_parameters(cu_nbparam_t *nbp,
208 const interaction_const_t *ic,
209 const PairlistParams &listParams)
211 nbp->ewald_beta = ic->ewaldcoeff_q;
212 nbp->sh_ewald = ic->sh_ewald;
213 nbp->epsfac = ic->epsfac;
214 nbp->two_k_rf = 2.0 * ic->k_rf;
215 nbp->c_rf = ic->c_rf;
216 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
217 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
218 nbp->rlistOuter_sq = listParams.rlistOuter * listParams.rlistOuter;
219 nbp->rlistInner_sq = listParams.rlistInner * listParams.rlistInner;
220 nbp->useDynamicPruning = listParams.useDynamicPruning;
222 nbp->sh_lj_ewald = ic->sh_lj_ewald;
223 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
225 nbp->rvdw_switch = ic->rvdw_switch;
226 nbp->dispersion_shift = ic->dispersion_shift;
227 nbp->repulsion_shift = ic->repulsion_shift;
228 nbp->vdw_switch = ic->vdw_switch;
231 /*! Initializes the nonbonded parameter data structure. */
232 static void init_nbparam(cu_nbparam_t *nbp,
233 const interaction_const_t *ic,
234 const PairlistParams &listParams,
235 const nbnxn_atomdata_t::Params &nbatParams)
239 ntypes = nbatParams.numTypes;
241 set_cutoff_parameters(nbp, ic, listParams);
243 /* The kernel code supports LJ combination rules (geometric and LB) for
244 * all kernel types, but we only generate useful combination rule kernels.
245 * We currently only use LJ combination rule (geometric and LB) kernels
246 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
247 * with PME and 20% with RF, the other kernels speed up about half as much.
248 * For LJ force-switch the geometric rule would give 7% speed-up, but this
249 * combination is rarely used. LJ force-switch with LB rule is more common,
250 * but gives only 1% speed-up.
252 if (ic->vdwtype == evdwCUT)
254 switch (ic->vdw_modifier)
257 case eintmodPOTSHIFT:
258 switch (nbatParams.comb_rule)
261 nbp->vdwtype = evdwCuCUT;
264 nbp->vdwtype = evdwCuCUTCOMBGEOM;
267 nbp->vdwtype = evdwCuCUTCOMBLB;
270 gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
273 case eintmodFORCESWITCH:
274 nbp->vdwtype = evdwCuFSWITCH;
276 case eintmodPOTSWITCH:
277 nbp->vdwtype = evdwCuPSWITCH;
280 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
283 else if (ic->vdwtype == evdwPME)
285 if (ic->ljpme_comb_rule == ljcrGEOM)
287 assert(nbatParams.comb_rule == ljcrGEOM);
288 nbp->vdwtype = evdwCuEWALDGEOM;
292 assert(nbatParams.comb_rule == ljcrLB);
293 nbp->vdwtype = evdwCuEWALDLB;
298 gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
301 if (ic->eeltype == eelCUT)
303 nbp->eeltype = eelCuCUT;
305 else if (EEL_RF(ic->eeltype))
307 nbp->eeltype = eelCuRF;
309 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
311 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
312 nbp->eeltype = pick_ewald_kernel_type(false);
316 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
317 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
320 /* generate table for PME */
321 nbp->coulomb_tab = nullptr;
322 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
324 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
325 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
328 /* set up LJ parameter lookup table */
329 if (!useLjCombRule(nbp))
331 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
332 nbatParams.nbfp.data(), 2*ntypes*ntypes);
335 /* set up LJ-PME parameter lookup table */
336 if (ic->vdwtype == evdwPME)
338 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
339 nbatParams.nbfp_comb.data(), 2*ntypes);
343 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
344 * electrostatic type switching to twin cut-off (or back) if needed. */
345 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t *nbv,
346 const interaction_const_t *ic)
348 if (!nbv || !nbv->useGpu())
352 cu_nbparam_t *nbp = nbv->gpu_nbv->nbparam;
354 set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
356 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw);
358 GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
359 init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
362 /*! Initializes the pair list data structure. */
363 static void init_plist(cu_plist_t *pl)
365 /* initialize to nullptr pointers to data that is not allocated here and will
366 need reallocation in nbnxn_gpu_init_pairlist */
372 /* size -1 indicates that the respective array hasn't been initialized yet */
379 pl->imask_nalloc = -1;
381 pl->excl_nalloc = -1;
382 pl->haveFreshList = false;
385 /*! Initializes the timings data structure. */
386 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
395 for (i = 0; i < 2; i++)
397 for (j = 0; j < 2; j++)
399 t->ktime[i][j].t = 0.0;
400 t->ktime[i][j].c = 0;
404 t->pruneTime.t = 0.0;
405 t->dynamicPruneTime.c = 0;
406 t->dynamicPruneTime.t = 0.0;
409 /*! Initializes simulation constant data. */
410 static void cuda_init_const(gmx_nbnxn_cuda_t *nb,
411 const interaction_const_t *ic,
412 const PairlistParams &listParams,
413 const nbnxn_atomdata_t::Params &nbatParams)
415 init_atomdata_first(nb->atdat, nbatParams.numTypes);
416 init_nbparam(nb->nbparam, ic, listParams, nbatParams);
418 /* clear energy and shift force outputs */
419 nbnxn_cuda_clear_e_fshift(nb);
423 gpu_init(const gmx_device_info_t *deviceInfo,
424 const interaction_const_t *ic,
425 const PairlistParams &listParams,
426 const nbnxn_atomdata_t *nbat,
428 gmx_bool bLocalAndNonlocal)
432 gmx_nbnxn_cuda_t *nb;
435 snew(nb->nbparam, 1);
436 snew(nb->plist[InteractionLocality::Local], 1);
437 if (bLocalAndNonlocal)
439 snew(nb->plist[InteractionLocality::NonLocal], 1);
442 nb->bUseTwoStreams = bLocalAndNonlocal;
444 nb->timers = new cu_timers_t();
445 snew(nb->timings, 1);
448 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
449 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
450 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
452 init_plist(nb->plist[InteractionLocality::Local]);
454 /* set device info, just point it to the right GPU among the detected ones */
455 nb->dev_info = deviceInfo;
457 /* local/non-local GPU streams */
458 stat = cudaStreamCreate(&nb->stream[InteractionLocality::Local]);
459 CU_RET_ERR(stat, "cudaStreamCreate on stream[InterationLocality::Local] failed");
460 if (nb->bUseTwoStreams)
462 init_plist(nb->plist[InteractionLocality::NonLocal]);
464 /* Note that the device we're running on does not have to support
465 * priorities, because we are querying the priority range which in this
466 * case will be a single value.
468 int highest_priority;
469 stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
470 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
472 stat = cudaStreamCreateWithPriority(&nb->stream[InteractionLocality::NonLocal],
475 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[InteractionLocality::NonLocal] failed");
478 /* init events for sychronization (timing disabled for performance reasons!) */
479 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
480 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
481 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
482 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
484 /* WARNING: CUDA timings are incorrect with multiple streams.
485 * This is the main reason why they are disabled by default.
487 // TODO: Consider turning on by default when we can detect nr of streams.
488 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
492 init_timings(nb->timings);
495 /* set the kernel type for the current GPU */
496 /* pick L1 cache configuration */
497 cuda_set_cacheconfig();
499 cuda_init_const(nb, ic, listParams, nbat->params());
502 nb->natoms_alloc = 0;
503 nb->atomIndicesSize = 0;
504 nb->atomIndicesSize_alloc = 0;
506 nb->ncxy_na_alloc = 0;
508 nb->ncxy_ind_alloc = 0;
510 nb->nfrvec_alloc = 0;
516 fprintf(debug, "Initialized CUDA data structures.\n");
522 void gpu_init_pairlist(gmx_nbnxn_cuda_t *nb,
523 const NbnxnPairlistGpu *h_plist,
524 const InteractionLocality iloc)
527 bool bDoTime = (nb->bDoTime && !h_plist->sci.empty());
528 cudaStream_t stream = nb->stream[iloc];
529 cu_plist_t *d_plist = nb->plist[iloc];
531 if (d_plist->na_c < 0)
533 d_plist->na_c = h_plist->na_ci;
537 if (d_plist->na_c != h_plist->na_ci)
539 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
540 d_plist->na_c, h_plist->na_ci);
545 gpu_timers_t::Interaction &iTimers = nb->timers->interaction[iloc];
549 iTimers.pl_h2d.openTimingRegion(stream);
550 iTimers.didPairlistH2D = true;
553 Context context = nullptr;
555 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(),
556 &d_plist->nsci, &d_plist->sci_nalloc, context);
557 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(),
558 stream, GpuApiCallBehavior::Async,
559 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
561 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(),
562 &d_plist->ncj4, &d_plist->cj4_nalloc, context);
563 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(),
564 stream, GpuApiCallBehavior::Async,
565 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
567 reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size()*c_nbnxnGpuClusterpairSplit,
568 &d_plist->nimask, &d_plist->imask_nalloc, context);
570 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(),
571 &d_plist->nexcl, &d_plist->excl_nalloc, context);
572 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(),
573 stream, GpuApiCallBehavior::Async,
574 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
578 iTimers.pl_h2d.closeTimingRegion(stream);
581 /* the next use of thist list we be the first one, so we need to prune */
582 d_plist->haveFreshList = true;
585 void gpu_upload_shiftvec(gmx_nbnxn_cuda_t *nb,
586 const nbnxn_atomdata_t *nbatom)
588 cu_atomdata_t *adat = nb->atdat;
589 cudaStream_t ls = nb->stream[InteractionLocality::Local];
591 /* only if we have a dynamic box */
592 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
594 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(),
595 SHIFTS * sizeof(*adat->shift_vec), ls);
596 adat->bShiftVecUploaded = true;
600 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
601 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
604 cu_atomdata_t *adat = nb->atdat;
605 cudaStream_t ls = nb->stream[InteractionLocality::Local];
607 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
608 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
611 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
612 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
615 cu_atomdata_t *adat = nb->atdat;
616 cudaStream_t ls = nb->stream[InteractionLocality::Local];
618 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
619 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
620 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
621 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
622 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
623 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
626 void gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
628 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
629 /* clear shift force array and energies if the outputs were
630 used in the current step */
631 if (flags & GMX_FORCE_VIRIAL)
633 nbnxn_cuda_clear_e_fshift(nb);
637 void gpu_init_atomdata(gmx_nbnxn_cuda_t *nb,
638 const nbnxn_atomdata_t *nbat)
643 bool bDoTime = nb->bDoTime;
644 cu_timers_t *timers = nb->timers;
645 cu_atomdata_t *d_atdat = nb->atdat;
646 cudaStream_t ls = nb->stream[InteractionLocality::Local];
648 natoms = nbat->numAtoms();
653 /* time async copy */
654 timers->atdat.openTimingRegion(ls);
657 /* need to reallocate if we have to copy more atoms than the amount of space
658 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
659 if (natoms > d_atdat->nalloc)
661 nalloc = over_alloc_small(natoms);
663 /* free up first if the arrays have already been initialized */
664 if (d_atdat->nalloc != -1)
666 freeDeviceBuffer(&d_atdat->f);
667 freeDeviceBuffer(&d_atdat->xq);
668 freeDeviceBuffer(&d_atdat->atom_types);
669 freeDeviceBuffer(&d_atdat->lj_comb);
672 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
673 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
674 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
675 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
676 if (useLjCombRule(nb->nbparam))
678 stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
679 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
683 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
684 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
687 d_atdat->nalloc = nalloc;
691 d_atdat->natoms = natoms;
692 d_atdat->natoms_local = nbat->natoms_local;
694 /* need to clear GPU f output if realloc happened */
697 nbnxn_cuda_clear_f(nb, nalloc);
700 if (useLjCombRule(nb->nbparam))
702 cu_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(),
703 natoms*sizeof(*d_atdat->lj_comb), ls);
707 cu_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(),
708 natoms*sizeof(*d_atdat->atom_types), ls);
713 timers->atdat.closeTimingRegion(ls);
717 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam)
719 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
721 destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
725 void gpu_free(gmx_nbnxn_cuda_t *nb)
728 cu_atomdata_t *atdat;
729 cu_nbparam_t *nbparam;
737 nbparam = nb->nbparam;
739 nbnxn_cuda_free_nbparam_table(nbparam);
741 stat = cudaEventDestroy(nb->nonlocal_done);
742 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
743 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
744 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
749 /* The non-local counters/stream (second in the array) are needed only with DD. */
750 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
752 stat = cudaStreamDestroy(nb->stream[i]);
753 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
757 if (!useLjCombRule(nb->nbparam))
759 destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj);
763 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
765 destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
768 stat = cudaFree(atdat->shift_vec);
769 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
770 stat = cudaFree(atdat->fshift);
771 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
773 stat = cudaFree(atdat->e_lj);
774 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
775 stat = cudaFree(atdat->e_el);
776 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
778 freeDeviceBuffer(&atdat->f);
779 freeDeviceBuffer(&atdat->xq);
780 freeDeviceBuffer(&atdat->atom_types);
781 freeDeviceBuffer(&atdat->lj_comb);
784 auto *plist = nb->plist[InteractionLocality::Local];
785 freeDeviceBuffer(&plist->sci);
786 freeDeviceBuffer(&plist->cj4);
787 freeDeviceBuffer(&plist->imask);
788 freeDeviceBuffer(&plist->excl);
790 if (nb->bUseTwoStreams)
792 auto *plist_nl = nb->plist[InteractionLocality::NonLocal];
793 freeDeviceBuffer(&plist_nl->sci);
794 freeDeviceBuffer(&plist_nl->cj4);
795 freeDeviceBuffer(&plist_nl->imask);
796 freeDeviceBuffer(&plist_nl->excl);
801 pfree(nb->nbst.e_lj);
802 nb->nbst.e_lj = nullptr;
804 pfree(nb->nbst.e_el);
805 nb->nbst.e_el = nullptr;
807 pfree(nb->nbst.fshift);
808 nb->nbst.fshift = nullptr;
817 fprintf(debug, "Cleaned up CUDA data structures.\n");
821 //! This function is documented in the header file
822 gmx_wallclock_gpu_nbnxn_t *gpu_get_timings(gmx_nbnxn_cuda_t *nb)
824 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
827 void gpu_reset_timings(nonbonded_verlet_t* nbv)
829 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
831 init_timings(nbv->gpu_nbv->timings);
835 int gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
837 return nb != nullptr ?
838 gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
842 gmx_bool gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
844 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
845 (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
848 void *gpu_get_command_stream(gmx_nbnxn_gpu_t *nb,
849 const InteractionLocality iloc)
853 return static_cast<void *>(&nb->stream[iloc]);
856 void *gpu_get_xq(gmx_nbnxn_gpu_t *nb)
860 return static_cast<void *>(nb->atdat->xq);
863 void *gpu_get_f(gmx_nbnxn_gpu_t *nb)
867 return static_cast<void *>(nb->atdat->f);
870 rvec *gpu_get_fshift(gmx_nbnxn_gpu_t *nb)
874 return reinterpret_cast<rvec *>(nb->atdat->fshift);
877 /* Initialization for X buffer operations on GPU. */
878 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
879 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet &gridSet,
880 gmx_nbnxn_gpu_t *gpu_nbv)
882 cudaStream_t stream = gpu_nbv->stream[InteractionLocality::Local];
883 bool bDoTime = gpu_nbv->bDoTime;
884 const int maxNumColumns = gridSet.numColumnsMax();
886 reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns*gridSet.grids().size(),
887 &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, nullptr);
888 reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns*gridSet.grids().size(),
889 &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, nullptr);
891 for (unsigned int g = 0; g < gridSet.grids().size(); g++)
894 const Nbnxm::Grid &grid = gridSet.grids()[g];
896 const int numColumns = grid.numColumns();
897 const int *atomIndices = gridSet.atomIndices().data();
898 const int atomIndicesSize = gridSet.atomIndices().size();
899 const int *cxy_na = grid.cxy_na().data();
900 const int *cxy_ind = grid.cxy_ind().data();
901 const int numRealAtomsTotal = gridSet.numRealAtomsTotal();
903 reallocateDeviceBuffer(&gpu_nbv->xrvec, numRealAtomsTotal, &gpu_nbv->natoms, &gpu_nbv->natoms_alloc, nullptr);
904 reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize, &gpu_nbv->atomIndicesSize_alloc, nullptr);
906 if (atomIndicesSize > 0)
911 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
914 copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, stream, GpuApiCallBehavior::Async, nullptr);
918 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
927 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
930 int* destPtr = &gpu_nbv->cxy_na[maxNumColumns*g];
931 copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
935 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
940 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
943 destPtr = &gpu_nbv->cxy_ind[maxNumColumns*g];
944 copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
948 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
954 // The above data is transferred on the local stream but is a
955 // dependency of the nonlocal stream (specifically the nonlocal X
956 // buf ops kernel). We therefore set a dependency to ensure
957 // that the nonlocal stream waits on the local stream here.
958 // This call records an event in the local stream:
959 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
960 // ...and this call instructs the nonlocal stream to wait on that event:
961 nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
966 /* Initialization for F buffer operations on GPU. */
967 void nbnxn_gpu_init_add_nbat_f_to_f(const int *cell,
968 gmx_nbnxn_gpu_t *gpu_nbv,
972 cudaStream_t stream = gpu_nbv->stream[InteractionLocality::Local];
974 reallocateDeviceBuffer(&gpu_nbv->frvec, natoms_total, &gpu_nbv->nfrvec, &gpu_nbv->nfrvec_alloc, nullptr);
976 if (natoms_total > 0)
978 reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc, nullptr);
979 copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, stream, GpuApiCallBehavior::Async, nullptr);