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);
90 /*! \brief Return whether texture objects are used on this device.
92 * \param[in] pointer to the GPU device info structure to inspect for texture objects support
93 * \return true if texture objects are used on this device
95 static bool use_texobj(const gmx_device_info_t *dev_info)
97 assert(!c_disableCudaTextures);
98 /* Only device CC >= 3.0 (Kepler and later) support texture objects */
99 return (dev_info->prop.major >= 3);
102 /*! \brief Return whether combination rules are used.
104 * \param[in] pointer to nonbonded paramter struct
105 * \return true if combination rules are used in this run, false otherwise
107 static inline bool useLjCombRule(const cu_nbparam_t *nbparam)
109 return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
110 nbparam->vdwtype == evdwCuCUTCOMBLB);
113 /*! \brief Set up texture object for an array of type T.
115 * Set up texture object for an array of type T and bind it to the device memory
116 * \p d_ptr points to.
118 * \tparam[in] T Raw data type
119 * \param[out] texObj texture object to initialize
120 * \param[in] d_ptr pointer to device global memory to bind \p texObj to
121 * \param[in] sizeInBytes size of memory area to bind \p texObj to
123 template <typename T>
124 static void setup1DTexture(cudaTextureObject_t &texObj,
128 assert(!c_disableCudaTextures);
134 memset(&rd, 0, sizeof(rd));
135 rd.resType = cudaResourceTypeLinear;
136 rd.res.linear.devPtr = d_ptr;
137 rd.res.linear.desc = cudaCreateChannelDesc<T>();
138 rd.res.linear.sizeInBytes = sizeInBytes;
140 memset(&td, 0, sizeof(td));
141 td.readMode = cudaReadModeElementType;
142 stat = cudaCreateTextureObject(&texObj, &rd, &td, NULL);
143 CU_RET_ERR(stat, "cudaCreateTextureObject failed");
146 /*! \brief Set up texture reference for an array of type T.
148 * Set up texture object for an array of type T and bind it to the device memory
149 * \p d_ptr points to.
151 * \tparam[in] T Raw data type
152 * \param[out] texObj texture reference to initialize
153 * \param[in] d_ptr pointer to device global memory to bind \p texObj to
154 * \param[in] sizeInBytes size of memory area to bind \p texObj to
156 template <typename T>
157 static void setup1DTexture(const struct texture<T, 1, cudaReadModeElementType> *texRef,
161 assert(!c_disableCudaTextures);
164 cudaChannelFormatDesc cd;
166 cd = cudaCreateChannelDesc<T>();
167 stat = cudaBindTexture(nullptr, texRef, d_ptr, &cd, sizeInBytes);
168 CU_RET_ERR(stat, "cudaBindTexture failed");
171 /*! \brief Initialize parameter lookup table.
173 * Initializes device memory, copies data from host and binds
174 * a texture to allocated device memory to be used for LJ/Ewald/... parameter
177 * \tparam[in] T Raw data type
178 * \param[out] d_ptr device pointer to the memory to be allocated
179 * \param[out] texObj texture object to be initialized
180 * \param[out] texRef texture reference to be initialized
181 * \param[in] h_ptr pointer to the host memory to be uploaded to the device
182 * \param[in] numElem number of elements in the h_ptr
183 * \param[in] devInfo pointer to the info struct of the device in use
185 template <typename T>
186 static void initParamLookupTable(T * &d_ptr,
187 cudaTextureObject_t &texObj,
188 const struct texture<T, 1, cudaReadModeElementType> *texRef,
191 const gmx_device_info_t *devInfo)
193 const size_t sizeInBytes = numElem * sizeof(*d_ptr);
194 cudaError_t stat = cudaMalloc((void **)&d_ptr, sizeInBytes);
195 CU_RET_ERR(stat, "cudaMalloc failed in initParamLookupTable");
196 cu_copy_H2D(d_ptr, (void *)h_ptr, sizeInBytes);
198 if (!c_disableCudaTextures)
200 if (use_texobj(devInfo))
202 setup1DTexture<T>(texObj, d_ptr, sizeInBytes);
206 setup1DTexture<T>(texRef, d_ptr, sizeInBytes);
211 /*! \brief Initialized the Ewald Coulomb correction GPU table.
213 Tabulates the Ewald Coulomb force and initializes the size/scale
214 and the table GPU array. If called with an already allocated table,
215 it just re-uploads the table.
217 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
219 const gmx_device_info_t *dev_info)
221 if (nbp->coulomb_tab != NULL)
223 nbnxn_cuda_free_nbparam_table(nbp, dev_info);
226 nbp->coulomb_tab_scale = ic->tabq_scale;
227 initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
228 &nbnxn_cuda_get_coulomb_tab_texref(),
229 ic->tabq_coul_F, ic->tabq_size, dev_info);
233 /*! Initializes the atomdata structure first time, it only gets filled at
235 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
240 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
241 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
242 ad->bShiftVecUploaded = false;
244 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
245 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
247 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
248 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
249 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
250 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
252 /* initialize to NULL poiters to data that is not allocated here and will
253 need reallocation in nbnxn_cuda_init_atomdata */
257 /* size -1 indicates that the respective array hasn't been initialized yet */
262 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
263 earlier GPUs, single or twin cut-off. */
264 static int pick_ewald_kernel_type(bool bTwinCut,
265 const gmx_device_info_t *dev_info)
267 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
270 /* Benchmarking/development environment variables to force the use of
271 analytical or tabulated Ewald kernel. */
272 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
273 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
275 if (bForceAnalyticalEwald && bForceTabulatedEwald)
277 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
278 "requested through environment variables.");
281 /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
282 if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
284 bUseAnalyticalEwald = true;
288 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
293 bUseAnalyticalEwald = false;
297 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
301 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
302 forces it (use it for debugging/benchmarking only). */
303 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
305 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
309 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
315 /*! Copies all parameters related to the cut-off from ic to nbp */
316 static void set_cutoff_parameters(cu_nbparam_t *nbp,
317 const interaction_const_t *ic,
318 const NbnxnListParameters *listParams)
320 nbp->ewald_beta = ic->ewaldcoeff_q;
321 nbp->sh_ewald = ic->sh_ewald;
322 nbp->epsfac = ic->epsfac;
323 nbp->two_k_rf = 2.0 * ic->k_rf;
324 nbp->c_rf = ic->c_rf;
325 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
326 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
327 nbp->rlistOuter_sq = listParams->rlistOuter * listParams->rlistOuter;
328 nbp->rlistInner_sq = listParams->rlistInner * listParams->rlistInner;
329 nbp->useDynamicPruning = listParams->useDynamicPruning;
331 nbp->sh_lj_ewald = ic->sh_lj_ewald;
332 nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
334 nbp->rvdw_switch = ic->rvdw_switch;
335 nbp->dispersion_shift = ic->dispersion_shift;
336 nbp->repulsion_shift = ic->repulsion_shift;
337 nbp->vdw_switch = ic->vdw_switch;
340 /*! Initializes the nonbonded parameter data structure. */
341 static void init_nbparam(cu_nbparam_t *nbp,
342 const interaction_const_t *ic,
343 const NbnxnListParameters *listParams,
344 const nbnxn_atomdata_t *nbat,
345 const gmx_device_info_t *dev_info)
349 ntypes = nbat->ntype;
351 set_cutoff_parameters(nbp, ic, listParams);
353 /* The kernel code supports LJ combination rules (geometric and LB) for
354 * all kernel types, but we only generate useful combination rule kernels.
355 * We currently only use LJ combination rule (geometric and LB) kernels
356 * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
357 * with PME and 20% with RF, the other kernels speed up about half as much.
358 * For LJ force-switch the geometric rule would give 7% speed-up, but this
359 * combination is rarely used. LJ force-switch with LB rule is more common,
360 * but gives only 1% speed-up.
362 if (ic->vdwtype == evdwCUT)
364 switch (ic->vdw_modifier)
367 case eintmodPOTSHIFT:
368 switch (nbat->comb_rule)
371 nbp->vdwtype = evdwCuCUT;
374 nbp->vdwtype = evdwCuCUTCOMBGEOM;
377 nbp->vdwtype = evdwCuCUTCOMBLB;
380 gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
384 case eintmodFORCESWITCH:
385 nbp->vdwtype = evdwCuFSWITCH;
387 case eintmodPOTSWITCH:
388 nbp->vdwtype = evdwCuPSWITCH;
391 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
395 else if (ic->vdwtype == evdwPME)
397 if (ic->ljpme_comb_rule == ljcrGEOM)
399 assert(nbat->comb_rule == ljcrGEOM);
400 nbp->vdwtype = evdwCuEWALDGEOM;
404 assert(nbat->comb_rule == ljcrLB);
405 nbp->vdwtype = evdwCuEWALDLB;
410 gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
413 if (ic->eeltype == eelCUT)
415 nbp->eeltype = eelCuCUT;
417 else if (EEL_RF(ic->eeltype))
419 nbp->eeltype = eelCuRF;
421 else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
423 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
424 nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
428 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
429 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
432 /* generate table for PME */
433 nbp->coulomb_tab = NULL;
434 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
436 init_ewald_coulomb_force_table(ic, nbp, dev_info);
439 /* set up LJ parameter lookup table */
440 if (!useLjCombRule(nbp))
442 initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
443 &nbnxn_cuda_get_nbfp_texref(),
444 nbat->nbfp, 2*ntypes*ntypes, dev_info);
447 /* set up LJ-PME parameter lookup table */
448 if (ic->vdwtype == evdwPME)
450 initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
451 &nbnxn_cuda_get_nbfp_comb_texref(),
452 nbat->nbfp_comb, 2*ntypes, dev_info);
456 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
457 * electrostatic type switching to twin cut-off (or back) if needed. */
458 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t *nbv,
459 const interaction_const_t *ic,
460 const NbnxnListParameters *listParams)
462 if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
466 gmx_nbnxn_cuda_t *nb = nbv->gpu_nbv;
467 cu_nbparam_t *nbp = nb->nbparam;
469 set_cutoff_parameters(nbp, ic, listParams);
471 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
474 init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
477 /*! Initializes the pair list data structure. */
478 static void init_plist(cu_plist_t *pl)
480 /* initialize to NULL pointers to data that is not allocated here and will
481 need reallocation in nbnxn_gpu_init_pairlist */
487 /* size -1 indicates that the respective array hasn't been initialized yet */
494 pl->imask_nalloc = -1;
496 pl->excl_nalloc = -1;
497 pl->haveFreshList = false;
500 /*! Initializes the timer data structure. */
501 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
504 int eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
506 stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
507 CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
508 stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
509 CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
511 /* The non-local counters/stream (second in the array) are needed only with DD. */
512 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
514 stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
515 CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
516 stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
517 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
519 stat = cudaEventCreateWithFlags(&(t->start_prune_k[i]), eventflags);
520 CU_RET_ERR(stat, "cudaEventCreate on start_prune_k failed");
521 stat = cudaEventCreateWithFlags(&(t->stop_prune_k[i]), eventflags);
522 CU_RET_ERR(stat, "cudaEventCreate on stop_prune_k failed");
524 stat = cudaEventCreateWithFlags(&(t->start_rollingPrune_k[i]), eventflags);
525 CU_RET_ERR(stat, "cudaEventCreate on start_rollingPrune_k failed");
526 stat = cudaEventCreateWithFlags(&(t->stop_rollingPrune_k[i]), eventflags);
527 CU_RET_ERR(stat, "cudaEventCreate on stop_rollingPrune_k failed");
529 stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
530 CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
531 stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
532 CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
534 stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
535 CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
536 stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
537 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
539 stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
540 CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
541 stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
542 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
544 t->didPairlistH2D[i] = false;
545 t->didPrune[i] = false;
546 t->didRollingPrune[i] = false;
550 /*! Initializes the timings data structure. */
551 static void init_timings(gmx_wallclock_gpu_t *t)
560 for (i = 0; i < 2; i++)
562 for (j = 0; j < 2; j++)
564 t->ktime[i][j].t = 0.0;
565 t->ktime[i][j].c = 0;
569 t->pruneTime.t = 0.0;
570 t->dynamicPruneTime.c = 0;
571 t->dynamicPruneTime.t = 0.0;
574 /*! Initializes simulation constant data. */
575 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t *nb,
576 const interaction_const_t *ic,
577 const NbnxnListParameters *listParams,
578 const nonbonded_verlet_group_t *nbv_group)
580 init_atomdata_first(nb->atdat, nbv_group[0].nbat->ntype);
581 init_nbparam(nb->nbparam, ic, listParams, nbv_group[0].nbat, nb->dev_info);
583 /* clear energy and shift force outputs */
584 nbnxn_cuda_clear_e_fshift(nb);
587 void nbnxn_gpu_init(gmx_nbnxn_cuda_t **p_nb,
588 const gmx_device_info_t *deviceInfo,
589 const interaction_const_t *ic,
590 const NbnxnListParameters *listParams,
591 nonbonded_verlet_group_t *nbv_grp,
593 gmx_bool bLocalAndNonlocal)
596 gmx_nbnxn_cuda_t *nb;
605 snew(nb->nbparam, 1);
606 snew(nb->plist[eintLocal], 1);
607 if (bLocalAndNonlocal)
609 snew(nb->plist[eintNonlocal], 1);
612 nb->bUseTwoStreams = bLocalAndNonlocal;
615 snew(nb->timings, 1);
618 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
619 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
620 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
622 init_plist(nb->plist[eintLocal]);
624 /* set device info, just point it to the right GPU among the detected ones */
625 nb->dev_info = deviceInfo;
627 /* local/non-local GPU streams */
628 stat = cudaStreamCreate(&nb->stream[eintLocal]);
629 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
630 if (nb->bUseTwoStreams)
632 init_plist(nb->plist[eintNonlocal]);
634 /* Note that the device we're running on does not have to support
635 * priorities, because we are querying the priority range which in this
636 * case will be a single value.
638 int highest_priority;
639 stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
640 CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
642 stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
645 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
648 /* init events for sychronization (timing disabled for performance reasons!) */
649 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
650 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
651 stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
652 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
654 /* CUDA timing disabled as event timers don't work:
655 - with multiple streams = domain-decomposition;
656 - when turned off by GMX_DISABLE_CUDA_TIMING/GMX_DISABLE_GPU_TIMING.
658 nb->bDoTime = (!nb->bUseTwoStreams &&
659 (getenv("GMX_DISABLE_CUDA_TIMING") == NULL) &&
660 (getenv("GMX_DISABLE_GPU_TIMING") == NULL));
664 init_timers(nb->timers, nb->bUseTwoStreams);
665 init_timings(nb->timings);
668 /* set the kernel type for the current GPU */
669 /* pick L1 cache configuration */
670 nbnxn_cuda_set_cacheconfig(nb->dev_info);
672 nbnxn_cuda_init_const(nb, ic, listParams, nbv_grp);
678 fprintf(debug, "Initialized CUDA data structures.\n");
682 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t *nb,
683 const nbnxn_pairlist_t *h_plist,
686 if (canSkipWork(nb, iloc))
693 bool bDoTime = nb->bDoTime;
694 cudaStream_t stream = nb->stream[iloc];
695 cu_plist_t *d_plist = nb->plist[iloc];
697 if (d_plist->na_c < 0)
699 d_plist->na_c = h_plist->na_ci;
703 if (d_plist->na_c != h_plist->na_ci)
705 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
706 d_plist->na_c, h_plist->na_ci);
713 stat = cudaEventRecord(nb->timers->start_pl_h2d[iloc], stream);
714 CU_RET_ERR(stat, "cudaEventRecord failed");
715 nb->timers->didPairlistH2D[iloc] = true;
718 cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
719 &d_plist->nsci, &d_plist->sci_nalloc,
723 cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
724 &d_plist->ncj4, &d_plist->cj4_nalloc,
728 /* this call only allocates space on the device (no data is transferred) */
729 cu_realloc_buffered((void **)&d_plist->imask, NULL, sizeof(*d_plist->imask),
730 &d_plist->nimask, &d_plist->imask_nalloc,
731 h_plist->ncj4*c_nbnxnGpuClusterpairSplit,
734 cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
735 &d_plist->nexcl, &d_plist->excl_nalloc,
741 stat = cudaEventRecord(nb->timers->stop_pl_h2d[iloc], stream);
742 CU_RET_ERR(stat, "cudaEventRecord failed");
745 /* the next use of thist list we be the first one, so we need to prune */
746 d_plist->haveFreshList = true;
749 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t *nb,
750 const nbnxn_atomdata_t *nbatom)
752 cu_atomdata_t *adat = nb->atdat;
753 cudaStream_t ls = nb->stream[eintLocal];
755 /* only if we have a dynamic box */
756 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
758 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
759 SHIFTS * sizeof(*adat->shift_vec), ls);
760 adat->bShiftVecUploaded = true;
764 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
765 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
768 cu_atomdata_t *adat = nb->atdat;
769 cudaStream_t ls = nb->stream[eintLocal];
771 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
772 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
775 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
776 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
779 cu_atomdata_t *adat = nb->atdat;
780 cudaStream_t ls = nb->stream[eintLocal];
782 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
783 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
784 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
785 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
786 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
787 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
790 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
792 nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
793 /* clear shift force array and energies if the outputs were
794 used in the current step */
795 if (flags & GMX_FORCE_VIRIAL)
797 nbnxn_cuda_clear_e_fshift(nb);
801 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t *nb,
802 const struct nbnxn_atomdata_t *nbat)
807 bool bDoTime = nb->bDoTime;
808 cu_timers_t *timers = nb->timers;
809 cu_atomdata_t *d_atdat = nb->atdat;
810 cudaStream_t ls = nb->stream[eintLocal];
812 natoms = nbat->natoms;
817 /* time async copy */
818 stat = cudaEventRecord(timers->start_atdat, ls);
819 CU_RET_ERR(stat, "cudaEventRecord failed");
822 /* need to reallocate if we have to copy more atoms than the amount of space
823 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
824 if (natoms > d_atdat->nalloc)
826 nalloc = over_alloc_small(natoms);
828 /* free up first if the arrays have already been initialized */
829 if (d_atdat->nalloc != -1)
831 cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
832 cu_free_buffered(d_atdat->xq);
833 cu_free_buffered(d_atdat->atom_types);
834 cu_free_buffered(d_atdat->lj_comb);
837 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
838 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
839 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
840 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
841 if (useLjCombRule(nb->nbparam))
843 stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
844 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
848 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
849 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
852 d_atdat->nalloc = nalloc;
856 d_atdat->natoms = natoms;
857 d_atdat->natoms_local = nbat->natoms_local;
859 /* need to clear GPU f output if realloc happened */
862 nbnxn_cuda_clear_f(nb, nalloc);
865 if (useLjCombRule(nb->nbparam))
867 cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
868 natoms*sizeof(*d_atdat->lj_comb), ls);
872 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
873 natoms*sizeof(*d_atdat->atom_types), ls);
878 stat = cudaEventRecord(timers->stop_atdat, ls);
879 CU_RET_ERR(stat, "cudaEventRecord failed");
883 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t *nbparam,
884 const gmx_device_info_t *dev_info)
888 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
890 if (!c_disableCudaTextures)
892 /* Only device CC >= 3.0 (Kepler and later) support texture objects */
893 if (use_texobj(dev_info))
895 stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
896 CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
900 GMX_UNUSED_VALUE(dev_info);
901 stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
902 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
905 cu_free_buffered(nbparam->coulomb_tab);
909 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
912 cu_atomdata_t *atdat;
913 cu_nbparam_t *nbparam;
914 cu_plist_t *plist, *plist_nl;
923 nbparam = nb->nbparam;
924 plist = nb->plist[eintLocal];
925 plist_nl = nb->plist[eintNonlocal];
928 nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
930 stat = cudaEventDestroy(nb->nonlocal_done);
931 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
932 stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
933 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
937 stat = cudaEventDestroy(timers->start_atdat);
938 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
939 stat = cudaEventDestroy(timers->stop_atdat);
940 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
942 /* The non-local counters/stream (second in the array) are needed only with DD. */
943 for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
945 stat = cudaEventDestroy(timers->start_nb_k[i]);
946 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
947 stat = cudaEventDestroy(timers->stop_nb_k[i]);
948 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
950 stat = cudaEventDestroy(timers->start_prune_k[i]);
951 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_prune_k");
952 stat = cudaEventDestroy(timers->stop_prune_k[i]);
953 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_prune_k");
955 stat = cudaEventDestroy(timers->start_rollingPrune_k[i]);
956 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_rollingPrune_k");
957 stat = cudaEventDestroy(timers->stop_rollingPrune_k[i]);
958 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_rollingPrune_k");
960 stat = cudaEventDestroy(timers->start_pl_h2d[i]);
961 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
962 stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
963 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
965 stat = cudaStreamDestroy(nb->stream[i]);
966 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
968 stat = cudaEventDestroy(timers->start_nb_h2d[i]);
969 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
970 stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
971 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
973 stat = cudaEventDestroy(timers->start_nb_d2h[i]);
974 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
975 stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
976 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
980 if (!useLjCombRule(nb->nbparam))
982 if (!c_disableCudaTextures)
984 /* Only device CC >= 3.0 (Kepler and later) support texture objects */
985 if (use_texobj(nb->dev_info))
987 stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
988 CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
992 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
993 CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
996 cu_free_buffered(nbparam->nbfp);
999 if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
1001 if (!c_disableCudaTextures)
1003 /* Only device CC >= 3.0 (Kepler and later) support texture objects */
1004 if (use_texobj(nb->dev_info))
1006 stat = cudaDestroyTextureObject(nbparam->nbfp_comb_texobj);
1007 CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_comb_texobj failed");
1011 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_comb_texref());
1012 CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_comb_texref failed");
1015 cu_free_buffered(nbparam->nbfp_comb);
1018 stat = cudaFree(atdat->shift_vec);
1019 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
1020 stat = cudaFree(atdat->fshift);
1021 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
1023 stat = cudaFree(atdat->e_lj);
1024 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
1025 stat = cudaFree(atdat->e_el);
1026 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
1028 cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
1029 cu_free_buffered(atdat->xq);
1030 cu_free_buffered(atdat->atom_types, &atdat->ntypes);
1031 cu_free_buffered(atdat->lj_comb);
1033 cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
1034 cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
1035 cu_free_buffered(plist->imask, &plist->nimask, &plist->imask_nalloc);
1036 cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
1037 if (nb->bUseTwoStreams)
1039 cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
1040 cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
1041 cu_free_buffered(plist_nl->imask, &plist_nl->nimask, &plist_nl->imask_nalloc);
1042 cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
1048 if (nb->bUseTwoStreams)
1058 fprintf(debug, "Cleaned up CUDA data structures.\n");
1062 gmx_wallclock_gpu_t * nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
1064 return (nb != NULL && nb->bDoTime) ? nb->timings : NULL;
1067 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
1069 if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
1071 init_timings(nbv->gpu_nbv->timings);
1075 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
1078 gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
1082 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
1084 return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1085 (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));