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
43 #include "gmx_fatal.h"
47 #include "types/nb_verlet.h"
48 #include "types/interaction_const.h"
49 #include "types/force_flags.h"
50 #include "../nbnxn_consts.h"
52 #include "nbnxn_cuda_types.h"
53 #include "../../gmxlib/cuda_tools/cudautils.cuh"
54 #include "nbnxn_cuda_data_mgmt.h"
55 #include "pmalloc_cuda.h"
56 #include "gpu_utils.h"
58 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
60 /* This is a heuristically determined parameter for the Fermi architecture for
61 * the minimum size of ci lists by multiplying this constant with the # of
62 * multiprocessors on the current device.
64 static unsigned int gpu_min_ci_balanced_factor = 40;
66 /* Functions from nbnxn_cuda.cu */
67 extern void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo);
68 extern const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_nbfp_texref();
69 extern const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_coulomb_tab_texref();
71 /* We should actually be using md_print_warn in md_logging.c,
72 * but we can't include mpi.h in CUDA code.
74 static void md_print_warn(FILE *fplog, const char *buf)
78 /* We should only print to stderr on the master node,
79 * in most cases fplog is only set on the master node, so this works.
81 fprintf(stderr, "\n%s\n", buf);
82 fprintf(fplog, "\n%s\n", buf);
87 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb);
90 /*! Tabulates the Ewald Coulomb force and initializes the size/scale
91 and the table GPU array. If called with an already allocated table,
92 it just re-uploads the table.
94 static void init_ewald_coulomb_force_table(cu_nbparam_t *nbp)
96 float *ftmp, *coul_tab;
101 tabsize = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
102 /* Subtract 2 iso 1 to avoid access out of range due to rounding */
103 tabscale = (tabsize - 2) / sqrt(nbp->rcoulomb_sq);
105 pmalloc((void**)&ftmp, tabsize*sizeof(*ftmp));
107 table_spline3_fill_ewald_lr(ftmp, NULL, NULL, tabsize,
108 1/tabscale, nbp->ewald_beta);
110 /* If the table pointer == NULL the table is generated the first time =>
111 the array pointer will be saved to nbparam and the texture is bound.
113 coul_tab = nbp->coulomb_tab;
114 if (coul_tab == NULL)
116 stat = cudaMalloc((void **)&coul_tab, tabsize*sizeof(*coul_tab));
117 CU_RET_ERR(stat, "cudaMalloc failed on coul_tab");
119 nbp->coulomb_tab = coul_tab;
121 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
122 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_coulomb_tab_texref(),
123 coul_tab, &cd, tabsize*sizeof(*coul_tab));
124 CU_RET_ERR(stat, "cudaBindTexture on coul_tab failed");
127 cu_copy_H2D(coul_tab, ftmp, tabsize*sizeof(*coul_tab));
129 nbp->coulomb_tab_size = tabsize;
130 nbp->coulomb_tab_scale = tabscale;
136 /*! Initializes the atomdata structure first time, it only gets filled at
138 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
143 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
144 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
145 ad->bShiftVecUploaded = false;
147 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
148 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
150 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
151 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
152 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
153 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
155 /* initialize to NULL poiters to data that is not allocated here and will
156 need reallocation in nbnxn_cuda_init_atomdata */
160 /* size -1 indicates that the respective array hasn't been initialized yet */
165 /*! Initializes the nonbonded parameter data structure. */
166 static void init_nbparam(cu_nbparam_t *nbp,
167 const interaction_const_t *ic,
168 const nonbonded_verlet_t *nbv)
173 ntypes = nbv->grp[0].nbat->ntype;
175 nbp->ewald_beta = ic->ewaldcoeff;
176 nbp->sh_ewald = ic->sh_ewald;
177 nbp->epsfac = ic->epsfac;
178 nbp->two_k_rf = 2.0 * ic->k_rf;
179 nbp->c_rf = ic->c_rf;
180 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
181 nbp->rcoulomb_sq= ic->rcoulomb * ic->rcoulomb;
182 nbp->rlist_sq = ic->rlist * ic->rlist;
183 nbp->sh_invrc6 = ic->sh_invrc6;
185 if (ic->eeltype == eelCUT)
187 nbp->eeltype = eelCuCUT;
189 else if (EEL_RF(ic->eeltype))
191 nbp->eeltype = eelCuRF;
193 else if ((EEL_PME(ic->eeltype) || ic->eeltype==eelEWALD))
195 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off, unless
196 forced by the env. var. (used only for benchmarking). */
197 if (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL)
199 nbp->eeltype = eelCuEWALD;
203 nbp->eeltype = eelCuEWALD_TWIN;
208 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
209 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
212 /* generate table for PME */
213 if (nbp->eeltype == eelCuEWALD)
215 nbp->coulomb_tab = NULL;
216 init_ewald_coulomb_force_table(nbp);
219 nnbfp = 2*ntypes*ntypes;
220 stat = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
221 CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
222 cu_copy_H2D(nbp->nbfp, nbv->grp[0].nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
224 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
225 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
226 nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
227 CU_RET_ERR(stat, "cudaBindTexture on nbfp failed");
230 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
231 * electrostatic type switching to twin cut-off (or back) if needed. */
232 void nbnxn_cuda_pme_loadbal_update_param(nbnxn_cuda_ptr_t cu_nb,
233 const interaction_const_t *ic)
235 cu_nbparam_t *nbp = cu_nb->nbparam;
237 nbp->rlist_sq = ic->rlist * ic->rlist;
238 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
239 nbp->ewald_beta = ic->ewaldcoeff;
241 /* When switching to/from twin cut-off, the electrostatics type needs updating.
242 (The env. var. that forces twin cut-off is for benchmarking only!) */
243 if (ic->rcoulomb == ic->rvdw &&
244 getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL)
246 nbp->eeltype = eelCuEWALD;
250 nbp->eeltype = eelCuEWALD_TWIN;
253 init_ewald_coulomb_force_table(cu_nb->nbparam);
256 /*! Initializes the pair list data structure. */
257 static void init_plist(cu_plist_t *pl)
259 /* initialize to NULL pointers to data that is not allocated here and will
260 need reallocation in nbnxn_cuda_init_pairlist */
265 /* size -1 indicates that the respective array hasn't been initialized yet */
272 pl->excl_nalloc = -1;
273 pl->bDoPrune = false;
276 /*! Initializes the timer data structure. */
277 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
280 int eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync: cudaEventDefault );
282 stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
283 CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
284 stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
285 CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
287 /* The non-local counters/stream (second in the array) are needed only with DD. */
288 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
290 stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
291 CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
292 stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
293 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
296 stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
297 CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
298 stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
299 CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
301 stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
302 CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
303 stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
304 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
306 stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
307 CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
308 stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
309 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
313 /*! Initializes the timings data structure. */
314 static void init_timings(wallclock_gpu_t *t)
323 for (i = 0; i < 2; i++)
325 for(j = 0; j < 2; j++)
327 t->ktime[i][j].t = 0.0;
328 t->ktime[i][j].c = 0;
333 /* Decide which kernel version to use (default or legacy) based on:
335 * - non-bonded kernel selector environment variables
336 * - GPU SM version TODO ???
338 static int pick_nbnxn_kernel_version()
340 bool bLegacyKernel, bDefaultKernel, bCUDA40, bCUDA32;
344 /* legacy kernel (former k2), kept for now for backward compatibility,
345 faster than the default with CUDA 3.2/4.0 (TODO: on Kepler?). */
346 bLegacyKernel = (getenv("GMX_CUDA_NB_LEGACY") != NULL);
347 /* default kernel (former k3). */
348 bDefaultKernel = (getenv("GMX_CUDA_NB_DEFAULT") != NULL);
350 if ((unsigned)(bLegacyKernel + bDefaultKernel) > 1)
352 gmx_fatal(FARGS, "Multiple CUDA non-bonded kernels requested; to manually pick a kernel set only one \n"
353 "of the following environment variables: \n"
354 "GMX_CUDA_NB_DEFAULT, GMX_CUDA_NB_LEGACY");
357 bCUDA32 = bCUDA40 = false;
358 #if CUDA_VERSION == 3200
360 sprintf(sbuf, "3.2");
361 #elif CUDA_VERSION == 4000
363 sprintf(sbuf, "4.0");
366 /* default is default ;) */
367 kver = eNbnxnCuKDefault;
369 if (bCUDA32 || bCUDA40)
371 /* use legacy kernel unless something else is forced by an env. var */
375 "\nNOTE: CUDA %s compilation detected; with this compiler version the legacy\n"
376 " non-bonded kernels perform best. However, the default kernels were\n"
377 " selected by the GMX_CUDA_NB_DEFAULT environment variable.\n"
378 " For best performance upgrade your CUDA toolkit.",
383 kver = eNbnxnCuKLegacy;
388 /* issue not if the non-default kernel is forced by an env. var */
392 "\nNOTE: Legacy non-bonded CUDA kernels were selected by the GMX_CUDA_NB_LEGACY\n"
393 " env. var. Consider using using the default kernels which should be faster!\n");
395 kver = eNbnxnCuKLegacy;
402 void nbnxn_cuda_init(FILE *fplog,
403 nbnxn_cuda_ptr_t *p_cu_nb,
404 gmx_gpu_info_t *gpu_info, int my_gpu_index,
405 gmx_bool bLocalAndNonlocal)
410 bool bStreamSync, bNoStreamSync, bTMPIAtomics, bX86, bOldDriver;
415 if (p_cu_nb == NULL) return;
419 snew(nb->nbparam, 1);
420 snew(nb->plist[eintLocal], 1);
421 if (bLocalAndNonlocal)
423 snew(nb->plist[eintNonlocal], 1);
426 nb->bUseTwoStreams = bLocalAndNonlocal;
429 snew(nb->timings, 1);
432 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
433 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
434 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
436 init_plist(nb->plist[eintLocal]);
438 /* local/non-local GPU streams */
439 stat = cudaStreamCreate(&nb->stream[eintLocal]);
440 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
441 if (nb->bUseTwoStreams)
443 init_plist(nb->plist[eintNonlocal]);
444 stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
445 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
448 /* init events for sychronization (timing disabled for performance reasons!) */
449 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
450 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
451 stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
452 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
454 /* set device info, just point it to the right GPU among the detected ones */
455 nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, my_gpu_index)];
457 /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
458 * (which increases with shorter time/step) caused by a known CUDA driver bug.
459 * To work around the issue we'll use an (admittedly fragile) memory polling
460 * waiting to preserve performance. This requires support for atomic
461 * operations and only works on x86/x86_64.
462 * With polling wait event-timing also needs to be disabled.
464 * The overhead is greatly reduced in API v5.0 drivers and the improvement
465 $ is independent of runtime version. Hence, with API v5.0 drivers and later
466 * we won't switch to polling.
468 * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
469 * ranks so we will also disable it in that case.
472 bStreamSync = getenv("GMX_CUDA_STREAMSYNC") != NULL;
473 bNoStreamSync = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
478 bTMPIAtomics = false;
481 #if defined(i386) || defined(__x86_64__)
487 if (bStreamSync && bNoStreamSync)
489 gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
492 stat = cudaDriverGetVersion(&cuda_drv_ver);
493 CU_RET_ERR(stat, "cudaDriverGetVersion failed");
494 bOldDriver = (cuda_drv_ver < 5000);
496 if (nb->dev_info->prop.ECCEnabled == 1)
500 nb->bUseStreamSync = true;
502 /* only warn if polling should be used */
503 if (bOldDriver && !gpu_info->bDevShare)
506 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
507 " cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
512 /* Can/should turn of cudaStreamSynchronize wait only if
513 * - we're on x86/x86_64
514 * - atomics are available
515 * - GPUs are not being shared
516 * - and driver is old. */
518 (bX86 && bTMPIAtomics && !gpu_info->bDevShare && bOldDriver) ?
521 if (nb->bUseStreamSync)
524 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
525 " cause performance loss. Switching to the alternative polling GPU waiting.\n"
526 " If you encounter issues, switch back to standard GPU waiting by setting\n"
527 " the GMX_CUDA_STREAMSYNC environment variable.\n");
531 /* Tell the user that the ECC+old driver combination can be bad */
533 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0. A bug in this\n"
534 " driver can cause performance loss.\n"
535 " However, the polling waiting workaround can not be used because\n%s\n"
536 " Consider updating the driver or turning ECC off.",
537 (!bX86 || !bTMPIAtomics) ?
538 " atomic operations are not supported by the platform/CPU+compiler." :
539 " GPU(s) are being oversubscribed.");
540 md_print_warn(fplog, sbuf);
548 nb->bUseStreamSync = false;
551 "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
555 /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
556 nb->bUseStreamSync = true;
560 /* CUDA timing disabled as event timers don't work:
561 - with multiple streams = domain-decomposition;
562 - with the polling waiting hack (without cudaStreamSynchronize);
563 - when turned off by GMX_DISABLE_CUDA_TIMING.
565 nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
566 (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
570 init_timers(nb->timers, nb->bUseTwoStreams);
571 init_timings(nb->timings);
574 /* set the kernel type for the current GPU */
575 nb->kernel_ver = pick_nbnxn_kernel_version();
576 /* pick L1 cache configuration */
577 nbnxn_cuda_set_cacheconfig(nb->dev_info);
583 fprintf(debug, "Initialized CUDA data structures.\n");
587 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t cu_nb,
588 const interaction_const_t *ic,
589 const nonbonded_verlet_t *nbv)
591 init_atomdata_first(cu_nb->atdat, nbv->grp[0].nbat->ntype);
592 init_nbparam(cu_nb->nbparam, ic, nbv);
594 /* clear energy and shift force outputs */
595 nbnxn_cuda_clear_e_fshift(cu_nb);
598 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t cu_nb,
599 const nbnxn_pairlist_t *h_plist,
604 bool bDoTime = cu_nb->bDoTime;
605 cudaStream_t stream = cu_nb->stream[iloc];
606 cu_plist_t *d_plist = cu_nb->plist[iloc];
608 if (d_plist->na_c < 0)
610 d_plist->na_c = h_plist->na_ci;
614 if (d_plist->na_c != h_plist->na_ci)
616 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
617 d_plist->na_c, h_plist->na_ci);
624 stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
625 CU_RET_ERR(stat, "cudaEventRecord failed");
628 cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
629 &d_plist->nsci, &d_plist->sci_nalloc,
633 cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
634 &d_plist->ncj4, &d_plist->cj4_nalloc,
638 cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
639 &d_plist->nexcl, &d_plist->excl_nalloc,
645 stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
646 CU_RET_ERR(stat, "cudaEventRecord failed");
649 /* need to prune the pair list during the next step */
650 d_plist->bDoPrune = true;
653 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t cu_nb,
654 const nbnxn_atomdata_t *nbatom)
656 cu_atomdata_t *adat = cu_nb->atdat;
657 cudaStream_t ls = cu_nb->stream[eintLocal];
659 /* only if we have a dynamic box */
660 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
662 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
663 SHIFTS * sizeof(*adat->shift_vec), ls);
664 adat->bShiftVecUploaded = true;
668 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
669 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
672 cu_atomdata_t *adat = cu_nb->atdat;
673 cudaStream_t ls = cu_nb->stream[eintLocal];
675 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
676 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
679 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
680 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
683 cu_atomdata_t *adat = cu_nb->atdat;
684 cudaStream_t ls = cu_nb->stream[eintLocal];
686 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
687 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
688 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
689 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
690 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
691 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
694 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
696 nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
697 /* clear shift force array and energies if the outputs were
698 used in the current step */
699 if (flags & GMX_FORCE_VIRIAL)
701 nbnxn_cuda_clear_e_fshift(cu_nb);
705 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t cu_nb,
706 const nbnxn_atomdata_t *nbat)
711 bool bDoTime = cu_nb->bDoTime;
712 cu_timers_t *timers = cu_nb->timers;
713 cu_atomdata_t *d_atdat = cu_nb->atdat;
714 cudaStream_t ls = cu_nb->stream[eintLocal];
716 natoms = nbat->natoms;
721 /* time async copy */
722 stat = cudaEventRecord(timers->start_atdat, ls);
723 CU_RET_ERR(stat, "cudaEventRecord failed");
726 /* need to reallocate if we have to copy more atoms than the amount of space
727 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
728 if (natoms > d_atdat->nalloc)
730 nalloc = over_alloc_small(natoms);
732 /* free up first if the arrays have already been initialized */
733 if (d_atdat->nalloc != -1)
735 cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
736 cu_free_buffered(d_atdat->xq);
737 cu_free_buffered(d_atdat->atom_types);
740 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
741 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
742 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
743 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
745 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
746 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
748 d_atdat->nalloc = nalloc;
752 d_atdat->natoms = natoms;
753 d_atdat->natoms_local = nbat->natoms_local;
755 /* need to clear GPU f output if realloc happened */
758 nbnxn_cuda_clear_f(cu_nb, nalloc);
761 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
762 natoms*sizeof(*d_atdat->atom_types), ls);
766 stat = cudaEventRecord(timers->stop_atdat, ls);
767 CU_RET_ERR(stat, "cudaEventRecord failed");
771 void nbnxn_cuda_free(FILE *fplog, nbnxn_cuda_ptr_t cu_nb)
774 cu_atomdata_t *atdat;
775 cu_nbparam_t *nbparam;
776 cu_plist_t *plist, *plist_nl;
779 if (cu_nb == NULL) return;
781 atdat = cu_nb->atdat;
782 nbparam = cu_nb->nbparam;
783 plist = cu_nb->plist[eintLocal];
784 plist_nl = cu_nb->plist[eintNonlocal];
785 timers = cu_nb->timers;
787 if (nbparam->eeltype == eelCuEWALD || nbparam->eeltype == eelCuEWALD_TWIN)
789 stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
790 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
791 cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
794 stat = cudaEventDestroy(cu_nb->nonlocal_done);
795 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
796 stat = cudaEventDestroy(cu_nb->misc_ops_done);
797 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
801 stat = cudaEventDestroy(timers->start_atdat);
802 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
803 stat = cudaEventDestroy(timers->stop_atdat);
804 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
806 /* The non-local counters/stream (second in the array) are needed only with DD. */
807 for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
809 stat = cudaEventDestroy(timers->start_nb_k[i]);
810 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
811 stat = cudaEventDestroy(timers->stop_nb_k[i]);
812 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
814 stat = cudaEventDestroy(timers->start_pl_h2d[i]);
815 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
816 stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
817 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
819 stat = cudaStreamDestroy(cu_nb->stream[i]);
820 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
822 stat = cudaEventDestroy(timers->start_nb_h2d[i]);
823 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
824 stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
825 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
827 stat = cudaEventDestroy(timers->start_nb_d2h[i]);
828 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
829 stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
830 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
834 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
835 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
836 cu_free_buffered(nbparam->nbfp);
838 stat = cudaFree(atdat->shift_vec);
839 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
840 stat = cudaFree(atdat->fshift);
841 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
843 stat = cudaFree(atdat->e_lj);
844 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
845 stat = cudaFree(atdat->e_el);
846 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
848 cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
849 cu_free_buffered(atdat->xq);
850 cu_free_buffered(atdat->atom_types, &atdat->ntypes);
852 cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
853 cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
854 cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
855 if (cu_nb->bUseTwoStreams)
857 cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
858 cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
859 cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
864 fprintf(debug, "Cleaned up CUDA data structures.\n");
868 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
871 cudaStream_t stream = cu_nb->stream[iloc];
873 stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
874 CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
877 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
879 return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
882 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
886 init_timings(cu_nb->timings);
890 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
892 return cu_nb != NULL ?
893 gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;