2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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 #include "gromacs/gpu_utils/cudautils.cuh"
48 #include "gromacs/gpu_utils/gpu_utils.h"
49 #include "gromacs/gpu_utils/pmalloc_cuda.h"
50 #include "gromacs/hardware/gpu_hw_info.h"
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/mdlib/force_flags.h"
53 #include "gromacs/mdlib/nb_verlet.h"
54 #include "gromacs/mdlib/nbnxn_consts.h"
55 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/timing/gpu_timing.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/real.h"
64 #include "gromacs/utility/smalloc.h"
66 #include "nbnxn_cuda.h"
67 #include "nbnxn_cuda_types.h"
69 /* This is a heuristically determined parameter for the Fermi, Kepler
70 * and Maxwell architectures for the minimum size of ci lists by multiplying
71 * this constant with the # of multiprocessors on the current device.
72 * Since the maximum number of blocks per multiprocessor is 16, the ideal
73 * count for small systems is 32 or 48 blocks per multiprocessor. Because
74 * there is a bit of fluctuations in the generated block counts, we use
75 * a target of 44 instead of the ideal value of 48.
77 static unsigned int gpu_min_ci_balanced_factor = 44;
80 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
83 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam,
84 const gmx_device_info_t *dev_info);
86 /*! \brief Return whether combination rules are used.
88 * \param[in] pointer to nonbonded paramter struct
89 * \return true if combination rules are used in this run, false otherwise
91 static inline bool useLjCombRule(const cu_nbparam_t *nbparam)
93 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
94 nbparam->vdwtype == evdwCuCUTCOMBLB);
97 /*! \brief Initialized the Ewald Coulomb correction GPU table.
99 Tabulates the Ewald Coulomb force and initializes the size/scale
100 and the table GPU array. If called with an already allocated table,
101 it just re-uploads the table.
103 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
105 const gmx_device_info_t *dev_info)
107 if (nbp->coulomb_tab != nullptr)
109 nbnxn_cuda_free_nbparam_table(nbp, dev_info);
112 nbp->coulomb_tab_scale = ic->tabq_scale;
113 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
114 ic->tabq_coul_F, ic->tabq_size, dev_info);
118 /*! Initializes the atomdata structure first time, it only gets filled at
120 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
125 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
126 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
127 ad->bShiftVecUploaded = false;
129 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
130 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
132 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
133 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
134 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
135 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
137 /* initialize to nullptr poiters to data that is not allocated here and will
138 need reallocation in nbnxn_cuda_init_atomdata */
142 /* size -1 indicates that the respective array hasn't been initialized yet */
147 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
148 earlier GPUs, single or twin cut-off. */
149 static int pick_ewald_kernel_type(bool bTwinCut,
150 const gmx_device_info_t *dev_info)
152 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
155 /* Benchmarking/development environment variables to force the use of
156 analytical or tabulated Ewald kernel. */
157 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
158 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
160 if (bForceAnalyticalEwald && bForceTabulatedEwald)
162 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
163 "requested through environment variables.");
166 /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
167 if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
169 bUseAnalyticalEwald = true;
173 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
178 bUseAnalyticalEwald = false;
182 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
186 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
187 forces it (use it for debugging/benchmarking only). */
188 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
190 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
194 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
200 /*! Copies all parameters related to the cut-off from ic to nbp */
201 static void set_cutoff_parameters(cu_nbparam_t *nbp,
202 const interaction_const_t *ic,
203 const NbnxnListParameters *listParams)
205 nbp->ewald_beta = ic->ewaldcoeff_q;
206 nbp->sh_ewald = ic->sh_ewald;
207 nbp->epsfac = ic->epsfac;
208 nbp->two_k_rf = 2.0 * ic->k_rf;
209 nbp->c_rf = ic->c_rf;
210 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
211 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
212 nbp->rlistOuter_sq = listParams->rlistOuter * listParams->rlistOuter;
213 nbp->rlistInner_sq = listParams->rlistInner * listParams->rlistInner;
214 nbp->useDynamicPruning = listParams->useDynamicPruning;
216 nbp->sh_lj_ewald = ic->sh_lj_ewald;
217 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
219 nbp->rvdw_switch = ic->rvdw_switch;
220 nbp->dispersion_shift = ic->dispersion_shift;
221 nbp->repulsion_shift = ic->repulsion_shift;
222 nbp->vdw_switch = ic->vdw_switch;
225 /*! Initializes the nonbonded parameter data structure. */
226 static void init_nbparam(cu_nbparam_t *nbp,
227 const interaction_const_t *ic,
228 const NbnxnListParameters *listParams,
229 const nbnxn_atomdata_t *nbat,
230 const gmx_device_info_t *dev_info)
234 ntypes = nbat->ntype;
236 set_cutoff_parameters(nbp, ic, listParams);
238 /* The kernel code supports LJ combination rules (geometric and LB) for
239 * all kernel types, but we only generate useful combination rule kernels.
240 * We currently only use LJ combination rule (geometric and LB) kernels
241 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
242 * with PME and 20% with RF, the other kernels speed up about half as much.
243 * For LJ force-switch the geometric rule would give 7% speed-up, but this
244 * combination is rarely used. LJ force-switch with LB rule is more common,
245 * but gives only 1% speed-up.
247 if (ic->vdwtype == evdwCUT)
249 switch (ic->vdw_modifier)
252 case eintmodPOTSHIFT:
253 switch (nbat->comb_rule)
256 nbp->vdwtype = evdwCuCUT;
259 nbp->vdwtype = evdwCuCUTCOMBGEOM;
262 nbp->vdwtype = evdwCuCUTCOMBLB;
265 gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
268 case eintmodFORCESWITCH:
269 nbp->vdwtype = evdwCuFSWITCH;
271 case eintmodPOTSWITCH:
272 nbp->vdwtype = evdwCuPSWITCH;
275 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
278 else if (ic->vdwtype == evdwPME)
280 if (ic->ljpme_comb_rule == ljcrGEOM)
282 assert(nbat->comb_rule == ljcrGEOM);
283 nbp->vdwtype = evdwCuEWALDGEOM;
287 assert(nbat->comb_rule == ljcrLB);
288 nbp->vdwtype = evdwCuEWALDLB;
293 gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
296 if (ic->eeltype == eelCUT)
298 nbp->eeltype = eelCuCUT;
300 else if (EEL_RF(ic->eeltype))
302 nbp->eeltype = eelCuRF;
304 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
306 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
307 nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
311 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
312 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
315 /* generate table for PME */
316 nbp->coulomb_tab = nullptr;
317 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
319 init_ewald_coulomb_force_table(ic, nbp, dev_info);
322 /* set up LJ parameter lookup table */
323 if (!useLjCombRule(nbp))
325 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
326 nbat->nbfp, 2*ntypes*ntypes, dev_info);
329 /* set up LJ-PME parameter lookup table */
330 if (ic->vdwtype == evdwPME)
332 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
333 nbat->nbfp_comb, 2*ntypes, dev_info);
337 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
338 * electrostatic type switching to twin cut-off (or back) if needed. */
339 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t *nbv,
340 const interaction_const_t *ic,
341 const NbnxnListParameters *listParams)
343 if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
347 gmx_nbnxn_cuda_t *nb = nbv->gpu_nbv;
348 cu_nbparam_t *nbp = nb->nbparam;
350 set_cutoff_parameters(nbp, ic, listParams);
352 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
355 init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
358 /*! Initializes the pair list data structure. */
359 static void init_plist(cu_plist_t *pl)
361 /* initialize to nullptr pointers to data that is not allocated here and will
362 need reallocation in nbnxn_gpu_init_pairlist */
368 /* size -1 indicates that the respective array hasn't been initialized yet */
375 pl->imask_nalloc = -1;
377 pl->excl_nalloc = -1;
378 pl->haveFreshList = false;
381 /*! Initializes the timer data structure. */
382 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
384 /* The non-local counters/stream (second in the array) are needed only with DD. */
385 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
387 t->didPairlistH2D[i] = false;
388 t->didPrune[i] = false;
389 t->didRollingPrune[i] = false;
393 /*! Initializes the timings data structure. */
394 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
403 for (i = 0; i < 2; i++)
405 for (j = 0; j < 2; j++)
407 t->ktime[i][j].t = 0.0;
408 t->ktime[i][j].c = 0;
412 t->pruneTime.t = 0.0;
413 t->dynamicPruneTime.c = 0;
414 t->dynamicPruneTime.t = 0.0;
417 /*! Initializes simulation constant data. */
418 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t *nb,
419 const interaction_const_t *ic,
420 const NbnxnListParameters *listParams,
421 const nbnxn_atomdata_t *nbat)
423 init_atomdata_first(nb->atdat, nbat->ntype);
424 init_nbparam(nb->nbparam, ic, listParams, nbat, nb->dev_info);
426 /* clear energy and shift force outputs */
427 nbnxn_cuda_clear_e_fshift(nb);
430 void nbnxn_gpu_init(gmx_nbnxn_cuda_t **p_nb,
431 const gmx_device_info_t *deviceInfo,
432 const interaction_const_t *ic,
433 const NbnxnListParameters *listParams,
434 const nbnxn_atomdata_t *nbat,
436 gmx_bool bLocalAndNonlocal)
439 gmx_nbnxn_cuda_t *nb;
448 snew(nb->nbparam, 1);
449 snew(nb->plist[eintLocal], 1);
450 if (bLocalAndNonlocal)
452 snew(nb->plist[eintNonlocal], 1);
455 nb->bUseTwoStreams = bLocalAndNonlocal;
457 nb->timers = new cu_timers_t();
458 snew(nb->timings, 1);
461 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
462 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
463 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
465 init_plist(nb->plist[eintLocal]);
467 /* set device info, just point it to the right GPU among the detected ones */
468 nb->dev_info = deviceInfo;
470 /* local/non-local GPU streams */
471 stat = cudaStreamCreate(&nb->stream[eintLocal]);
472 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
473 if (nb->bUseTwoStreams)
475 init_plist(nb->plist[eintNonlocal]);
477 /* Note that the device we're running on does not have to support
478 * priorities, because we are querying the priority range which in this
479 * case will be a single value.
481 int highest_priority;
482 stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
483 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
485 stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
488 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
491 /* init events for sychronization (timing disabled for performance reasons!) */
492 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
493 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
494 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
495 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
497 /* WARNING: CUDA timings are incorrect with multiple streams.
498 * This is the main reason why they are disabled by default.
500 // TODO: Consider turning on by default when we can detect nr of streams.
501 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
505 init_timers(nb->timers, nb->bUseTwoStreams);
506 init_timings(nb->timings);
509 /* set the kernel type for the current GPU */
510 /* pick L1 cache configuration */
511 nbnxn_cuda_set_cacheconfig(nb->dev_info);
513 nbnxn_cuda_init_const(nb, ic, listParams, nbat);
519 fprintf(debug, "Initialized CUDA data structures.\n");
523 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t *nb,
524 const nbnxn_pairlist_t *h_plist,
528 bool bDoTime = (nb->bDoTime && h_plist->nsci > 0);
529 cudaStream_t stream = nb->stream[iloc];
530 cu_plist_t *d_plist = nb->plist[iloc];
532 if (d_plist->na_c < 0)
534 d_plist->na_c = h_plist->na_ci;
538 if (d_plist->na_c != h_plist->na_ci)
540 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
541 d_plist->na_c, h_plist->na_ci);
548 nb->timers->pl_h2d[iloc].openTimingRegion(stream);
549 nb->timers->didPairlistH2D[iloc] = true;
552 Context context = nullptr;
554 reallocateDeviceBuffer(&d_plist->sci, h_plist->nsci,
555 &d_plist->nsci, &d_plist->sci_nalloc, context);
556 copyToDeviceBuffer(&d_plist->sci, h_plist->sci, 0, h_plist->nsci,
557 stream, GpuApiCallBehavior::Async,
558 bDoTime ? nb->timers->pl_h2d[iloc].fetchNextEvent() : nullptr);
560 reallocateDeviceBuffer(&d_plist->cj4, h_plist->ncj4,
561 &d_plist->ncj4, &d_plist->cj4_nalloc, context);
562 copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4, 0, h_plist->ncj4,
563 stream, GpuApiCallBehavior::Async,
564 bDoTime ? nb->timers->pl_h2d[iloc].fetchNextEvent() : nullptr);
566 reallocateDeviceBuffer(&d_plist->imask, h_plist->ncj4*c_nbnxnGpuClusterpairSplit,
567 &d_plist->nimask, &d_plist->imask_nalloc, context);
569 reallocateDeviceBuffer(&d_plist->excl, h_plist->nexcl,
570 &d_plist->nexcl, &d_plist->excl_nalloc, context);
571 copyToDeviceBuffer(&d_plist->excl, h_plist->excl, 0, h_plist->nexcl,
572 stream, GpuApiCallBehavior::Async,
573 bDoTime ? nb->timers->pl_h2d[iloc].fetchNextEvent() : nullptr);
577 nb->timers->pl_h2d[iloc].closeTimingRegion(stream);
580 /* the next use of thist list we be the first one, so we need to prune */
581 d_plist->haveFreshList = true;
584 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t *nb,
585 const nbnxn_atomdata_t *nbatom)
587 cu_atomdata_t *adat = nb->atdat;
588 cudaStream_t ls = nb->stream[eintLocal];
590 /* only if we have a dynamic box */
591 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
593 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
594 SHIFTS * sizeof(*adat->shift_vec), ls);
595 adat->bShiftVecUploaded = true;
599 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
600 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
603 cu_atomdata_t *adat = nb->atdat;
604 cudaStream_t ls = nb->stream[eintLocal];
606 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
607 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
610 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
611 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
614 cu_atomdata_t *adat = nb->atdat;
615 cudaStream_t ls = nb->stream[eintLocal];
617 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
618 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
619 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
620 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
621 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
622 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
625 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
627 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
628 /* clear shift force array and energies if the outputs were
629 used in the current step */
630 if (flags & GMX_FORCE_VIRIAL)
632 nbnxn_cuda_clear_e_fshift(nb);
636 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t *nb,
637 const nbnxn_atomdata_t *nbat)
642 bool bDoTime = nb->bDoTime;
643 cu_timers_t *timers = nb->timers;
644 cu_atomdata_t *d_atdat = nb->atdat;
645 cudaStream_t ls = nb->stream[eintLocal];
647 natoms = nbat->natoms;
652 /* time async copy */
653 timers->atdat.openTimingRegion(ls);
656 /* need to reallocate if we have to copy more atoms than the amount of space
657 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
658 if (natoms > d_atdat->nalloc)
660 nalloc = over_alloc_small(natoms);
662 /* free up first if the arrays have already been initialized */
663 if (d_atdat->nalloc != -1)
665 freeDeviceBuffer(&d_atdat->f);
666 freeDeviceBuffer(&d_atdat->xq);
667 freeDeviceBuffer(&d_atdat->atom_types);
668 freeDeviceBuffer(&d_atdat->lj_comb);
671 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
672 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
673 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
674 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
675 if (useLjCombRule(nb->nbparam))
677 stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
678 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
682 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
683 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
686 d_atdat->nalloc = nalloc;
690 d_atdat->natoms = natoms;
691 d_atdat->natoms_local = nbat->natoms_local;
693 /* need to clear GPU f output if realloc happened */
696 nbnxn_cuda_clear_f(nb, nalloc);
699 if (useLjCombRule(nb->nbparam))
701 cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
702 natoms*sizeof(*d_atdat->lj_comb), ls);
706 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
707 natoms*sizeof(*d_atdat->atom_types), ls);
712 timers->atdat.closeTimingRegion(ls);
716 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam,
717 const gmx_device_info_t *dev_info)
719 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
721 destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj,
726 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
729 cu_atomdata_t *atdat;
730 cu_nbparam_t *nbparam;
738 nbparam = nb->nbparam;
740 nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
742 stat = cudaEventDestroy(nb->nonlocal_done);
743 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
744 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
745 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
750 /* The non-local counters/stream (second in the array) are needed only with DD. */
751 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
753 stat = cudaStreamDestroy(nb->stream[i]);
754 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
758 if (!useLjCombRule(nb->nbparam))
760 destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj,
765 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
767 destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj,
771 stat = cudaFree(atdat->shift_vec);
772 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
773 stat = cudaFree(atdat->fshift);
774 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
776 stat = cudaFree(atdat->e_lj);
777 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
778 stat = cudaFree(atdat->e_el);
779 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
781 freeDeviceBuffer(&atdat->f);
782 freeDeviceBuffer(&atdat->xq);
783 freeDeviceBuffer(&atdat->atom_types);
784 freeDeviceBuffer(&atdat->lj_comb);
787 auto *plist = nb->plist[eintLocal];
788 freeDeviceBuffer(&plist->sci);
789 freeDeviceBuffer(&plist->cj4);
790 freeDeviceBuffer(&plist->imask);
791 freeDeviceBuffer(&plist->excl);
793 if (nb->bUseTwoStreams)
795 auto *plist_nl = nb->plist[eintNonlocal];
796 freeDeviceBuffer(&plist_nl->sci);
797 freeDeviceBuffer(&plist_nl->cj4);
798 freeDeviceBuffer(&plist_nl->imask);
799 freeDeviceBuffer(&plist_nl->excl);
804 pfree(nb->nbst.e_lj);
805 nb->nbst.e_lj = nullptr;
807 pfree(nb->nbst.e_el);
808 nb->nbst.e_el = nullptr;
810 pfree(nb->nbst.fshift);
811 nb->nbst.fshift = nullptr;
820 fprintf(debug, "Cleaned up CUDA data structures.\n");
824 //! This function is documented in the header file
825 gmx_wallclock_gpu_nbnxn_t *nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
827 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
830 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
832 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
834 init_timings(nbv->gpu_nbv->timings);
838 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
840 return nb != nullptr ?
841 gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
845 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
847 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
848 (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));