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/gpu_data_mgmt.h"
61 #include "gromacs/nbnxm/nbnxm.h"
62 #include "gromacs/pbcutil/ishift.h"
63 #include "gromacs/timing/gpu_timing.h"
64 #include "gromacs/utility/basedefinitions.h"
65 #include "gromacs/utility/cstringutil.h"
66 #include "gromacs/utility/fatalerror.h"
67 #include "gromacs/utility/real.h"
68 #include "gromacs/utility/smalloc.h"
70 #include "nbnxm_cuda.h"
75 /* This is a heuristically determined parameter for the Kepler
76 * and Maxwell architectures for the minimum size of ci lists by multiplying
77 * this constant with the # of multiprocessors on the current device.
78 * Since the maximum number of blocks per multiprocessor is 16, the ideal
79 * count for small systems is 32 or 48 blocks per multiprocessor. Because
80 * there is a bit of fluctuations in the generated block counts, we use
81 * a target of 44 instead of the ideal value of 48.
83 static unsigned int gpu_min_ci_balanced_factor = 44;
86 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
89 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam);
91 /*! \brief Return whether combination rules are used.
93 * \param[in] pointer to nonbonded paramter struct
94 * \return true if combination rules are used in this run, false otherwise
96 static inline bool useLjCombRule(const cu_nbparam_t *nbparam)
98 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
99 nbparam->vdwtype == evdwCuCUTCOMBLB);
102 /*! \brief Initialized the Ewald Coulomb correction GPU table.
104 Tabulates the Ewald Coulomb force and initializes the size/scale
105 and the table GPU array. If called with an already allocated table,
106 it just re-uploads the table.
108 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
111 if (nbp->coulomb_tab != nullptr)
113 nbnxn_cuda_free_nbparam_table(nbp);
116 nbp->coulomb_tab_scale = ic->tabq_scale;
117 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
118 ic->tabq_coul_F, ic->tabq_size);
122 /*! Initializes the atomdata structure first time, it only gets filled at
124 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
129 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
130 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
131 ad->bShiftVecUploaded = false;
133 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
134 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
136 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
137 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
138 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
139 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
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(bool bTwinCut)
155 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
158 /* Benchmarking/development environment variables to force the use of
159 analytical or tabulated Ewald kernel. */
160 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
161 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
163 if (bForceAnalyticalEwald && bForceTabulatedEwald)
165 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
166 "requested through environment variables.");
169 /* By default use analytical Ewald. */
170 bUseAnalyticalEwald = true;
171 if (bForceAnalyticalEwald)
175 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
178 else if (bForceTabulatedEwald)
180 bUseAnalyticalEwald = false;
184 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
188 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
189 forces it (use it for debugging/benchmarking only). */
190 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
192 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
196 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
202 /*! Copies all parameters related to the cut-off from ic to nbp */
203 static void set_cutoff_parameters(cu_nbparam_t *nbp,
204 const interaction_const_t *ic,
205 const NbnxnListParameters *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 NbnxnListParameters *listParams,
231 const nbnxn_atomdata_t::Params &nbatParams)
235 ntypes = nbatParams.numTypes;
237 set_cutoff_parameters(nbp, ic, listParams);
239 /* The kernel code supports LJ combination rules (geometric and LB) for
240 * all kernel types, but we only generate useful combination rule kernels.
241 * We currently only use LJ combination rule (geometric and LB) kernels
242 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
243 * with PME and 20% with RF, the other kernels speed up about half as much.
244 * For LJ force-switch the geometric rule would give 7% speed-up, but this
245 * combination is rarely used. LJ force-switch with LB rule is more common,
246 * but gives only 1% speed-up.
248 if (ic->vdwtype == evdwCUT)
250 switch (ic->vdw_modifier)
253 case eintmodPOTSHIFT:
254 switch (nbatParams.comb_rule)
257 nbp->vdwtype = evdwCuCUT;
260 nbp->vdwtype = evdwCuCUTCOMBGEOM;
263 nbp->vdwtype = evdwCuCUTCOMBLB;
266 gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
269 case eintmodFORCESWITCH:
270 nbp->vdwtype = evdwCuFSWITCH;
272 case eintmodPOTSWITCH:
273 nbp->vdwtype = evdwCuPSWITCH;
276 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
279 else if (ic->vdwtype == evdwPME)
281 if (ic->ljpme_comb_rule == ljcrGEOM)
283 assert(nbatParams.comb_rule == ljcrGEOM);
284 nbp->vdwtype = evdwCuEWALDGEOM;
288 assert(nbatParams.comb_rule == ljcrLB);
289 nbp->vdwtype = evdwCuEWALDLB;
294 gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
297 if (ic->eeltype == eelCUT)
299 nbp->eeltype = eelCuCUT;
301 else if (EEL_RF(ic->eeltype))
303 nbp->eeltype = eelCuRF;
305 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
307 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
308 nbp->eeltype = pick_ewald_kernel_type(false);
312 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
313 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
316 /* generate table for PME */
317 nbp->coulomb_tab = nullptr;
318 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
320 init_ewald_coulomb_force_table(ic, nbp);
323 /* set up LJ parameter lookup table */
324 if (!useLjCombRule(nbp))
326 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
327 nbatParams.nbfp.data(), 2*ntypes*ntypes);
330 /* set up LJ-PME parameter lookup table */
331 if (ic->vdwtype == evdwPME)
333 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
334 nbatParams.nbfp_comb.data(), 2*ntypes);
338 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
339 * electrostatic type switching to twin cut-off (or back) if needed. */
340 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t *nbv,
341 const interaction_const_t *ic,
342 const NbnxnListParameters *listParams)
344 if (!nbv || !nbv->useGpu())
348 cu_nbparam_t *nbp = nbv->gpu_nbv->nbparam;
350 set_cutoff_parameters(nbp, ic, listParams);
352 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw);
354 init_ewald_coulomb_force_table(ic, nbp);
357 /*! Initializes the pair list data structure. */
358 static void init_plist(cu_plist_t *pl)
360 /* initialize to nullptr pointers to data that is not allocated here and will
361 need reallocation in nbnxn_gpu_init_pairlist */
367 /* size -1 indicates that the respective array hasn't been initialized yet */
374 pl->imask_nalloc = -1;
376 pl->excl_nalloc = -1;
377 pl->haveFreshList = false;
380 /*! Initializes the timings data structure. */
381 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
390 for (i = 0; i < 2; i++)
392 for (j = 0; j < 2; j++)
394 t->ktime[i][j].t = 0.0;
395 t->ktime[i][j].c = 0;
399 t->pruneTime.t = 0.0;
400 t->dynamicPruneTime.c = 0;
401 t->dynamicPruneTime.t = 0.0;
404 /*! Initializes simulation constant data. */
405 static void cuda_init_const(gmx_nbnxn_cuda_t *nb,
406 const interaction_const_t *ic,
407 const NbnxnListParameters *listParams,
408 const nbnxn_atomdata_t::Params &nbatParams)
410 init_atomdata_first(nb->atdat, nbatParams.numTypes);
411 init_nbparam(nb->nbparam, ic, listParams, nbatParams);
413 /* clear energy and shift force outputs */
414 nbnxn_cuda_clear_e_fshift(nb);
418 gpu_init(const gmx_device_info_t *deviceInfo,
419 const interaction_const_t *ic,
420 const NbnxnListParameters *listParams,
421 const nbnxn_atomdata_t *nbat,
423 gmx_bool bLocalAndNonlocal)
427 gmx_nbnxn_cuda_t *nb;
430 snew(nb->nbparam, 1);
431 snew(nb->plist[InteractionLocality::Local], 1);
432 if (bLocalAndNonlocal)
434 snew(nb->plist[InteractionLocality::NonLocal], 1);
437 nb->bUseTwoStreams = bLocalAndNonlocal;
439 nb->timers = new cu_timers_t();
440 snew(nb->timings, 1);
443 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
444 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
445 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
447 init_plist(nb->plist[InteractionLocality::Local]);
449 /* set device info, just point it to the right GPU among the detected ones */
450 nb->dev_info = deviceInfo;
452 /* local/non-local GPU streams */
453 stat = cudaStreamCreate(&nb->stream[InteractionLocality::Local]);
454 CU_RET_ERR(stat, "cudaStreamCreate on stream[InterationLocality::Local] failed");
455 if (nb->bUseTwoStreams)
457 init_plist(nb->plist[InteractionLocality::NonLocal]);
459 /* Note that the device we're running on does not have to support
460 * priorities, because we are querying the priority range which in this
461 * case will be a single value.
463 int highest_priority;
464 stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
465 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
467 stat = cudaStreamCreateWithPriority(&nb->stream[InteractionLocality::NonLocal],
470 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[InteractionLocality::NonLocal] failed");
473 /* init events for sychronization (timing disabled for performance reasons!) */
474 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
475 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
476 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
477 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
479 /* WARNING: CUDA timings are incorrect with multiple streams.
480 * This is the main reason why they are disabled by default.
482 // TODO: Consider turning on by default when we can detect nr of streams.
483 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
487 init_timings(nb->timings);
490 /* set the kernel type for the current GPU */
491 /* pick L1 cache configuration */
492 cuda_set_cacheconfig();
494 cuda_init_const(nb, ic, listParams, nbat->params());
498 fprintf(debug, "Initialized CUDA data structures.\n");
504 void gpu_init_pairlist(gmx_nbnxn_cuda_t *nb,
505 const NbnxnPairlistGpu *h_plist,
506 const InteractionLocality iloc)
509 bool bDoTime = (nb->bDoTime && !h_plist->sci.empty());
510 cudaStream_t stream = nb->stream[iloc];
511 cu_plist_t *d_plist = nb->plist[iloc];
513 if (d_plist->na_c < 0)
515 d_plist->na_c = h_plist->na_ci;
519 if (d_plist->na_c != h_plist->na_ci)
521 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
522 d_plist->na_c, h_plist->na_ci);
527 gpu_timers_t::Interaction &iTimers = nb->timers->interaction[iloc];
531 iTimers.pl_h2d.openTimingRegion(stream);
532 iTimers.didPairlistH2D = true;
535 Context context = nullptr;
537 reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(),
538 &d_plist->nsci, &d_plist->sci_nalloc, context);
539 copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(),
540 stream, GpuApiCallBehavior::Async,
541 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
543 reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(),
544 &d_plist->ncj4, &d_plist->cj4_nalloc, context);
545 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(),
546 stream, GpuApiCallBehavior::Async,
547 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, context);
552 reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(),
553 &d_plist->nexcl, &d_plist->excl_nalloc, context);
554 copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(),
555 stream, GpuApiCallBehavior::Async,
556 bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
560 iTimers.pl_h2d.closeTimingRegion(stream);
563 /* the next use of thist list we be the first one, so we need to prune */
564 d_plist->haveFreshList = true;
567 void gpu_upload_shiftvec(gmx_nbnxn_cuda_t *nb,
568 const nbnxn_atomdata_t *nbatom)
570 cu_atomdata_t *adat = nb->atdat;
571 cudaStream_t ls = nb->stream[InteractionLocality::Local];
573 /* only if we have a dynamic box */
574 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
576 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(),
577 SHIFTS * sizeof(*adat->shift_vec), ls);
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(gmx_nbnxn_cuda_t *nb, int natoms_clear)
586 cu_atomdata_t *adat = nb->atdat;
587 cudaStream_t ls = nb->stream[InteractionLocality::Local];
589 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
590 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
593 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
594 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
597 cu_atomdata_t *adat = nb->atdat;
598 cudaStream_t ls = nb->stream[InteractionLocality::Local];
600 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
601 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
602 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
603 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
604 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
605 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
608 void gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
610 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
611 /* clear shift force array and energies if the outputs were
612 used in the current step */
613 if (flags & GMX_FORCE_VIRIAL)
615 nbnxn_cuda_clear_e_fshift(nb);
619 void gpu_init_atomdata(gmx_nbnxn_cuda_t *nb,
620 const nbnxn_atomdata_t *nbat)
625 bool bDoTime = nb->bDoTime;
626 cu_timers_t *timers = nb->timers;
627 cu_atomdata_t *d_atdat = nb->atdat;
628 cudaStream_t ls = nb->stream[InteractionLocality::Local];
630 natoms = nbat->numAtoms();
635 /* time async copy */
636 timers->atdat.openTimingRegion(ls);
639 /* need to reallocate if we have to copy more atoms than the amount of space
640 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
641 if (natoms > d_atdat->nalloc)
643 nalloc = over_alloc_small(natoms);
645 /* free up first if the arrays have already been initialized */
646 if (d_atdat->nalloc != -1)
648 freeDeviceBuffer(&d_atdat->f);
649 freeDeviceBuffer(&d_atdat->xq);
650 freeDeviceBuffer(&d_atdat->atom_types);
651 freeDeviceBuffer(&d_atdat->lj_comb);
654 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
655 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
656 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
657 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
658 if (useLjCombRule(nb->nbparam))
660 stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
661 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
665 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
666 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
669 d_atdat->nalloc = nalloc;
673 d_atdat->natoms = natoms;
674 d_atdat->natoms_local = nbat->natoms_local;
676 /* need to clear GPU f output if realloc happened */
679 nbnxn_cuda_clear_f(nb, nalloc);
682 if (useLjCombRule(nb->nbparam))
684 cu_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(),
685 natoms*sizeof(*d_atdat->lj_comb), ls);
689 cu_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(),
690 natoms*sizeof(*d_atdat->atom_types), ls);
695 timers->atdat.closeTimingRegion(ls);
699 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam)
701 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
703 destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
707 void gpu_free(gmx_nbnxn_cuda_t *nb)
710 cu_atomdata_t *atdat;
711 cu_nbparam_t *nbparam;
719 nbparam = nb->nbparam;
721 nbnxn_cuda_free_nbparam_table(nbparam);
723 stat = cudaEventDestroy(nb->nonlocal_done);
724 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
725 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
726 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
731 /* The non-local counters/stream (second in the array) are needed only with DD. */
732 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
734 stat = cudaStreamDestroy(nb->stream[i]);
735 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
739 if (!useLjCombRule(nb->nbparam))
741 destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj);
745 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
747 destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
750 stat = cudaFree(atdat->shift_vec);
751 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
752 stat = cudaFree(atdat->fshift);
753 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
755 stat = cudaFree(atdat->e_lj);
756 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
757 stat = cudaFree(atdat->e_el);
758 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
760 freeDeviceBuffer(&atdat->f);
761 freeDeviceBuffer(&atdat->xq);
762 freeDeviceBuffer(&atdat->atom_types);
763 freeDeviceBuffer(&atdat->lj_comb);
766 auto *plist = nb->plist[InteractionLocality::Local];
767 freeDeviceBuffer(&plist->sci);
768 freeDeviceBuffer(&plist->cj4);
769 freeDeviceBuffer(&plist->imask);
770 freeDeviceBuffer(&plist->excl);
772 if (nb->bUseTwoStreams)
774 auto *plist_nl = nb->plist[InteractionLocality::NonLocal];
775 freeDeviceBuffer(&plist_nl->sci);
776 freeDeviceBuffer(&plist_nl->cj4);
777 freeDeviceBuffer(&plist_nl->imask);
778 freeDeviceBuffer(&plist_nl->excl);
783 pfree(nb->nbst.e_lj);
784 nb->nbst.e_lj = nullptr;
786 pfree(nb->nbst.e_el);
787 nb->nbst.e_el = nullptr;
789 pfree(nb->nbst.fshift);
790 nb->nbst.fshift = nullptr;
799 fprintf(debug, "Cleaned up CUDA data structures.\n");
803 //! This function is documented in the header file
804 gmx_wallclock_gpu_nbnxn_t *gpu_get_timings(gmx_nbnxn_cuda_t *nb)
806 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
809 void gpu_reset_timings(nonbonded_verlet_t* nbv)
811 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
813 init_timings(nbv->gpu_nbv->timings);
817 int gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
819 return nb != nullptr ?
820 gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
824 gmx_bool gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
826 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
827 (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
830 void *gpu_get_command_stream(gmx_nbnxn_gpu_t *nb,
831 const InteractionLocality iloc)
835 return static_cast<void *>(&nb->stream[iloc]);
838 void *gpu_get_xq(gmx_nbnxn_gpu_t *nb)
842 return static_cast<void *>(nb->atdat->xq);
845 void *gpu_get_f(gmx_nbnxn_gpu_t *nb)
849 return static_cast<void *>(nb->atdat->f);
852 rvec *gpu_get_fshift(gmx_nbnxn_gpu_t *nb)
856 return reinterpret_cast<rvec *>(nb->atdat->fshift);