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/pairlistsets.h"
65 #include "gromacs/pbcutil/ishift.h"
66 #include "gromacs/timing/gpu_timing.h"
67 #include "gromacs/utility/basedefinitions.h"
68 #include "gromacs/utility/cstringutil.h"
69 #include "gromacs/utility/fatalerror.h"
70 #include "gromacs/utility/real.h"
71 #include "gromacs/utility/smalloc.h"
73 #include "nbnxm_cuda.h"
78 /* This is a heuristically determined parameter for the Kepler
79 * and Maxwell architectures for the minimum size of ci lists by multiplying
80 * this constant with the # of multiprocessors on the current device.
81 * Since the maximum number of blocks per multiprocessor is 16, the ideal
82 * count for small systems is 32 or 48 blocks per multiprocessor. Because
83 * there is a bit of fluctuations in the generated block counts, we use
84 * a target of 44 instead of the ideal value of 48.
86 static unsigned int gpu_min_ci_balanced_factor = 44;
89 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
92 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam);
94 /*! \brief Return whether combination rules are used.
96 * \param[in] pointer to nonbonded paramter struct
97 * \return true if combination rules are used in this run, false otherwise
99 static inline bool useLjCombRule(const cu_nbparam_t *nbparam)
101 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
102 nbparam->vdwtype == evdwCuCUTCOMBLB);
105 /*! \brief Initialized the Ewald Coulomb correction GPU table.
107 Tabulates the Ewald Coulomb force and initializes the size/scale
108 and the table GPU array. If called with an already allocated table,
109 it just re-uploads the table.
111 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
114 if (nbp->coulomb_tab != nullptr)
116 nbnxn_cuda_free_nbparam_table(nbp);
119 nbp->coulomb_tab_scale = ic->tabq_scale;
120 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
121 ic->tabq_coul_F, ic->tabq_size);
125 /*! Initializes the atomdata structure first time, it only gets filled at
127 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
132 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
133 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
134 ad->bShiftVecUploaded = false;
136 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
137 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
139 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
140 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
141 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
142 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
144 /* initialize to nullptr poiters to data that is not allocated here and will
145 need reallocation in nbnxn_cuda_init_atomdata */
149 /* size -1 indicates that the respective array hasn't been initialized yet */
154 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
155 earlier GPUs, single or twin cut-off. */
156 static int pick_ewald_kernel_type(bool bTwinCut)
158 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
161 /* Benchmarking/development environment variables to force the use of
162 analytical or tabulated Ewald kernel. */
163 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
164 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
166 if (bForceAnalyticalEwald && bForceTabulatedEwald)
168 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
169 "requested through environment variables.");
172 /* By default use analytical Ewald. */
173 bUseAnalyticalEwald = true;
174 if (bForceAnalyticalEwald)
178 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
181 else if (bForceTabulatedEwald)
183 bUseAnalyticalEwald = false;
187 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
191 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
192 forces it (use it for debugging/benchmarking only). */
193 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
195 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
199 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
205 /*! Copies all parameters related to the cut-off from ic to nbp */
206 static void set_cutoff_parameters(cu_nbparam_t *nbp,
207 const interaction_const_t *ic,
208 const PairlistParams &listParams)
210 nbp->ewald_beta = ic->ewaldcoeff_q;
211 nbp->sh_ewald = ic->sh_ewald;
212 nbp->epsfac = ic->epsfac;
213 nbp->two_k_rf = 2.0 * ic->k_rf;
214 nbp->c_rf = ic->c_rf;
215 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
216 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
217 nbp->rlistOuter_sq = listParams.rlistOuter * listParams.rlistOuter;
218 nbp->rlistInner_sq = listParams.rlistInner * listParams.rlistInner;
219 nbp->useDynamicPruning = listParams.useDynamicPruning;
221 nbp->sh_lj_ewald = ic->sh_lj_ewald;
222 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
224 nbp->rvdw_switch = ic->rvdw_switch;
225 nbp->dispersion_shift = ic->dispersion_shift;
226 nbp->repulsion_shift = ic->repulsion_shift;
227 nbp->vdw_switch = ic->vdw_switch;
230 /*! Initializes the nonbonded parameter data structure. */
231 static void init_nbparam(cu_nbparam_t *nbp,
232 const interaction_const_t *ic,
233 const PairlistParams &listParams,
234 const nbnxn_atomdata_t::Params &nbatParams)
238 ntypes = nbatParams.numTypes;
240 set_cutoff_parameters(nbp, ic, listParams);
242 /* The kernel code supports LJ combination rules (geometric and LB) for
243 * all kernel types, but we only generate useful combination rule kernels.
244 * We currently only use LJ combination rule (geometric and LB) kernels
245 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
246 * with PME and 20% with RF, the other kernels speed up about half as much.
247 * For LJ force-switch the geometric rule would give 7% speed-up, but this
248 * combination is rarely used. LJ force-switch with LB rule is more common,
249 * but gives only 1% speed-up.
251 if (ic->vdwtype == evdwCUT)
253 switch (ic->vdw_modifier)
256 case eintmodPOTSHIFT:
257 switch (nbatParams.comb_rule)
260 nbp->vdwtype = evdwCuCUT;
263 nbp->vdwtype = evdwCuCUTCOMBGEOM;
266 nbp->vdwtype = evdwCuCUTCOMBLB;
269 gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
272 case eintmodFORCESWITCH:
273 nbp->vdwtype = evdwCuFSWITCH;
275 case eintmodPOTSWITCH:
276 nbp->vdwtype = evdwCuPSWITCH;
279 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
282 else if (ic->vdwtype == evdwPME)
284 if (ic->ljpme_comb_rule == ljcrGEOM)
286 assert(nbatParams.comb_rule == ljcrGEOM);
287 nbp->vdwtype = evdwCuEWALDGEOM;
291 assert(nbatParams.comb_rule == ljcrLB);
292 nbp->vdwtype = evdwCuEWALDLB;
297 gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
300 if (ic->eeltype == eelCUT)
302 nbp->eeltype = eelCuCUT;
304 else if (EEL_RF(ic->eeltype))
306 nbp->eeltype = eelCuRF;
308 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
310 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
311 nbp->eeltype = pick_ewald_kernel_type(false);
315 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
316 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
319 /* generate table for PME */
320 nbp->coulomb_tab = nullptr;
321 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
323 init_ewald_coulomb_force_table(ic, nbp);
326 /* set up LJ parameter lookup table */
327 if (!useLjCombRule(nbp))
329 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
330 nbatParams.nbfp.data(), 2*ntypes*ntypes);
333 /* set up LJ-PME parameter lookup table */
334 if (ic->vdwtype == evdwPME)
336 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
337 nbatParams.nbfp_comb.data(), 2*ntypes);
341 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
342 * electrostatic type switching to twin cut-off (or back) if needed. */
343 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t *nbv,
344 const interaction_const_t *ic)
346 if (!nbv || !nbv->useGpu())
350 cu_nbparam_t *nbp = nbv->gpu_nbv->nbparam;
352 set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
354 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw);
356 init_ewald_coulomb_force_table(ic, nbp);
359 /*! Initializes the pair list data structure. */
360 static void init_plist(cu_plist_t *pl)
362 /* initialize to nullptr pointers to data that is not allocated here and will
363 need reallocation in nbnxn_gpu_init_pairlist */
369 /* size -1 indicates that the respective array hasn't been initialized yet */
376 pl->imask_nalloc = -1;
378 pl->excl_nalloc = -1;
379 pl->haveFreshList = false;
382 /*! Initializes the timings data structure. */
383 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
392 for (i = 0; i < 2; i++)
394 for (j = 0; j < 2; j++)
396 t->ktime[i][j].t = 0.0;
397 t->ktime[i][j].c = 0;
401 t->pruneTime.t = 0.0;
402 t->dynamicPruneTime.c = 0;
403 t->dynamicPruneTime.t = 0.0;
406 /*! Initializes simulation constant data. */
407 static void cuda_init_const(gmx_nbnxn_cuda_t *nb,
408 const interaction_const_t *ic,
409 const PairlistParams &listParams,
410 const nbnxn_atomdata_t::Params &nbatParams)
412 init_atomdata_first(nb->atdat, nbatParams.numTypes);
413 init_nbparam(nb->nbparam, ic, listParams, nbatParams);
415 /* clear energy and shift force outputs */
416 nbnxn_cuda_clear_e_fshift(nb);
420 gpu_init(const gmx_device_info_t *deviceInfo,
421 const interaction_const_t *ic,
422 const PairlistParams &listParams,
423 const nbnxn_atomdata_t *nbat,
425 gmx_bool bLocalAndNonlocal)
429 gmx_nbnxn_cuda_t *nb;
432 snew(nb->nbparam, 1);
433 snew(nb->plist[InteractionLocality::Local], 1);
434 if (bLocalAndNonlocal)
436 snew(nb->plist[InteractionLocality::NonLocal], 1);
439 nb->bUseTwoStreams = bLocalAndNonlocal;
441 nb->timers = new cu_timers_t();
442 snew(nb->timings, 1);
445 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
446 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
447 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
449 init_plist(nb->plist[InteractionLocality::Local]);
451 /* set device info, just point it to the right GPU among the detected ones */
452 nb->dev_info = deviceInfo;
454 /* local/non-local GPU streams */
455 stat = cudaStreamCreate(&nb->stream[InteractionLocality::Local]);
456 CU_RET_ERR(stat, "cudaStreamCreate on stream[InterationLocality::Local] failed");
457 if (nb->bUseTwoStreams)
459 init_plist(nb->plist[InteractionLocality::NonLocal]);
461 /* Note that the device we're running on does not have to support
462 * priorities, because we are querying the priority range which in this
463 * case will be a single value.
465 int highest_priority;
466 stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
467 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
469 stat = cudaStreamCreateWithPriority(&nb->stream[InteractionLocality::NonLocal],
472 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[InteractionLocality::NonLocal] failed");
475 /* init events for sychronization (timing disabled for performance reasons!) */
476 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
477 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
478 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
479 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
481 /* WARNING: CUDA timings are incorrect with multiple streams.
482 * This is the main reason why they are disabled by default.
484 // TODO: Consider turning on by default when we can detect nr of streams.
485 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
489 init_timings(nb->timings);
492 /* set the kernel type for the current GPU */
493 /* pick L1 cache configuration */
494 cuda_set_cacheconfig();
496 cuda_init_const(nb, ic, listParams, nbat->params());
500 fprintf(debug, "Initialized CUDA data structures.\n");
506 void gpu_init_pairlist(gmx_nbnxn_cuda_t *nb,
507 const NbnxnPairlistGpu *h_plist,
508 const InteractionLocality iloc)
511 bool bDoTime = (nb->bDoTime && !h_plist->sci.empty());
512 cudaStream_t stream = nb->stream[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(stream);
534 iTimers.didPairlistH2D = true;
537 Context context = nullptr;
539 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(),
540 &d_plist->nsci, &d_plist->sci_nalloc, context);
541 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(),
542 stream, GpuApiCallBehavior::Async,
543 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
545 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(),
546 &d_plist->ncj4, &d_plist->cj4_nalloc, context);
547 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(),
548 stream, GpuApiCallBehavior::Async,
549 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
551 reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size()*c_nbnxnGpuClusterpairSplit,
552 &d_plist->nimask, &d_plist->imask_nalloc, context);
554 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(),
555 &d_plist->nexcl, &d_plist->excl_nalloc, context);
556 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(),
557 stream, GpuApiCallBehavior::Async,
558 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
562 iTimers.pl_h2d.closeTimingRegion(stream);
565 /* the next use of thist list we be the first one, so we need to prune */
566 d_plist->haveFreshList = true;
569 void gpu_upload_shiftvec(gmx_nbnxn_cuda_t *nb,
570 const nbnxn_atomdata_t *nbatom)
572 cu_atomdata_t *adat = nb->atdat;
573 cudaStream_t ls = nb->stream[InteractionLocality::Local];
575 /* only if we have a dynamic box */
576 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
578 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(),
579 SHIFTS * sizeof(*adat->shift_vec), ls);
580 adat->bShiftVecUploaded = true;
584 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
585 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
588 cu_atomdata_t *adat = nb->atdat;
589 cudaStream_t ls = nb->stream[InteractionLocality::Local];
591 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
592 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
595 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
596 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
599 cu_atomdata_t *adat = nb->atdat;
600 cudaStream_t ls = nb->stream[InteractionLocality::Local];
602 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
603 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
604 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
605 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
606 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
607 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
610 void gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
612 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
613 /* clear shift force array and energies if the outputs were
614 used in the current step */
615 if (flags & GMX_FORCE_VIRIAL)
617 nbnxn_cuda_clear_e_fshift(nb);
621 void gpu_init_atomdata(gmx_nbnxn_cuda_t *nb,
622 const nbnxn_atomdata_t *nbat)
627 bool bDoTime = nb->bDoTime;
628 cu_timers_t *timers = nb->timers;
629 cu_atomdata_t *d_atdat = nb->atdat;
630 cudaStream_t ls = nb->stream[InteractionLocality::Local];
632 natoms = nbat->numAtoms();
637 /* time async copy */
638 timers->atdat.openTimingRegion(ls);
641 /* need to reallocate if we have to copy more atoms than the amount of space
642 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
643 if (natoms > d_atdat->nalloc)
645 nalloc = over_alloc_small(natoms);
647 /* free up first if the arrays have already been initialized */
648 if (d_atdat->nalloc != -1)
650 freeDeviceBuffer(&d_atdat->f);
651 freeDeviceBuffer(&d_atdat->xq);
652 freeDeviceBuffer(&d_atdat->atom_types);
653 freeDeviceBuffer(&d_atdat->lj_comb);
656 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
657 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
658 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
659 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
660 if (useLjCombRule(nb->nbparam))
662 stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
663 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
667 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
668 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
671 d_atdat->nalloc = nalloc;
675 d_atdat->natoms = natoms;
676 d_atdat->natoms_local = nbat->natoms_local;
678 /* need to clear GPU f output if realloc happened */
681 nbnxn_cuda_clear_f(nb, nalloc);
684 if (useLjCombRule(nb->nbparam))
686 cu_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(),
687 natoms*sizeof(*d_atdat->lj_comb), ls);
691 cu_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(),
692 natoms*sizeof(*d_atdat->atom_types), ls);
697 timers->atdat.closeTimingRegion(ls);
701 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam)
703 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
705 destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
709 void gpu_free(gmx_nbnxn_cuda_t *nb)
712 cu_atomdata_t *atdat;
713 cu_nbparam_t *nbparam;
721 nbparam = nb->nbparam;
723 nbnxn_cuda_free_nbparam_table(nbparam);
725 stat = cudaEventDestroy(nb->nonlocal_done);
726 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
727 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
728 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
733 /* The non-local counters/stream (second in the array) are needed only with DD. */
734 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
736 stat = cudaStreamDestroy(nb->stream[i]);
737 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
741 if (!useLjCombRule(nb->nbparam))
743 destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj);
747 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
749 destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
752 stat = cudaFree(atdat->shift_vec);
753 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
754 stat = cudaFree(atdat->fshift);
755 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
757 stat = cudaFree(atdat->e_lj);
758 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
759 stat = cudaFree(atdat->e_el);
760 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
762 freeDeviceBuffer(&atdat->f);
763 freeDeviceBuffer(&atdat->xq);
764 freeDeviceBuffer(&atdat->atom_types);
765 freeDeviceBuffer(&atdat->lj_comb);
768 auto *plist = nb->plist[InteractionLocality::Local];
769 freeDeviceBuffer(&plist->sci);
770 freeDeviceBuffer(&plist->cj4);
771 freeDeviceBuffer(&plist->imask);
772 freeDeviceBuffer(&plist->excl);
774 if (nb->bUseTwoStreams)
776 auto *plist_nl = nb->plist[InteractionLocality::NonLocal];
777 freeDeviceBuffer(&plist_nl->sci);
778 freeDeviceBuffer(&plist_nl->cj4);
779 freeDeviceBuffer(&plist_nl->imask);
780 freeDeviceBuffer(&plist_nl->excl);
785 pfree(nb->nbst.e_lj);
786 nb->nbst.e_lj = nullptr;
788 pfree(nb->nbst.e_el);
789 nb->nbst.e_el = nullptr;
791 pfree(nb->nbst.fshift);
792 nb->nbst.fshift = nullptr;
801 fprintf(debug, "Cleaned up CUDA data structures.\n");
805 //! This function is documented in the header file
806 gmx_wallclock_gpu_nbnxn_t *gpu_get_timings(gmx_nbnxn_cuda_t *nb)
808 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
811 void gpu_reset_timings(nonbonded_verlet_t* nbv)
813 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
815 init_timings(nbv->gpu_nbv->timings);
819 int gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
821 return nb != nullptr ?
822 gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
826 gmx_bool gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
828 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
829 (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
832 void *gpu_get_command_stream(gmx_nbnxn_gpu_t *nb,
833 const InteractionLocality iloc)
837 return static_cast<void *>(&nb->stream[iloc]);
840 void *gpu_get_xq(gmx_nbnxn_gpu_t *nb)
844 return static_cast<void *>(nb->atdat->xq);
847 void *gpu_get_f(gmx_nbnxn_gpu_t *nb)
851 return static_cast<void *>(nb->atdat->f);
854 rvec *gpu_get_fshift(gmx_nbnxn_gpu_t *nb)
858 return reinterpret_cast<rvec *>(nb->atdat->fshift);
861 /* Initialization for X buffer operations on GPU. */
862 /* TODO Remove explicit pinning from host arrays from here and manage in a more natural way*/
863 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet &gridSet,
864 gmx_nbnxn_gpu_t *gpu_nbv,
865 const Nbnxm::AtomLocality locality)
868 const Nbnxm::InteractionLocality iloc = ((locality == AtomLocality::Local) ?
869 InteractionLocality::Local : InteractionLocality::NonLocal);
870 cudaStream_t stream = gpu_nbv->stream[iloc];
871 bool bDoTime = gpu_nbv->bDoTime;
872 int gridBegin = 0, gridEnd = 0;
877 case Nbnxm::AtomLocality::All:
879 gridEnd = gridSet.grids().size();
881 case Nbnxm::AtomLocality::Local:
885 case Nbnxm::AtomLocality::NonLocal:
887 gridEnd = gridSet.grids().size();
889 case Nbnxm::AtomLocality::Count:
890 GMX_ASSERT(false, "Count is invalid locality specifier");
896 for (int g = gridBegin; g < gridEnd; g++)
899 const Nbnxm::Grid &grid = gridSet.grids()[g];
901 //TODO improve naming here. Either use the getters straight or
902 //using variables with about the same name as the getters
903 const int ncxy = grid.numColumns();
904 const int *a = gridSet.atomIndices().data();
905 const int a_nalloc = gridSet.atomIndices().size();
906 const int *na_all = grid.cxy_na().data();
907 const int *cxy_ind = grid.cxy_ind().data();
908 const int natoms_nonlocal = gridSet.numRealAtomsTotal();
910 if (iloc == Nbnxm::InteractionLocality::Local)
915 freeDeviceBuffer(&gpu_nbv->xrvec);
917 //TODO replace with reallocateDeviceBuffer
918 allocateDeviceBuffer(&gpu_nbv->xrvec, natoms_nonlocal, nullptr);
920 if (gpu_nbv->abufops)
922 freeDeviceBuffer(&gpu_nbv->abufops);
924 //TODO replace with reallocateDeviceBuffer
925 allocateDeviceBuffer(&gpu_nbv->abufops, a_nalloc, nullptr);
929 // source data must be pinned for H2D assertion. This should be moved into place where data is (re-)alloced.
930 stat = cudaHostRegister((void*) a, a_nalloc*sizeof(int), cudaHostRegisterDefault);
931 CU_RET_ERR(stat, "cudaHostRegister failed on a");
935 gpu_nbv->timers->xf[locality].nb_h2d.openTimingRegion(stream);
938 copyToDeviceBuffer(&gpu_nbv->abufops, a, 0, a_nalloc, stream, GpuApiCallBehavior::Async, nullptr);
942 gpu_nbv->timers->xf[locality].nb_h2d.closeTimingRegion(stream);
945 stat = cudaHostUnregister((void*) a);
946 CU_RET_ERR(stat, "cudaHostUnRegister failed on a");
951 if (gpu_nbv->nabufops[locality])
953 freeDeviceBuffer(&gpu_nbv->nabufops[locality]);
956 //TODO replace with reallocateDeviceBuffer
957 allocateDeviceBuffer(&gpu_nbv->nabufops[locality], ncxy, nullptr);
959 if (gpu_nbv->cxybufops[locality])
961 freeDeviceBuffer(&gpu_nbv->cxybufops[locality]);
963 //TODO replace with reallocateDeviceBuffer
964 allocateDeviceBuffer(&gpu_nbv->cxybufops[locality], ncxy, nullptr);
968 // source data must be pinned for H2D assertion. This should be moved into place where data is (re-)alloced.
969 stat = cudaHostRegister((void*) na_all, ncxy*sizeof(int), cudaHostRegisterDefault);
970 CU_RET_ERR(stat, "cudaHostRegister failed on na_all");
974 gpu_nbv->timers->xf[locality].nb_h2d.openTimingRegion(stream);
977 copyToDeviceBuffer(&gpu_nbv->nabufops[locality], na_all, 0, ncxy, stream, GpuApiCallBehavior::Async, nullptr);
981 gpu_nbv->timers->xf[locality].nb_h2d.closeTimingRegion(stream);
984 stat = cudaHostUnregister((void*) na_all);
985 CU_RET_ERR(stat, "cudaHostUnRegister failed on na_all");
987 // source data must be pinned for H2D assertion. This should be moved into place where data is (re-)alloced.
988 stat = cudaHostRegister((void*) cxy_ind, ncxy*sizeof(int), cudaHostRegisterDefault);
989 CU_RET_ERR(stat, "cudaHostRegister failed on cxy_ind");
993 gpu_nbv->timers->xf[locality].nb_h2d.openTimingRegion(stream);
996 copyToDeviceBuffer(&gpu_nbv->cxybufops[locality], cxy_ind, 0, ncxy, stream, GpuApiCallBehavior::Async, nullptr);
1000 gpu_nbv->timers->xf[locality].nb_h2d.closeTimingRegion(stream);
1003 stat = cudaHostUnregister((void*) cxy_ind);
1004 CU_RET_ERR(stat, "cudaHostUnRegister failed on cxy_ind");
1010 } // namespace Nbnxm