1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2012, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
45 #include "gmx_fatal.h"
49 #include "types/nb_verlet.h"
50 #include "types/interaction_const.h"
51 #include "types/force_flags.h"
52 #include "../nbnxn_consts.h"
54 #include "nbnxn_cuda_types.h"
55 #include "../../gmxlib/cuda_tools/cudautils.cuh"
56 #include "nbnxn_cuda_data_mgmt.h"
57 #include "pmalloc_cuda.h"
58 #include "gpu_utils.h"
60 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
62 /* This is a heuristically determined parameter for the Fermi architecture for
63 * the minimum size of ci lists by multiplying this constant with the # of
64 * multiprocessors on the current device.
66 static unsigned int gpu_min_ci_balanced_factor = 40;
68 /* Functions from nbnxn_cuda.cu */
69 extern void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo);
70 extern const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_nbfp_texref();
71 extern const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_coulomb_tab_texref();
73 /* We should actually be using md_print_warn in md_logging.c,
74 * but we can't include mpi.h in CUDA code.
76 static void md_print_warn(FILE *fplog,
83 /* We should only print to stderr on the master node,
84 * in most cases fplog is only set on the master node, so this works.
87 fprintf(stderr, "\n");
88 vfprintf(stderr, fmt, ap);
89 fprintf(stderr, "\n");
94 vfprintf(fplog, fmt, ap);
102 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb);
105 /*! Tabulates the Ewald Coulomb force and initializes the size/scale
106 and the table GPU array. If called with an already allocated table,
107 it just re-uploads the table.
109 static void init_ewald_coulomb_force_table(cu_nbparam_t *nbp)
111 float *ftmp, *coul_tab;
116 tabsize = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
117 /* Subtract 2 iso 1 to avoid access out of range due to rounding */
118 tabscale = (tabsize - 2) / sqrt(nbp->rcoulomb_sq);
120 pmalloc((void**)&ftmp, tabsize*sizeof(*ftmp));
122 table_spline3_fill_ewald_lr(ftmp, NULL, NULL, tabsize,
123 1/tabscale, nbp->ewald_beta);
125 /* If the table pointer == NULL the table is generated the first time =>
126 the array pointer will be saved to nbparam and the texture is bound.
128 coul_tab = nbp->coulomb_tab;
129 if (coul_tab == NULL)
131 stat = cudaMalloc((void **)&coul_tab, tabsize*sizeof(*coul_tab));
132 CU_RET_ERR(stat, "cudaMalloc failed on coul_tab");
134 nbp->coulomb_tab = coul_tab;
136 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
137 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_coulomb_tab_texref(),
138 coul_tab, &cd, tabsize*sizeof(*coul_tab));
139 CU_RET_ERR(stat, "cudaBindTexture on coul_tab failed");
142 cu_copy_H2D(coul_tab, ftmp, tabsize*sizeof(*coul_tab));
144 nbp->coulomb_tab_size = tabsize;
145 nbp->coulomb_tab_scale = tabscale;
151 /*! Initializes the atomdata structure first time, it only gets filled at
153 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
158 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
159 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
160 ad->bShiftVecUploaded = false;
162 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
163 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
165 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
166 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
167 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
168 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
170 /* initialize to NULL poiters to data that is not allocated here and will
171 need reallocation in nbnxn_cuda_init_atomdata */
175 /* size -1 indicates that the respective array hasn't been initialized yet */
180 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
181 earlier GPUs, single or twin cut-off. */
182 static int pick_ewald_kernel_type(bool bTwinCut,
183 const cuda_dev_info_t *dev_info)
185 bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
188 /* Benchmarking/development environment variables to force the use of
189 analytical or tabulated Ewald kernel. */
190 bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
191 bForceTabulatedEwald = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
193 if (bForceAnalyticalEwald && bForceTabulatedEwald)
195 gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
196 "requested through environment variables.");
199 /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
200 if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
202 bUseAnalyticalEwald = true;
206 fprintf(debug, "Using analytical Ewald CUDA kernels\n");
211 bUseAnalyticalEwald = false;
215 fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
219 /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
220 forces it (use it for debugging/benchmarking only). */
221 if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
223 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
227 kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
234 /*! Initializes the nonbonded parameter data structure. */
235 static void init_nbparam(cu_nbparam_t *nbp,
236 const interaction_const_t *ic,
237 const nbnxn_atomdata_t *nbat,
238 const cuda_dev_info_t *dev_info)
243 ntypes = nbat->ntype;
245 nbp->ewald_beta = ic->ewaldcoeff;
246 nbp->sh_ewald = ic->sh_ewald;
247 nbp->epsfac = ic->epsfac;
248 nbp->two_k_rf = 2.0 * ic->k_rf;
249 nbp->c_rf = ic->c_rf;
250 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
251 nbp->rcoulomb_sq= ic->rcoulomb * ic->rcoulomb;
252 nbp->rlist_sq = ic->rlist * ic->rlist;
253 nbp->sh_invrc6 = ic->sh_invrc6;
255 if (ic->eeltype == eelCUT)
257 nbp->eeltype = eelCuCUT;
259 else if (EEL_RF(ic->eeltype))
261 nbp->eeltype = eelCuRF;
263 else if ((EEL_PME(ic->eeltype) || ic->eeltype==eelEWALD))
265 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
266 nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
270 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
271 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
274 /* generate table for PME */
275 nbp->coulomb_tab = NULL;
276 if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
278 init_ewald_coulomb_force_table(nbp);
281 nnbfp = 2*ntypes*ntypes;
282 stat = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
283 CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
284 cu_copy_H2D(nbp->nbfp, nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
286 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
287 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
288 nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
289 CU_RET_ERR(stat, "cudaBindTexture on nbfp failed");
292 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
293 * electrostatic type switching to twin cut-off (or back) if needed. */
294 void nbnxn_cuda_pme_loadbal_update_param(nbnxn_cuda_ptr_t cu_nb,
295 const interaction_const_t *ic)
297 cu_nbparam_t *nbp = cu_nb->nbparam;
299 nbp->rlist_sq = ic->rlist * ic->rlist;
300 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
301 nbp->ewald_beta = ic->ewaldcoeff;
303 nbp->eeltype = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
306 init_ewald_coulomb_force_table(cu_nb->nbparam);
309 /*! Initializes the pair list data structure. */
310 static void init_plist(cu_plist_t *pl)
312 /* initialize to NULL pointers to data that is not allocated here and will
313 need reallocation in nbnxn_cuda_init_pairlist */
318 /* size -1 indicates that the respective array hasn't been initialized yet */
325 pl->excl_nalloc = -1;
326 pl->bDoPrune = false;
329 /*! Initializes the timer data structure. */
330 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
333 int eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync: cudaEventDefault );
335 stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
336 CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
337 stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
338 CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
340 /* The non-local counters/stream (second in the array) are needed only with DD. */
341 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
343 stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
344 CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
345 stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
346 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
349 stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
350 CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
351 stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
352 CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
354 stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
355 CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
356 stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
357 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
359 stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
360 CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
361 stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
362 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
366 /*! Initializes the timings data structure. */
367 static void init_timings(wallclock_gpu_t *t)
376 for (i = 0; i < 2; i++)
378 for(j = 0; j < 2; j++)
380 t->ktime[i][j].t = 0.0;
381 t->ktime[i][j].c = 0;
386 /* Decide which kernel version to use (default or legacy) based on:
387 * - CUDA version used for compilation
388 * - non-bonded kernel selector environment variables
389 * - GPU architecture version
391 static int pick_nbnxn_kernel_version(FILE *fplog,
392 cuda_dev_info_t *devinfo)
394 bool bForceLegacyKernel, bForceDefaultKernel, bCUDA40, bCUDA32;
398 /* Legacy kernel (former k2), kept for backward compatibility as it is
399 faster than the default with CUDA 3.2/4.0 on Fermi (not on Kepler). */
400 bForceLegacyKernel = (getenv("GMX_CUDA_NB_LEGACY") != NULL);
401 /* default kernel (former k3). */
402 bForceDefaultKernel = (getenv("GMX_CUDA_NB_DEFAULT") != NULL);
404 if ((unsigned)(bForceLegacyKernel + bForceDefaultKernel) > 1)
406 gmx_fatal(FARGS, "Multiple CUDA non-bonded kernels requested; to manually pick a kernel set only one \n"
407 "of the following environment variables: \n"
408 "GMX_CUDA_NB_DEFAULT, GMX_CUDA_NB_LEGACY");
411 bCUDA32 = bCUDA40 = false;
412 #if CUDA_VERSION == 3200
414 sprintf(sbuf, "3.2");
415 #elif CUDA_VERSION == 4000
417 sprintf(sbuf, "4.0");
420 /* default is default ;) */
421 kver = eNbnxnCuKDefault;
423 /* Consider switching to legacy kernels only on Fermi */
424 if (devinfo->prop.major < 3 && (bCUDA32 || bCUDA40))
426 /* use legacy kernel unless something else is forced by an env. var */
427 if (bForceDefaultKernel)
430 "NOTE: CUDA %s compilation detected; with this compiler version the legacy\n"
431 " non-bonded kernels perform best. However, the default kernels were\n"
432 " selected by the GMX_CUDA_NB_DEFAULT environment variable.\n"
433 " For best performance upgrade your CUDA toolkit.\n",
438 kver = eNbnxnCuKLegacy;
443 /* issue note if the non-default kernel is forced by an env. var */
444 if (bForceLegacyKernel)
447 "NOTE: Legacy non-bonded CUDA kernels selected by the GMX_CUDA_NB_LEGACY\n"
448 " env. var. Consider using using the default kernels which should be faster!\n");
450 kver = eNbnxnCuKLegacy;
457 void nbnxn_cuda_init(FILE *fplog,
458 nbnxn_cuda_ptr_t *p_cu_nb,
459 gmx_gpu_info_t *gpu_info, int my_gpu_index,
460 gmx_bool bLocalAndNonlocal)
465 bool bStreamSync, bNoStreamSync, bTMPIAtomics, bX86, bOldDriver;
470 if (p_cu_nb == NULL) return;
474 snew(nb->nbparam, 1);
475 snew(nb->plist[eintLocal], 1);
476 if (bLocalAndNonlocal)
478 snew(nb->plist[eintNonlocal], 1);
481 nb->bUseTwoStreams = bLocalAndNonlocal;
484 snew(nb->timings, 1);
487 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
488 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
489 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
491 init_plist(nb->plist[eintLocal]);
493 /* local/non-local GPU streams */
494 stat = cudaStreamCreate(&nb->stream[eintLocal]);
495 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
496 if (nb->bUseTwoStreams)
498 init_plist(nb->plist[eintNonlocal]);
499 stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
500 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
503 /* init events for sychronization (timing disabled for performance reasons!) */
504 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
505 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
506 stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
507 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
509 /* set device info, just point it to the right GPU among the detected ones */
510 nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, my_gpu_index)];
512 /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
513 * (which increases with shorter time/step) caused by a known CUDA driver bug.
514 * To work around the issue we'll use an (admittedly fragile) memory polling
515 * waiting to preserve performance. This requires support for atomic
516 * operations and only works on x86/x86_64.
517 * With polling wait event-timing also needs to be disabled.
519 * The overhead is greatly reduced in API v5.0 drivers and the improvement
520 $ is independent of runtime version. Hence, with API v5.0 drivers and later
521 * we won't switch to polling.
523 * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
524 * ranks so we will also disable it in that case.
527 bStreamSync = getenv("GMX_CUDA_STREAMSYNC") != NULL;
528 bNoStreamSync = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
533 bTMPIAtomics = false;
536 #if defined(i386) || defined(__x86_64__)
542 if (bStreamSync && bNoStreamSync)
544 gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
547 stat = cudaDriverGetVersion(&cuda_drv_ver);
548 CU_RET_ERR(stat, "cudaDriverGetVersion failed");
550 bOldDriver = (cuda_drv_ver < 5000);
552 if ((nb->dev_info->prop.ECCEnabled == 1) && bOldDriver)
554 /* Polling wait should be used instead of cudaStreamSynchronize only if:
555 * - ECC is ON & driver is old (checked above),
556 * - we're on x86/x86_64,
557 * - atomics are available, and
558 * - GPUs are not being shared.
560 bool bShouldUsePollSync = (bX86 && bTMPIAtomics && !gpu_info->bDevShare);
564 nb->bUseStreamSync = true;
566 /* only warn if polling should be used */
567 if (bShouldUsePollSync)
570 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
571 " cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
576 nb->bUseStreamSync = !bShouldUsePollSync;
578 if (bShouldUsePollSync)
581 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
582 " cause performance loss. Switching to the alternative polling GPU wait.\n"
583 " If you encounter issues, switch back to standard GPU waiting by setting\n"
584 " the GMX_CUDA_STREAMSYNC environment variable.\n");
588 /* Tell the user that the ECC+old driver combination can be bad */
590 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0.\n"
591 " A known bug in this driver version can cause performance loss.\n"
592 " However, the polling wait workaround can not be used because\n%s\n"
593 " Consider updating the driver or turning ECC off.",
594 (bX86 && bTMPIAtomics) ?
595 " GPU(s) are being oversubscribed." :
596 " atomic operations are not supported by the platform/CPU+compiler.");
597 md_print_warn(fplog, sbuf);
605 nb->bUseStreamSync = false;
608 "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
612 /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
613 nb->bUseStreamSync = true;
617 /* CUDA timing disabled as event timers don't work:
618 - with multiple streams = domain-decomposition;
619 - with the polling waiting hack (without cudaStreamSynchronize);
620 - when turned off by GMX_DISABLE_CUDA_TIMING.
622 nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
623 (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
627 init_timers(nb->timers, nb->bUseTwoStreams);
628 init_timings(nb->timings);
631 /* set the kernel type for the current GPU */
632 nb->kernel_ver = pick_nbnxn_kernel_version(fplog, nb->dev_info);
633 /* pick L1 cache configuration */
634 nbnxn_cuda_set_cacheconfig(nb->dev_info);
640 fprintf(debug, "Initialized CUDA data structures.\n");
644 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t cu_nb,
645 const interaction_const_t *ic,
646 const nonbonded_verlet_group_t *nbv_group)
648 init_atomdata_first(cu_nb->atdat, nbv_group[0].nbat->ntype);
649 init_nbparam(cu_nb->nbparam, ic, nbv_group[0].nbat, cu_nb->dev_info);
651 /* clear energy and shift force outputs */
652 nbnxn_cuda_clear_e_fshift(cu_nb);
655 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t cu_nb,
656 const nbnxn_pairlist_t *h_plist,
661 bool bDoTime = cu_nb->bDoTime;
662 cudaStream_t stream = cu_nb->stream[iloc];
663 cu_plist_t *d_plist = cu_nb->plist[iloc];
665 if (d_plist->na_c < 0)
667 d_plist->na_c = h_plist->na_ci;
671 if (d_plist->na_c != h_plist->na_ci)
673 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
674 d_plist->na_c, h_plist->na_ci);
681 stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
682 CU_RET_ERR(stat, "cudaEventRecord failed");
685 cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
686 &d_plist->nsci, &d_plist->sci_nalloc,
690 cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
691 &d_plist->ncj4, &d_plist->cj4_nalloc,
695 cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
696 &d_plist->nexcl, &d_plist->excl_nalloc,
702 stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
703 CU_RET_ERR(stat, "cudaEventRecord failed");
706 /* need to prune the pair list during the next step */
707 d_plist->bDoPrune = true;
710 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t cu_nb,
711 const nbnxn_atomdata_t *nbatom)
713 cu_atomdata_t *adat = cu_nb->atdat;
714 cudaStream_t ls = cu_nb->stream[eintLocal];
716 /* only if we have a dynamic box */
717 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
719 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
720 SHIFTS * sizeof(*adat->shift_vec), ls);
721 adat->bShiftVecUploaded = true;
725 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
726 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
729 cu_atomdata_t *adat = cu_nb->atdat;
730 cudaStream_t ls = cu_nb->stream[eintLocal];
732 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
733 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
736 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
737 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
740 cu_atomdata_t *adat = cu_nb->atdat;
741 cudaStream_t ls = cu_nb->stream[eintLocal];
743 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
744 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
745 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
746 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
747 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
748 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
751 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
753 nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
754 /* clear shift force array and energies if the outputs were
755 used in the current step */
756 if (flags & GMX_FORCE_VIRIAL)
758 nbnxn_cuda_clear_e_fshift(cu_nb);
762 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t cu_nb,
763 const nbnxn_atomdata_t *nbat)
768 bool bDoTime = cu_nb->bDoTime;
769 cu_timers_t *timers = cu_nb->timers;
770 cu_atomdata_t *d_atdat = cu_nb->atdat;
771 cudaStream_t ls = cu_nb->stream[eintLocal];
773 natoms = nbat->natoms;
778 /* time async copy */
779 stat = cudaEventRecord(timers->start_atdat, ls);
780 CU_RET_ERR(stat, "cudaEventRecord failed");
783 /* need to reallocate if we have to copy more atoms than the amount of space
784 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
785 if (natoms > d_atdat->nalloc)
787 nalloc = over_alloc_small(natoms);
789 /* free up first if the arrays have already been initialized */
790 if (d_atdat->nalloc != -1)
792 cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
793 cu_free_buffered(d_atdat->xq);
794 cu_free_buffered(d_atdat->atom_types);
797 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
798 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
799 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
800 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
802 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
803 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
805 d_atdat->nalloc = nalloc;
809 d_atdat->natoms = natoms;
810 d_atdat->natoms_local = nbat->natoms_local;
812 /* need to clear GPU f output if realloc happened */
815 nbnxn_cuda_clear_f(cu_nb, nalloc);
818 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
819 natoms*sizeof(*d_atdat->atom_types), ls);
823 stat = cudaEventRecord(timers->stop_atdat, ls);
824 CU_RET_ERR(stat, "cudaEventRecord failed");
828 void nbnxn_cuda_free(FILE *fplog, nbnxn_cuda_ptr_t cu_nb)
831 cu_atomdata_t *atdat;
832 cu_nbparam_t *nbparam;
833 cu_plist_t *plist, *plist_nl;
836 if (cu_nb == NULL) return;
838 atdat = cu_nb->atdat;
839 nbparam = cu_nb->nbparam;
840 plist = cu_nb->plist[eintLocal];
841 plist_nl = cu_nb->plist[eintNonlocal];
842 timers = cu_nb->timers;
844 if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
846 stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
847 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
848 cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
851 stat = cudaEventDestroy(cu_nb->nonlocal_done);
852 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
853 stat = cudaEventDestroy(cu_nb->misc_ops_done);
854 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
858 stat = cudaEventDestroy(timers->start_atdat);
859 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
860 stat = cudaEventDestroy(timers->stop_atdat);
861 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
863 /* The non-local counters/stream (second in the array) are needed only with DD. */
864 for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
866 stat = cudaEventDestroy(timers->start_nb_k[i]);
867 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
868 stat = cudaEventDestroy(timers->stop_nb_k[i]);
869 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
871 stat = cudaEventDestroy(timers->start_pl_h2d[i]);
872 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
873 stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
874 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
876 stat = cudaStreamDestroy(cu_nb->stream[i]);
877 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
879 stat = cudaEventDestroy(timers->start_nb_h2d[i]);
880 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
881 stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
882 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
884 stat = cudaEventDestroy(timers->start_nb_d2h[i]);
885 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
886 stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
887 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
891 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
892 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
893 cu_free_buffered(nbparam->nbfp);
895 stat = cudaFree(atdat->shift_vec);
896 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
897 stat = cudaFree(atdat->fshift);
898 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
900 stat = cudaFree(atdat->e_lj);
901 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
902 stat = cudaFree(atdat->e_el);
903 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
905 cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
906 cu_free_buffered(atdat->xq);
907 cu_free_buffered(atdat->atom_types, &atdat->ntypes);
909 cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
910 cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
911 cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
912 if (cu_nb->bUseTwoStreams)
914 cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
915 cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
916 cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
922 if (cu_nb->bUseTwoStreams)
927 sfree(cu_nb->timings);
932 fprintf(debug, "Cleaned up CUDA data structures.\n");
936 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
939 cudaStream_t stream = cu_nb->stream[iloc];
941 stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
942 CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
945 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
947 return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
950 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
954 init_timings(cu_nb->timings);
958 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
960 return cu_nb != NULL ?
961 gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;