2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017, 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.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 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
72 /* This is a heuristically determined parameter for the Fermi, Kepler
73 * and Maxwell architectures for the minimum size of ci lists by multiplying
74 * this constant with the # of multiprocessors on the current device.
75 * Since the maximum number of blocks per multiprocessor is 16, the ideal
76 * count for small systems is 32 or 48 blocks per multiprocessor. Because
77 * there is a bit of fluctuations in the generated block counts, we use
78 * a target of 44 instead of the ideal value of 48.
80 static unsigned int gpu_min_ci_balanced_factor = 44;
83 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
86 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam,
87 const gmx_device_info_t *dev_info);
89 /*! \brief Return whether combination rules are used.
91 * \param[in] pointer to nonbonded paramter struct
92 * \return true if combination rules are used in this run, false otherwise
94 static inline bool useLjCombRule(const cu_nbparam_t *nbparam)
96 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
97 nbparam->vdwtype == evdwCuCUTCOMBLB);
100 /*! \brief Initialized the Ewald Coulomb correction GPU table.
102 Tabulates the Ewald Coulomb force and initializes the size/scale
103 and the table GPU array. If called with an already allocated table,
104 it just re-uploads the table.
106 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
108 const gmx_device_info_t *dev_info)
110 if (nbp->coulomb_tab != NULL)
112 nbnxn_cuda_free_nbparam_table(nbp, dev_info);
115 nbp->coulomb_tab_scale = ic->tabq_scale;
116 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
117 &nbnxn_cuda_get_coulomb_tab_texref(),
118 ic->tabq_coul_F, ic->tabq_size, dev_info);
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 NULL 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,
154 const gmx_device_info_t *dev_info)
156 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
159 /* Benchmarking/development environment variables to force the use of
160 analytical or tabulated Ewald kernel. */
161 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
162 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
164 if (bForceAnalyticalEwald && bForceTabulatedEwald)
166 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
167 "requested through environment variables.");
170 /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
171 if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
173 bUseAnalyticalEwald = true;
177 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
182 bUseAnalyticalEwald = false;
186 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
190 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
191 forces it (use it for debugging/benchmarking only). */
192 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
194 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
198 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
204 /*! Copies all parameters related to the cut-off from ic to nbp */
205 static void set_cutoff_parameters(cu_nbparam_t *nbp,
206 const interaction_const_t *ic,
207 const NbnxnListParameters *listParams)
209 nbp->ewald_beta = ic->ewaldcoeff_q;
210 nbp->sh_ewald = ic->sh_ewald;
211 nbp->epsfac = ic->epsfac;
212 nbp->two_k_rf = 2.0 * ic->k_rf;
213 nbp->c_rf = ic->c_rf;
214 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
215 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
216 nbp->rlistOuter_sq = listParams->rlistOuter * listParams->rlistOuter;
217 nbp->rlistInner_sq = listParams->rlistInner * listParams->rlistInner;
218 nbp->useDynamicPruning = listParams->useDynamicPruning;
220 nbp->sh_lj_ewald = ic->sh_lj_ewald;
221 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
223 nbp->rvdw_switch = ic->rvdw_switch;
224 nbp->dispersion_shift = ic->dispersion_shift;
225 nbp->repulsion_shift = ic->repulsion_shift;
226 nbp->vdw_switch = ic->vdw_switch;
229 /*! Initializes the nonbonded parameter data structure. */
230 static void init_nbparam(cu_nbparam_t *nbp,
231 const interaction_const_t *ic,
232 const NbnxnListParameters *listParams,
233 const nbnxn_atomdata_t *nbat,
234 const gmx_device_info_t *dev_info)
238 ntypes = nbat->ntype;
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 (nbat->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!");
273 case eintmodFORCESWITCH:
274 nbp->vdwtype = evdwCuFSWITCH;
276 case eintmodPOTSWITCH:
277 nbp->vdwtype = evdwCuPSWITCH;
280 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
284 else if (ic->vdwtype == evdwPME)
286 if (ic->ljpme_comb_rule == ljcrGEOM)
288 assert(nbat->comb_rule == ljcrGEOM);
289 nbp->vdwtype = evdwCuEWALDGEOM;
293 assert(nbat->comb_rule == ljcrLB);
294 nbp->vdwtype = evdwCuEWALDLB;
299 gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
302 if (ic->eeltype == eelCUT)
304 nbp->eeltype = eelCuCUT;
306 else if (EEL_RF(ic->eeltype))
308 nbp->eeltype = eelCuRF;
310 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
312 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
313 nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
317 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
318 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
321 /* generate table for PME */
322 nbp->coulomb_tab = NULL;
323 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
325 init_ewald_coulomb_force_table(ic, nbp, dev_info);
328 /* set up LJ parameter lookup table */
329 if (!useLjCombRule(nbp))
331 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
332 &nbnxn_cuda_get_nbfp_texref(),
333 nbat->nbfp, 2*ntypes*ntypes, dev_info);
336 /* set up LJ-PME parameter lookup table */
337 if (ic->vdwtype == evdwPME)
339 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
340 &nbnxn_cuda_get_nbfp_comb_texref(),
341 nbat->nbfp_comb, 2*ntypes, dev_info);
345 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
346 * electrostatic type switching to twin cut-off (or back) if needed. */
347 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t *nbv,
348 const interaction_const_t *ic,
349 const NbnxnListParameters *listParams)
351 if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
355 gmx_nbnxn_cuda_t *nb = nbv->gpu_nbv;
356 cu_nbparam_t *nbp = nb->nbparam;
358 set_cutoff_parameters(nbp, ic, listParams);
360 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
363 init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
366 /*! Initializes the pair list data structure. */
367 static void init_plist(cu_plist_t *pl)
369 /* initialize to NULL pointers to data that is not allocated here and will
370 need reallocation in nbnxn_gpu_init_pairlist */
376 /* size -1 indicates that the respective array hasn't been initialized yet */
383 pl->imask_nalloc = -1;
385 pl->excl_nalloc = -1;
386 pl->haveFreshList = false;
389 /*! Initializes the timer data structure. */
390 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
393 int eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
395 stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
396 CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
397 stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
398 CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
400 /* The non-local counters/stream (second in the array) are needed only with DD. */
401 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
403 stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
404 CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
405 stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
406 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
408 stat = cudaEventCreateWithFlags(&(t->start_prune_k[i]), eventflags);
409 CU_RET_ERR(stat, "cudaEventCreate on start_prune_k failed");
410 stat = cudaEventCreateWithFlags(&(t->stop_prune_k[i]), eventflags);
411 CU_RET_ERR(stat, "cudaEventCreate on stop_prune_k failed");
413 stat = cudaEventCreateWithFlags(&(t->start_rollingPrune_k[i]), eventflags);
414 CU_RET_ERR(stat, "cudaEventCreate on start_rollingPrune_k failed");
415 stat = cudaEventCreateWithFlags(&(t->stop_rollingPrune_k[i]), eventflags);
416 CU_RET_ERR(stat, "cudaEventCreate on stop_rollingPrune_k failed");
418 stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
419 CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
420 stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
421 CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
423 stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
424 CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
425 stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
426 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
428 stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
429 CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
430 stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
431 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
433 t->didPairlistH2D[i] = false;
434 t->didPrune[i] = false;
435 t->didRollingPrune[i] = false;
439 /*! Initializes the timings data structure. */
440 static void init_timings(gmx_wallclock_gpu_t *t)
449 for (i = 0; i < 2; i++)
451 for (j = 0; j < 2; j++)
453 t->ktime[i][j].t = 0.0;
454 t->ktime[i][j].c = 0;
458 t->pruneTime.t = 0.0;
459 t->dynamicPruneTime.c = 0;
460 t->dynamicPruneTime.t = 0.0;
463 /*! Initializes simulation constant data. */
464 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t *nb,
465 const interaction_const_t *ic,
466 const NbnxnListParameters *listParams,
467 const nonbonded_verlet_group_t *nbv_group)
469 init_atomdata_first(nb->atdat, nbv_group[0].nbat->ntype);
470 init_nbparam(nb->nbparam, ic, listParams, nbv_group[0].nbat, nb->dev_info);
472 /* clear energy and shift force outputs */
473 nbnxn_cuda_clear_e_fshift(nb);
476 void nbnxn_gpu_init(gmx_nbnxn_cuda_t **p_nb,
477 const gmx_device_info_t *deviceInfo,
478 const interaction_const_t *ic,
479 const NbnxnListParameters *listParams,
480 nonbonded_verlet_group_t *nbv_grp,
482 gmx_bool bLocalAndNonlocal)
485 gmx_nbnxn_cuda_t *nb;
494 snew(nb->nbparam, 1);
495 snew(nb->plist[eintLocal], 1);
496 if (bLocalAndNonlocal)
498 snew(nb->plist[eintNonlocal], 1);
501 nb->bUseTwoStreams = bLocalAndNonlocal;
504 snew(nb->timings, 1);
507 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
508 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
509 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
511 init_plist(nb->plist[eintLocal]);
513 /* set device info, just point it to the right GPU among the detected ones */
514 nb->dev_info = deviceInfo;
516 /* local/non-local GPU streams */
517 stat = cudaStreamCreate(&nb->stream[eintLocal]);
518 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
519 if (nb->bUseTwoStreams)
521 init_plist(nb->plist[eintNonlocal]);
523 /* Note that the device we're running on does not have to support
524 * priorities, because we are querying the priority range which in this
525 * case will be a single value.
527 int highest_priority;
528 stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
529 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
531 stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
534 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
537 /* init events for sychronization (timing disabled for performance reasons!) */
538 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
539 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
540 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
541 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
543 /* CUDA timing disabled as event timers don't work:
544 - with multiple streams = domain-decomposition;
545 - when turned off by GMX_DISABLE_CUDA_TIMING/GMX_DISABLE_GPU_TIMING.
547 nb->bDoTime = (!nb->bUseTwoStreams &&
548 (getenv("GMX_DISABLE_CUDA_TIMING") == NULL) &&
549 (getenv("GMX_DISABLE_GPU_TIMING") == NULL));
553 init_timers(nb->timers, nb->bUseTwoStreams);
554 init_timings(nb->timings);
557 /* set the kernel type for the current GPU */
558 /* pick L1 cache configuration */
559 nbnxn_cuda_set_cacheconfig(nb->dev_info);
561 nbnxn_cuda_init_const(nb, ic, listParams, nbv_grp);
567 fprintf(debug, "Initialized CUDA data structures.\n");
571 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t *nb,
572 const nbnxn_pairlist_t *h_plist,
575 if (canSkipWork(nb, iloc))
582 bool bDoTime = nb->bDoTime;
583 cudaStream_t stream = nb->stream[iloc];
584 cu_plist_t *d_plist = nb->plist[iloc];
586 if (d_plist->na_c < 0)
588 d_plist->na_c = h_plist->na_ci;
592 if (d_plist->na_c != h_plist->na_ci)
594 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
595 d_plist->na_c, h_plist->na_ci);
602 stat = cudaEventRecord(nb->timers->start_pl_h2d[iloc], stream);
603 CU_RET_ERR(stat, "cudaEventRecord failed");
604 nb->timers->didPairlistH2D[iloc] = true;
607 cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
608 &d_plist->nsci, &d_plist->sci_nalloc,
612 cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
613 &d_plist->ncj4, &d_plist->cj4_nalloc,
617 /* this call only allocates space on the device (no data is transferred) */
618 cu_realloc_buffered((void **)&d_plist->imask, NULL, sizeof(*d_plist->imask),
619 &d_plist->nimask, &d_plist->imask_nalloc,
620 h_plist->ncj4*c_nbnxnGpuClusterpairSplit,
623 cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
624 &d_plist->nexcl, &d_plist->excl_nalloc,
630 stat = cudaEventRecord(nb->timers->stop_pl_h2d[iloc], stream);
631 CU_RET_ERR(stat, "cudaEventRecord failed");
634 /* the next use of thist list we be the first one, so we need to prune */
635 d_plist->haveFreshList = true;
638 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t *nb,
639 const nbnxn_atomdata_t *nbatom)
641 cu_atomdata_t *adat = nb->atdat;
642 cudaStream_t ls = nb->stream[eintLocal];
644 /* only if we have a dynamic box */
645 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
647 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
648 SHIFTS * sizeof(*adat->shift_vec), ls);
649 adat->bShiftVecUploaded = true;
653 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
654 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
657 cu_atomdata_t *adat = nb->atdat;
658 cudaStream_t ls = nb->stream[eintLocal];
660 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
661 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
664 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
665 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
668 cu_atomdata_t *adat = nb->atdat;
669 cudaStream_t ls = nb->stream[eintLocal];
671 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
672 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
673 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
674 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
675 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
676 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
679 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
681 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
682 /* clear shift force array and energies if the outputs were
683 used in the current step */
684 if (flags & GMX_FORCE_VIRIAL)
686 nbnxn_cuda_clear_e_fshift(nb);
690 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t *nb,
691 const struct nbnxn_atomdata_t *nbat)
696 bool bDoTime = nb->bDoTime;
697 cu_timers_t *timers = nb->timers;
698 cu_atomdata_t *d_atdat = nb->atdat;
699 cudaStream_t ls = nb->stream[eintLocal];
701 natoms = nbat->natoms;
706 /* time async copy */
707 stat = cudaEventRecord(timers->start_atdat, ls);
708 CU_RET_ERR(stat, "cudaEventRecord failed");
711 /* need to reallocate if we have to copy more atoms than the amount of space
712 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
713 if (natoms > d_atdat->nalloc)
715 nalloc = over_alloc_small(natoms);
717 /* free up first if the arrays have already been initialized */
718 if (d_atdat->nalloc != -1)
720 cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
721 cu_free_buffered(d_atdat->xq);
722 cu_free_buffered(d_atdat->atom_types);
723 cu_free_buffered(d_atdat->lj_comb);
726 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
727 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
728 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
729 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
730 if (useLjCombRule(nb->nbparam))
732 stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
733 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
737 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
738 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
741 d_atdat->nalloc = nalloc;
745 d_atdat->natoms = natoms;
746 d_atdat->natoms_local = nbat->natoms_local;
748 /* need to clear GPU f output if realloc happened */
751 nbnxn_cuda_clear_f(nb, nalloc);
754 if (useLjCombRule(nb->nbparam))
756 cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
757 natoms*sizeof(*d_atdat->lj_comb), ls);
761 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
762 natoms*sizeof(*d_atdat->atom_types), ls);
767 stat = cudaEventRecord(timers->stop_atdat, ls);
768 CU_RET_ERR(stat, "cudaEventRecord failed");
772 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam,
773 const gmx_device_info_t *dev_info)
775 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
777 destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj,
778 &nbnxn_cuda_get_coulomb_tab_texref(), dev_info);
782 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
785 cu_atomdata_t *atdat;
786 cu_nbparam_t *nbparam;
787 cu_plist_t *plist, *plist_nl;
796 nbparam = nb->nbparam;
797 plist = nb->plist[eintLocal];
798 plist_nl = nb->plist[eintNonlocal];
801 nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
803 stat = cudaEventDestroy(nb->nonlocal_done);
804 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
805 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
806 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
810 stat = cudaEventDestroy(timers->start_atdat);
811 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
812 stat = cudaEventDestroy(timers->stop_atdat);
813 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
815 /* The non-local counters/stream (second in the array) are needed only with DD. */
816 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
818 stat = cudaEventDestroy(timers->start_nb_k[i]);
819 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
820 stat = cudaEventDestroy(timers->stop_nb_k[i]);
821 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
823 stat = cudaEventDestroy(timers->start_prune_k[i]);
824 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_prune_k");
825 stat = cudaEventDestroy(timers->stop_prune_k[i]);
826 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_prune_k");
828 stat = cudaEventDestroy(timers->start_rollingPrune_k[i]);
829 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_rollingPrune_k");
830 stat = cudaEventDestroy(timers->stop_rollingPrune_k[i]);
831 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_rollingPrune_k");
833 stat = cudaEventDestroy(timers->start_pl_h2d[i]);
834 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
835 stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
836 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
838 stat = cudaStreamDestroy(nb->stream[i]);
839 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
841 stat = cudaEventDestroy(timers->start_nb_h2d[i]);
842 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
843 stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
844 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
846 stat = cudaEventDestroy(timers->start_nb_d2h[i]);
847 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
848 stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
849 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
853 if (!useLjCombRule(nb->nbparam))
855 destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj,
856 &nbnxn_cuda_get_nbfp_texref(), nb->dev_info);
860 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
862 destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj,
863 &nbnxn_cuda_get_nbfp_comb_texref(), nb->dev_info);
866 stat = cudaFree(atdat->shift_vec);
867 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
868 stat = cudaFree(atdat->fshift);
869 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
871 stat = cudaFree(atdat->e_lj);
872 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
873 stat = cudaFree(atdat->e_el);
874 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
876 cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
877 cu_free_buffered(atdat->xq);
878 cu_free_buffered(atdat->atom_types, &atdat->ntypes);
879 cu_free_buffered(atdat->lj_comb);
881 cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
882 cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
883 cu_free_buffered(plist->imask, &plist->nimask, &plist->imask_nalloc);
884 cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
885 if (nb->bUseTwoStreams)
887 cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
888 cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
889 cu_free_buffered(plist_nl->imask, &plist_nl->nimask, &plist_nl->imask_nalloc);
890 cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
896 if (nb->bUseTwoStreams)
906 fprintf(debug, "Cleaned up CUDA data structures.\n");
910 gmx_wallclock_gpu_t * nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
912 return (nb != NULL && nb->bDoTime) ? nb->timings : NULL;
915 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
917 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
919 init_timings(nbv->gpu_nbv->timings);
923 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
926 gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
930 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
932 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
933 (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));