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_common_utils.h"
56 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
57 #include "gromacs/mdtypes/interaction_const.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/pbcutil/ishift.h"
60 #include "gromacs/timing/gpu_timing.h"
61 #include "gromacs/utility/basedefinitions.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/real.h"
65 #include "gromacs/utility/smalloc.h"
67 #include "nbnxn_cuda.h"
68 #include "nbnxn_cuda_types.h"
70 /* This is a heuristically determined parameter for the Fermi, Kepler
71 * and Maxwell architectures for the minimum size of ci lists by multiplying
72 * this constant with the # of multiprocessors on the current device.
73 * Since the maximum number of blocks per multiprocessor is 16, the ideal
74 * count for small systems is 32 or 48 blocks per multiprocessor. Because
75 * there is a bit of fluctuations in the generated block counts, we use
76 * a target of 44 instead of the ideal value of 48.
78 static unsigned int gpu_min_ci_balanced_factor = 44;
81 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
84 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam,
85 const gmx_device_info_t *dev_info);
87 /*! \brief Return whether combination rules are used.
89 * \param[in] pointer to nonbonded paramter struct
90 * \return true if combination rules are used in this run, false otherwise
92 static inline bool useLjCombRule(const cu_nbparam_t *nbparam)
94 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
95 nbparam->vdwtype == evdwCuCUTCOMBLB);
98 /*! \brief Initialized the Ewald Coulomb correction GPU table.
100 Tabulates the Ewald Coulomb force and initializes the size/scale
101 and the table GPU array. If called with an already allocated table,
102 it just re-uploads the table.
104 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
106 const gmx_device_info_t *dev_info)
108 if (nbp->coulomb_tab != NULL)
110 nbnxn_cuda_free_nbparam_table(nbp, dev_info);
113 nbp->coulomb_tab_scale = ic->tabq_scale;
114 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
115 ic->tabq_coul_F, ic->tabq_size, dev_info);
119 /*! Initializes the atomdata structure first time, it only gets filled at
121 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
126 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
127 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
128 ad->bShiftVecUploaded = false;
130 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
131 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
133 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
134 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
135 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
136 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
138 /* initialize to NULL poiters to data that is not allocated here and will
139 need reallocation in nbnxn_cuda_init_atomdata */
143 /* size -1 indicates that the respective array hasn't been initialized yet */
148 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
149 earlier GPUs, single or twin cut-off. */
150 static int pick_ewald_kernel_type(bool bTwinCut,
151 const gmx_device_info_t *dev_info)
153 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
156 /* Benchmarking/development environment variables to force the use of
157 analytical or tabulated Ewald kernel. */
158 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
159 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
161 if (bForceAnalyticalEwald && bForceTabulatedEwald)
163 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
164 "requested through environment variables.");
167 /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
168 if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
170 bUseAnalyticalEwald = true;
174 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
179 bUseAnalyticalEwald = false;
183 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
187 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
188 forces it (use it for debugging/benchmarking only). */
189 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
191 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
195 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
201 /*! Copies all parameters related to the cut-off from ic to nbp */
202 static void set_cutoff_parameters(cu_nbparam_t *nbp,
203 const interaction_const_t *ic,
204 const NbnxnListParameters *listParams)
206 nbp->ewald_beta = ic->ewaldcoeff_q;
207 nbp->sh_ewald = ic->sh_ewald;
208 nbp->epsfac = ic->epsfac;
209 nbp->two_k_rf = 2.0 * ic->k_rf;
210 nbp->c_rf = ic->c_rf;
211 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
212 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
213 nbp->rlistOuter_sq = listParams->rlistOuter * listParams->rlistOuter;
214 nbp->rlistInner_sq = listParams->rlistInner * listParams->rlistInner;
215 nbp->useDynamicPruning = listParams->useDynamicPruning;
217 nbp->sh_lj_ewald = ic->sh_lj_ewald;
218 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
220 nbp->rvdw_switch = ic->rvdw_switch;
221 nbp->dispersion_shift = ic->dispersion_shift;
222 nbp->repulsion_shift = ic->repulsion_shift;
223 nbp->vdw_switch = ic->vdw_switch;
226 /*! Initializes the nonbonded parameter data structure. */
227 static void init_nbparam(cu_nbparam_t *nbp,
228 const interaction_const_t *ic,
229 const NbnxnListParameters *listParams,
230 const nbnxn_atomdata_t *nbat,
231 const gmx_device_info_t *dev_info)
235 ntypes = nbat->ntype;
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 (nbat->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!");
270 case eintmodFORCESWITCH:
271 nbp->vdwtype = evdwCuFSWITCH;
273 case eintmodPOTSWITCH:
274 nbp->vdwtype = evdwCuPSWITCH;
277 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
281 else if (ic->vdwtype == evdwPME)
283 if (ic->ljpme_comb_rule == ljcrGEOM)
285 assert(nbat->comb_rule == ljcrGEOM);
286 nbp->vdwtype = evdwCuEWALDGEOM;
290 assert(nbat->comb_rule == ljcrLB);
291 nbp->vdwtype = evdwCuEWALDLB;
296 gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
299 if (ic->eeltype == eelCUT)
301 nbp->eeltype = eelCuCUT;
303 else if (EEL_RF(ic->eeltype))
305 nbp->eeltype = eelCuRF;
307 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
309 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
310 nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
314 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
315 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
318 /* generate table for PME */
319 nbp->coulomb_tab = NULL;
320 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
322 init_ewald_coulomb_force_table(ic, nbp, dev_info);
325 /* set up LJ parameter lookup table */
326 if (!useLjCombRule(nbp))
328 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
329 nbat->nbfp, 2*ntypes*ntypes, dev_info);
332 /* set up LJ-PME parameter lookup table */
333 if (ic->vdwtype == evdwPME)
335 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
336 nbat->nbfp_comb, 2*ntypes, dev_info);
340 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
341 * electrostatic type switching to twin cut-off (or back) if needed. */
342 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t *nbv,
343 const interaction_const_t *ic,
344 const NbnxnListParameters *listParams)
346 if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
350 gmx_nbnxn_cuda_t *nb = nbv->gpu_nbv;
351 cu_nbparam_t *nbp = nb->nbparam;
353 set_cutoff_parameters(nbp, ic, listParams);
355 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
358 init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
361 /*! Initializes the pair list data structure. */
362 static void init_plist(cu_plist_t *pl)
364 /* initialize to NULL pointers to data that is not allocated here and will
365 need reallocation in nbnxn_gpu_init_pairlist */
371 /* size -1 indicates that the respective array hasn't been initialized yet */
378 pl->imask_nalloc = -1;
380 pl->excl_nalloc = -1;
381 pl->haveFreshList = false;
384 /*! Initializes the timer data structure. */
385 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
387 /* The non-local counters/stream (second in the array) are needed only with DD. */
388 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
390 t->didPairlistH2D[i] = false;
391 t->didPrune[i] = false;
392 t->didRollingPrune[i] = false;
396 /*! Initializes the timings data structure. */
397 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
406 for (i = 0; i < 2; i++)
408 for (j = 0; j < 2; j++)
410 t->ktime[i][j].t = 0.0;
411 t->ktime[i][j].c = 0;
415 t->pruneTime.t = 0.0;
416 t->dynamicPruneTime.c = 0;
417 t->dynamicPruneTime.t = 0.0;
420 /*! Initializes simulation constant data. */
421 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t *nb,
422 const interaction_const_t *ic,
423 const NbnxnListParameters *listParams,
424 const nbnxn_atomdata_t *nbat)
426 init_atomdata_first(nb->atdat, nbat->ntype);
427 init_nbparam(nb->nbparam, ic, listParams, nbat, nb->dev_info);
429 /* clear energy and shift force outputs */
430 nbnxn_cuda_clear_e_fshift(nb);
433 void nbnxn_gpu_init(gmx_nbnxn_cuda_t **p_nb,
434 const gmx_device_info_t *deviceInfo,
435 const interaction_const_t *ic,
436 const NbnxnListParameters *listParams,
437 const nbnxn_atomdata_t *nbat,
439 gmx_bool bLocalAndNonlocal)
442 gmx_nbnxn_cuda_t *nb;
451 snew(nb->nbparam, 1);
452 snew(nb->plist[eintLocal], 1);
453 if (bLocalAndNonlocal)
455 snew(nb->plist[eintNonlocal], 1);
458 nb->bUseTwoStreams = bLocalAndNonlocal;
460 nb->timers = new cu_timers_t();
461 snew(nb->timings, 1);
464 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
465 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
466 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
468 init_plist(nb->plist[eintLocal]);
470 /* set device info, just point it to the right GPU among the detected ones */
471 nb->dev_info = deviceInfo;
473 /* local/non-local GPU streams */
474 stat = cudaStreamCreate(&nb->stream[eintLocal]);
475 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
476 if (nb->bUseTwoStreams)
478 init_plist(nb->plist[eintNonlocal]);
480 /* Note that the device we're running on does not have to support
481 * priorities, because we are querying the priority range which in this
482 * case will be a single value.
484 int highest_priority;
485 stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
486 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
488 stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
491 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
494 /* init events for sychronization (timing disabled for performance reasons!) */
495 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
496 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
497 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
498 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
500 /* WARNING: CUDA timings are incorrect with multiple streams.
501 * This is the main reason why they are disabled by default.
503 // TODO: Consider turning on by default when we can detect nr of streams.
504 nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != NULL);
508 init_timers(nb->timers, nb->bUseTwoStreams);
509 init_timings(nb->timings);
512 /* set the kernel type for the current GPU */
513 /* pick L1 cache configuration */
514 nbnxn_cuda_set_cacheconfig(nb->dev_info);
516 nbnxn_cuda_init_const(nb, ic, listParams, nbat);
522 fprintf(debug, "Initialized CUDA data structures.\n");
526 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t *nb,
527 const nbnxn_pairlist_t *h_plist,
530 if (canSkipWork(nb, iloc))
536 bool bDoTime = nb->bDoTime;
537 cudaStream_t stream = nb->stream[iloc];
538 cu_plist_t *d_plist = nb->plist[iloc];
540 if (d_plist->na_c < 0)
542 d_plist->na_c = h_plist->na_ci;
546 if (d_plist->na_c != h_plist->na_ci)
548 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
549 d_plist->na_c, h_plist->na_ci);
556 nb->timers->pl_h2d[iloc].openTimingRegion(stream);
557 nb->timers->didPairlistH2D[iloc] = true;
560 cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
561 &d_plist->nsci, &d_plist->sci_nalloc,
565 cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
566 &d_plist->ncj4, &d_plist->cj4_nalloc,
570 /* this call only allocates space on the device (no data is transferred) */
571 cu_realloc_buffered((void **)&d_plist->imask, NULL, sizeof(*d_plist->imask),
572 &d_plist->nimask, &d_plist->imask_nalloc,
573 h_plist->ncj4*c_nbnxnGpuClusterpairSplit,
576 cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
577 &d_plist->nexcl, &d_plist->excl_nalloc,
583 nb->timers->pl_h2d[iloc].closeTimingRegion(stream);
586 /* the next use of thist list we be the first one, so we need to prune */
587 d_plist->haveFreshList = true;
590 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t *nb,
591 const nbnxn_atomdata_t *nbatom)
593 cu_atomdata_t *adat = nb->atdat;
594 cudaStream_t ls = nb->stream[eintLocal];
596 /* only if we have a dynamic box */
597 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
599 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
600 SHIFTS * sizeof(*adat->shift_vec), ls);
601 adat->bShiftVecUploaded = true;
605 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
606 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
609 cu_atomdata_t *adat = nb->atdat;
610 cudaStream_t ls = nb->stream[eintLocal];
612 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
613 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
616 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
617 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
620 cu_atomdata_t *adat = nb->atdat;
621 cudaStream_t ls = nb->stream[eintLocal];
623 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
624 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
625 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
626 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
627 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
628 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
631 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
633 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
634 /* clear shift force array and energies if the outputs were
635 used in the current step */
636 if (flags & GMX_FORCE_VIRIAL)
638 nbnxn_cuda_clear_e_fshift(nb);
642 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t *nb,
643 const nbnxn_atomdata_t *nbat)
648 bool bDoTime = nb->bDoTime;
649 cu_timers_t *timers = nb->timers;
650 cu_atomdata_t *d_atdat = nb->atdat;
651 cudaStream_t ls = nb->stream[eintLocal];
653 natoms = nbat->natoms;
658 /* time async copy */
659 timers->atdat.openTimingRegion(ls);
662 /* need to reallocate if we have to copy more atoms than the amount of space
663 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
664 if (natoms > d_atdat->nalloc)
666 nalloc = over_alloc_small(natoms);
668 /* free up first if the arrays have already been initialized */
669 if (d_atdat->nalloc != -1)
671 cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
672 cu_free_buffered(d_atdat->xq);
673 cu_free_buffered(d_atdat->atom_types);
674 cu_free_buffered(d_atdat->lj_comb);
677 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
678 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
679 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
680 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
681 if (useLjCombRule(nb->nbparam))
683 stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
684 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
688 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
689 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
692 d_atdat->nalloc = nalloc;
696 d_atdat->natoms = natoms;
697 d_atdat->natoms_local = nbat->natoms_local;
699 /* need to clear GPU f output if realloc happened */
702 nbnxn_cuda_clear_f(nb, nalloc);
705 if (useLjCombRule(nb->nbparam))
707 cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
708 natoms*sizeof(*d_atdat->lj_comb), ls);
712 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
713 natoms*sizeof(*d_atdat->atom_types), ls);
718 timers->atdat.closeTimingRegion(ls);
722 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam,
723 const gmx_device_info_t *dev_info)
725 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
727 destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj,
732 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
735 cu_atomdata_t *atdat;
736 cu_nbparam_t *nbparam;
737 cu_plist_t *plist, *plist_nl;
745 nbparam = nb->nbparam;
746 plist = nb->plist[eintLocal];
747 plist_nl = nb->plist[eintNonlocal];
749 nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
751 stat = cudaEventDestroy(nb->nonlocal_done);
752 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
753 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
754 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
759 /* The non-local counters/stream (second in the array) are needed only with DD. */
760 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
762 stat = cudaStreamDestroy(nb->stream[i]);
763 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
767 if (!useLjCombRule(nb->nbparam))
769 destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj,
774 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
776 destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj,
780 stat = cudaFree(atdat->shift_vec);
781 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
782 stat = cudaFree(atdat->fshift);
783 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
785 stat = cudaFree(atdat->e_lj);
786 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
787 stat = cudaFree(atdat->e_el);
788 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
790 cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
791 cu_free_buffered(atdat->xq);
792 cu_free_buffered(atdat->atom_types, &atdat->ntypes);
793 cu_free_buffered(atdat->lj_comb);
795 cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
796 cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
797 cu_free_buffered(plist->imask, &plist->nimask, &plist->imask_nalloc);
798 cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
799 if (nb->bUseTwoStreams)
801 cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
802 cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
803 cu_free_buffered(plist_nl->imask, &plist_nl->nimask, &plist_nl->imask_nalloc);
804 cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
810 if (nb->bUseTwoStreams)
819 fprintf(debug, "Cleaned up CUDA data structures.\n");
823 //! This function is documented in the header file
824 gmx_wallclock_gpu_nbnxn_t *nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
826 return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
829 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
831 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
833 init_timings(nbv->gpu_nbv->timings);
837 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
840 gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
844 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
846 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
847 (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));