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();
72 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb);
75 /*! Tabulates the Ewald Coulomb force and initializes the size/scale
76 and the table GPU array. If called with an already allocated table,
77 it just re-uploads the table.
79 static void init_ewald_coulomb_force_table(cu_nbparam_t *nbp)
81 float *ftmp, *coul_tab;
86 tabsize = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
87 /* Subtract 2 iso 1 to avoid access out of range due to rounding */
88 tabscale = (tabsize - 2) / sqrt(nbp->rcoulomb_sq);
90 pmalloc((void**)&ftmp, tabsize*sizeof(*ftmp));
92 table_spline3_fill_ewald_lr(ftmp, NULL, NULL, tabsize,
93 1/tabscale, nbp->ewald_beta);
95 /* If the table pointer == NULL the table is generated the first time =>
96 the array pointer will be saved to nbparam and the texture is bound.
98 coul_tab = nbp->coulomb_tab;
101 stat = cudaMalloc((void **)&coul_tab, tabsize*sizeof(*coul_tab));
102 CU_RET_ERR(stat, "cudaMalloc failed on coul_tab");
104 nbp->coulomb_tab = coul_tab;
106 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
107 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_coulomb_tab_texref(),
108 coul_tab, &cd, tabsize*sizeof(*coul_tab));
109 CU_RET_ERR(stat, "cudaBindTexture on coul_tab failed");
112 cu_copy_H2D(coul_tab, ftmp, tabsize*sizeof(*coul_tab));
114 nbp->coulomb_tab_size = tabsize;
115 nbp->coulomb_tab_scale = tabscale;
121 /*! Initializes the atomdata structure first time, it only gets filled at
123 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
128 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
129 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
130 ad->bShiftVecUploaded = false;
132 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
133 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
135 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
136 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
137 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
138 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
140 /* initialize to NULL poiters to data that is not allocated here and will
141 need reallocation in nbnxn_cuda_init_atomdata */
145 /* size -1 indicates that the respective array hasn't been initialized yet */
150 /*! Initializes the nonbonded parameter data structure. */
151 static void init_nbparam(cu_nbparam_t *nbp,
152 const interaction_const_t *ic,
153 const nonbonded_verlet_t *nbv)
158 ntypes = nbv->grp[0].nbat->ntype;
160 nbp->ewald_beta = ic->ewaldcoeff;
161 nbp->sh_ewald = ic->sh_ewald;
162 nbp->epsfac = ic->epsfac;
163 nbp->two_k_rf = 2.0 * ic->k_rf;
164 nbp->c_rf = ic->c_rf;
165 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
166 nbp->rcoulomb_sq= ic->rcoulomb * ic->rcoulomb;
167 nbp->rlist_sq = ic->rlist * ic->rlist;
168 nbp->sh_invrc6 = ic->sh_invrc6;
170 if (ic->eeltype == eelCUT)
172 nbp->eeltype = eelCuCUT;
174 else if (EEL_RF(ic->eeltype))
176 nbp->eeltype = eelCuRF;
178 else if ((EEL_PME(ic->eeltype) || ic->eeltype==eelEWALD))
180 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off, unless
181 forced by the env. var. (used only for benchmarking). */
182 if (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL)
184 nbp->eeltype = eelCuEWALD;
188 nbp->eeltype = eelCuEWALD_TWIN;
193 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
194 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
197 /* generate table for PME */
198 if (nbp->eeltype == eelCuEWALD)
200 nbp->coulomb_tab = NULL;
201 init_ewald_coulomb_force_table(nbp);
204 nnbfp = 2*ntypes*ntypes;
205 stat = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
206 CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
207 cu_copy_H2D(nbp->nbfp, nbv->grp[0].nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
209 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
210 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
211 nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
212 CU_RET_ERR(stat, "cudaBindTexture on nbfp failed");
215 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
216 * electrostatic type switching to twin cut-off (or back) if needed. */
217 void nbnxn_cuda_pme_loadbal_update_param(nbnxn_cuda_ptr_t cu_nb,
218 const interaction_const_t *ic)
220 cu_nbparam_t *nbp = cu_nb->nbparam;
222 nbp->rlist_sq = ic->rlist * ic->rlist;
223 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
224 nbp->ewald_beta = ic->ewaldcoeff;
226 /* When switching to/from twin cut-off, the electrostatics type needs updating.
227 (The env. var. that forces twin cut-off is for benchmarking only!) */
228 if (ic->rcoulomb == ic->rvdw &&
229 getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL)
231 nbp->eeltype = eelCuEWALD;
235 nbp->eeltype = eelCuEWALD_TWIN;
238 init_ewald_coulomb_force_table(cu_nb->nbparam);
241 /*! Initializes the pair list data structure. */
242 static void init_plist(cu_plist_t *pl)
244 /* initialize to NULL pointers to data that is not allocated here and will
245 need reallocation in nbnxn_cuda_init_pairlist */
250 /* size -1 indicates that the respective array hasn't been initialized yet */
257 pl->excl_nalloc = -1;
258 pl->bDoPrune = false;
261 /*! Initializes the timer data structure. */
262 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
265 int eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync: cudaEventDefault );
267 stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
268 CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
269 stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
270 CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
272 /* The non-local counters/stream (second in the array) are needed only with DD. */
273 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
275 stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
276 CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
277 stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
278 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
281 stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
282 CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
283 stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
284 CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
286 stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
287 CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
288 stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
289 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
291 stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
292 CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
293 stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
294 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
298 /*! Initializes the timings data structure. */
299 static void init_timings(wallclock_gpu_t *t)
308 for (i = 0; i < 2; i++)
310 for(j = 0; j < 2; j++)
312 t->ktime[i][j].t = 0.0;
313 t->ktime[i][j].c = 0;
318 /* Decide which kernel version to use (default or legacy) based on:
320 * - non-bonded kernel selector environment variables
321 * - GPU SM version TODO ???
323 static int pick_nbnxn_kernel_version()
325 bool bLegacyKernel, bDefaultKernel, bCUDA40, bCUDA32;
329 /* legacy kernel (former k2), kept for now for backward compatibility,
330 faster than the default with CUDA 3.2/4.0 (TODO: on Kepler?). */
331 bLegacyKernel = (getenv("GMX_CUDA_NB_LEGACY") != NULL);
332 /* default kernel (former k3). */
333 bDefaultKernel = (getenv("GMX_CUDA_NB_DEFAULT") != NULL);
335 if ((unsigned)(bLegacyKernel + bDefaultKernel) > 1)
337 gmx_fatal(FARGS, "Multiple CUDA non-bonded kernels requested; to manually pick a kernel set only one \n"
338 "of the following environment variables: \n"
339 "GMX_CUDA_NB_DEFAULT, GMX_CUDA_NB_LEGACY");
342 bCUDA32 = bCUDA40 = false;
343 #if CUDA_VERSION == 3200
345 sprintf(sbuf, "3.2");
346 #elif CUDA_VERSION == 4000
348 sprintf(sbuf, "4.0");
351 /* default is default ;) */
352 kver = eNbnxnCuKDefault;
354 if (bCUDA32 || bCUDA40)
356 /* use legacy kernel unless something else is forced by an env. var */
360 "\nNOTE: CUDA %s compilation detected; with this compiler version the legacy\n"
361 " non-bonded kernels perform best. However, the default kernels were\n"
362 " selected by the GMX_CUDA_NB_DEFAULT environment variable.\n"
363 " For best performance upgrade your CUDA toolkit.",
368 kver = eNbnxnCuKLegacy;
373 /* issue not if the non-default kernel is forced by an env. var */
377 "\nNOTE: Legacy non-bonded CUDA kernels were selected by the GMX_CUDA_NB_LEGACY\n"
378 " env. var. Consider using using the default kernels which should be faster!\n");
380 kver = eNbnxnCuKLegacy;
387 void nbnxn_cuda_init(FILE *fplog,
388 nbnxn_cuda_ptr_t *p_cu_nb,
389 gmx_gpu_info_t *gpu_info, int my_gpu_index,
390 gmx_bool bLocalAndNonlocal)
395 bool bStreamSync, bNoStreamSync, bTMPIAtomics, bX86;
399 if (p_cu_nb == NULL) return;
403 snew(nb->nbparam, 1);
404 snew(nb->plist[eintLocal], 1);
405 if (bLocalAndNonlocal)
407 snew(nb->plist[eintNonlocal], 1);
410 nb->bUseTwoStreams = bLocalAndNonlocal;
413 snew(nb->timings, 1);
416 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
417 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
418 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
420 init_plist(nb->plist[eintLocal]);
422 /* local/non-local GPU streams */
423 stat = cudaStreamCreate(&nb->stream[eintLocal]);
424 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
425 if (nb->bUseTwoStreams)
427 init_plist(nb->plist[eintNonlocal]);
428 stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
429 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
432 /* init events for sychronization (timing disabled for performance reasons!) */
433 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
434 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
435 stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
436 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
438 /* set device info, just point it to the right GPU among the detected ones */
439 nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, my_gpu_index)];
441 /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
442 * (which increases with shorter time/step) caused by a known CUDA driver bug.
443 * To work around the issue we'll use an (admittedly fragile) memory polling
444 * waiting to preserve performance. This requires support for atomic
445 * operations and only works on x86/x86_64.
446 * With polling wait event-timing also needs to be disabled.
449 bStreamSync = getenv("GMX_CUDA_STREAMSYNC") != NULL;
450 bNoStreamSync = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
455 bTMPIAtomics = false;
458 #if defined(i386) || defined(__x86_64__)
464 if (bStreamSync && bNoStreamSync)
466 gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
469 if (nb->dev_info->prop.ECCEnabled == 1)
473 nb->bUseStreamSync = true;
476 "NOTE: Using a GPU with ECC enabled, but cudaStreamSynchronize-based waiting is\n"
477 " forced by the GMX_CUDA_STREAMSYNC env. var. Due to a CUDA bug, this \n"
478 " combination causes performance loss.");
479 fprintf(stderr, "\n%s\n", sbuf);
482 fprintf(fplog, "\n%s\n", sbuf);
487 /* can use polling wait only on x86/x86_64 *if* atomics are available */
488 nb->bUseStreamSync = ((bX86 && bTMPIAtomics) == false);
493 "Using a GPU with ECC on; the standard cudaStreamSynchronize waiting, due to a\n"
494 " CUDA bug, causes performance loss when used in combination with ECC.\n"
495 " However, the polling waiting workaround can not be used as it is only\n"
496 " supported on x86/x86_64, but not on the current architecture.");
497 gmx_warning("%s\n", sbuf);
500 fprintf(fplog, "\n%s\n", sbuf);
504 else if (bTMPIAtomics)
509 "NOTE: Using a GPU with ECC enabled; will use polling waiting.\n");
515 "Using a GPU with ECC on; the standard cudaStreamSynchronize waiting, due to a\n"
516 " CUDA bug, causes performance loss when used in combination with ECC.\n"
517 " However, the polling waiting workaround can not be used as atomic\n"
518 " operations are not supported by the current CPU+compiler combination.");
519 gmx_warning("%s\n", sbuf);
522 fprintf(fplog, "\n%s\n", sbuf);
531 nb->bUseStreamSync = false;
534 "NOTE: Using a GPU with no/disabled ECC, but cudaStreamSynchronize-based waiting\n"
535 " is turned off and polling turned on by the GMX_NO_CUDA_STREAMSYNC env. var.");
536 fprintf(stderr, "\n%s\n", sbuf);
539 fprintf(fplog, "\n%s\n", sbuf);
544 /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
545 nb->bUseStreamSync = true;
549 /* CUDA timing disabled as event timers don't work:
550 - with multiple streams = domain-decomposition;
551 - with the polling waiting hack (without cudaStreamSynchronize);
552 - when turned off by GMX_DISABLE_CUDA_TIMING.
554 nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
555 (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
559 init_timers(nb->timers, nb->bUseTwoStreams);
560 init_timings(nb->timings);
563 /* set the kernel type for the current GPU */
564 nb->kernel_ver = pick_nbnxn_kernel_version();
565 /* pick L1 cache configuration */
566 nbnxn_cuda_set_cacheconfig(nb->dev_info);
572 fprintf(debug, "Initialized CUDA data structures.\n");
576 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t cu_nb,
577 const interaction_const_t *ic,
578 const nonbonded_verlet_t *nbv)
580 init_atomdata_first(cu_nb->atdat, nbv->grp[0].nbat->ntype);
581 init_nbparam(cu_nb->nbparam, ic, nbv);
583 /* clear energy and shift force outputs */
584 nbnxn_cuda_clear_e_fshift(cu_nb);
587 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t cu_nb,
588 const nbnxn_pairlist_t *h_plist,
593 bool bDoTime = cu_nb->bDoTime;
594 cudaStream_t stream = cu_nb->stream[iloc];
595 cu_plist_t *d_plist = cu_nb->plist[iloc];
597 if (d_plist->na_c < 0)
599 d_plist->na_c = h_plist->na_ci;
603 if (d_plist->na_c != h_plist->na_ci)
605 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
606 d_plist->na_c, h_plist->na_ci);
613 stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
614 CU_RET_ERR(stat, "cudaEventRecord failed");
617 cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
618 &d_plist->nsci, &d_plist->sci_nalloc,
622 cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
623 &d_plist->ncj4, &d_plist->cj4_nalloc,
627 cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
628 &d_plist->nexcl, &d_plist->excl_nalloc,
634 stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
635 CU_RET_ERR(stat, "cudaEventRecord failed");
638 /* need to prune the pair list during the next step */
639 d_plist->bDoPrune = true;
642 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t cu_nb,
643 const nbnxn_atomdata_t *nbatom)
645 cu_atomdata_t *adat = cu_nb->atdat;
646 cudaStream_t ls = cu_nb->stream[eintLocal];
648 /* only if we have a dynamic box */
649 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
651 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
652 SHIFTS * sizeof(*adat->shift_vec), ls);
653 adat->bShiftVecUploaded = true;
657 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
658 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
661 cu_atomdata_t *adat = cu_nb->atdat;
662 cudaStream_t ls = cu_nb->stream[eintLocal];
664 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
665 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
668 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
669 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
672 cu_atomdata_t *adat = cu_nb->atdat;
673 cudaStream_t ls = cu_nb->stream[eintLocal];
675 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
676 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
677 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
678 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
679 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
680 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
683 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
685 nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
686 /* clear shift force array and energies if the outputs were
687 used in the current step */
688 if (flags & GMX_FORCE_VIRIAL)
690 nbnxn_cuda_clear_e_fshift(cu_nb);
694 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t cu_nb,
695 const nbnxn_atomdata_t *nbat)
700 bool bDoTime = cu_nb->bDoTime;
701 cu_timers_t *timers = cu_nb->timers;
702 cu_atomdata_t *d_atdat = cu_nb->atdat;
703 cudaStream_t ls = cu_nb->stream[eintLocal];
705 natoms = nbat->natoms;
710 /* time async copy */
711 stat = cudaEventRecord(timers->start_atdat, ls);
712 CU_RET_ERR(stat, "cudaEventRecord failed");
715 /* need to reallocate if we have to copy more atoms than the amount of space
716 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
717 if (natoms > d_atdat->nalloc)
719 nalloc = over_alloc_small(natoms);
721 /* free up first if the arrays have already been initialized */
722 if (d_atdat->nalloc != -1)
724 cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
725 cu_free_buffered(d_atdat->xq);
726 cu_free_buffered(d_atdat->atom_types);
729 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
730 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
731 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
732 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
734 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
735 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
737 d_atdat->nalloc = nalloc;
741 d_atdat->natoms = natoms;
742 d_atdat->natoms_local = nbat->natoms_local;
744 /* need to clear GPU f output if realloc happened */
747 nbnxn_cuda_clear_f(cu_nb, nalloc);
750 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
751 natoms*sizeof(*d_atdat->atom_types), ls);
755 stat = cudaEventRecord(timers->stop_atdat, ls);
756 CU_RET_ERR(stat, "cudaEventRecord failed");
760 void nbnxn_cuda_free(FILE *fplog, nbnxn_cuda_ptr_t cu_nb)
763 cu_atomdata_t *atdat;
764 cu_nbparam_t *nbparam;
765 cu_plist_t *plist, *plist_nl;
768 if (cu_nb == NULL) return;
770 atdat = cu_nb->atdat;
771 nbparam = cu_nb->nbparam;
772 plist = cu_nb->plist[eintLocal];
773 plist_nl = cu_nb->plist[eintNonlocal];
774 timers = cu_nb->timers;
776 if (nbparam->eeltype == eelCuEWALD || nbparam->eeltype == eelCuEWALD_TWIN)
778 stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
779 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
780 cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
783 stat = cudaEventDestroy(cu_nb->nonlocal_done);
784 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
785 stat = cudaEventDestroy(cu_nb->misc_ops_done);
786 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
790 stat = cudaEventDestroy(timers->start_atdat);
791 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
792 stat = cudaEventDestroy(timers->stop_atdat);
793 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
795 /* The non-local counters/stream (second in the array) are needed only with DD. */
796 for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
798 stat = cudaEventDestroy(timers->start_nb_k[i]);
799 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
800 stat = cudaEventDestroy(timers->stop_nb_k[i]);
801 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
803 stat = cudaEventDestroy(timers->start_pl_h2d[i]);
804 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
805 stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
806 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
808 stat = cudaStreamDestroy(cu_nb->stream[i]);
809 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
811 stat = cudaEventDestroy(timers->start_nb_h2d[i]);
812 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
813 stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
814 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
816 stat = cudaEventDestroy(timers->start_nb_d2h[i]);
817 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
818 stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
819 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
823 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
824 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
825 cu_free_buffered(nbparam->nbfp);
827 stat = cudaFree(atdat->shift_vec);
828 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
829 stat = cudaFree(atdat->fshift);
830 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
832 stat = cudaFree(atdat->e_lj);
833 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
834 stat = cudaFree(atdat->e_el);
835 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
837 cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
838 cu_free_buffered(atdat->xq);
839 cu_free_buffered(atdat->atom_types, &atdat->ntypes);
841 cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
842 cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
843 cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
844 if (cu_nb->bUseTwoStreams)
846 cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
847 cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
848 cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
853 fprintf(debug, "Cleaned up CUDA data structures.\n");
857 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
860 cudaStream_t stream = cu_nb->stream[iloc];
862 stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
863 CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
866 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
868 return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
871 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
875 init_timings(cu_nb->timings);
879 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
881 return cu_nb != NULL ?
882 gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;