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, const char *buf)
80 /* We should only print to stderr on the master node,
81 * in most cases fplog is only set on the master node, so this works.
83 fprintf(stderr, "\n%s\n", buf);
84 fprintf(fplog, "\n%s\n", buf);
89 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb);
92 /*! Tabulates the Ewald Coulomb force and initializes the size/scale
93 and the table GPU array. If called with an already allocated table,
94 it just re-uploads the table.
96 static void init_ewald_coulomb_force_table(cu_nbparam_t *nbp)
98 float *ftmp, *coul_tab;
103 tabsize = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
104 /* Subtract 2 iso 1 to avoid access out of range due to rounding */
105 tabscale = (tabsize - 2) / sqrt(nbp->rcoulomb_sq);
107 pmalloc((void**)&ftmp, tabsize*sizeof(*ftmp));
109 table_spline3_fill_ewald_lr(ftmp, NULL, NULL, tabsize,
110 1/tabscale, nbp->ewald_beta);
112 /* If the table pointer == NULL the table is generated the first time =>
113 the array pointer will be saved to nbparam and the texture is bound.
115 coul_tab = nbp->coulomb_tab;
116 if (coul_tab == NULL)
118 stat = cudaMalloc((void **)&coul_tab, tabsize*sizeof(*coul_tab));
119 CU_RET_ERR(stat, "cudaMalloc failed on coul_tab");
121 nbp->coulomb_tab = coul_tab;
123 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
124 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_coulomb_tab_texref(),
125 coul_tab, &cd, tabsize*sizeof(*coul_tab));
126 CU_RET_ERR(stat, "cudaBindTexture on coul_tab failed");
129 cu_copy_H2D(coul_tab, ftmp, tabsize*sizeof(*coul_tab));
131 nbp->coulomb_tab_size = tabsize;
132 nbp->coulomb_tab_scale = tabscale;
138 /*! Initializes the atomdata structure first time, it only gets filled at
140 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
145 stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
146 CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
147 ad->bShiftVecUploaded = false;
149 stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
150 CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
152 stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
153 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
154 stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
155 CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
157 /* initialize to NULL poiters to data that is not allocated here and will
158 need reallocation in nbnxn_cuda_init_atomdata */
162 /* size -1 indicates that the respective array hasn't been initialized yet */
167 /*! Initializes the nonbonded parameter data structure. */
168 static void init_nbparam(cu_nbparam_t *nbp,
169 const interaction_const_t *ic,
170 const nonbonded_verlet_t *nbv)
175 ntypes = nbv->grp[0].nbat->ntype;
177 nbp->ewald_beta = ic->ewaldcoeff;
178 nbp->sh_ewald = ic->sh_ewald;
179 nbp->epsfac = ic->epsfac;
180 nbp->two_k_rf = 2.0 * ic->k_rf;
181 nbp->c_rf = ic->c_rf;
182 nbp->rvdw_sq = ic->rvdw * ic->rvdw;
183 nbp->rcoulomb_sq= ic->rcoulomb * ic->rcoulomb;
184 nbp->rlist_sq = ic->rlist * ic->rlist;
185 nbp->sh_invrc6 = ic->sh_invrc6;
187 if (ic->eeltype == eelCUT)
189 nbp->eeltype = eelCuCUT;
191 else if (EEL_RF(ic->eeltype))
193 nbp->eeltype = eelCuRF;
195 else if ((EEL_PME(ic->eeltype) || ic->eeltype==eelEWALD))
197 /* Initially rcoulomb == rvdw, so it's surely not twin cut-off, unless
198 forced by the env. var. (used only for benchmarking). */
199 if (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL)
201 nbp->eeltype = eelCuEWALD;
205 nbp->eeltype = eelCuEWALD_TWIN;
210 /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
211 gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
214 /* generate table for PME */
215 if (nbp->eeltype == eelCuEWALD)
217 nbp->coulomb_tab = NULL;
218 init_ewald_coulomb_force_table(nbp);
221 nnbfp = 2*ntypes*ntypes;
222 stat = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
223 CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
224 cu_copy_H2D(nbp->nbfp, nbv->grp[0].nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
226 cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
227 stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
228 nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
229 CU_RET_ERR(stat, "cudaBindTexture on nbfp failed");
232 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
233 * electrostatic type switching to twin cut-off (or back) if needed. */
234 void nbnxn_cuda_pme_loadbal_update_param(nbnxn_cuda_ptr_t cu_nb,
235 const interaction_const_t *ic)
237 cu_nbparam_t *nbp = cu_nb->nbparam;
239 nbp->rlist_sq = ic->rlist * ic->rlist;
240 nbp->rcoulomb_sq = ic->rcoulomb * ic->rcoulomb;
241 nbp->ewald_beta = ic->ewaldcoeff;
243 /* When switching to/from twin cut-off, the electrostatics type needs updating.
244 (The env. var. that forces twin cut-off is for benchmarking only!) */
245 if (ic->rcoulomb == ic->rvdw &&
246 getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL)
248 nbp->eeltype = eelCuEWALD;
252 nbp->eeltype = eelCuEWALD_TWIN;
255 init_ewald_coulomb_force_table(cu_nb->nbparam);
258 /*! Initializes the pair list data structure. */
259 static void init_plist(cu_plist_t *pl)
261 /* initialize to NULL pointers to data that is not allocated here and will
262 need reallocation in nbnxn_cuda_init_pairlist */
267 /* size -1 indicates that the respective array hasn't been initialized yet */
274 pl->excl_nalloc = -1;
275 pl->bDoPrune = false;
278 /*! Initializes the timer data structure. */
279 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
282 int eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync: cudaEventDefault );
284 stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
285 CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
286 stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
287 CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
289 /* The non-local counters/stream (second in the array) are needed only with DD. */
290 for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
292 stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
293 CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
294 stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
295 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
298 stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
299 CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
300 stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
301 CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
303 stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
304 CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
305 stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
306 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
308 stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
309 CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
310 stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
311 CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
315 /*! Initializes the timings data structure. */
316 static void init_timings(wallclock_gpu_t *t)
325 for (i = 0; i < 2; i++)
327 for(j = 0; j < 2; j++)
329 t->ktime[i][j].t = 0.0;
330 t->ktime[i][j].c = 0;
335 /* Decide which kernel version to use (default or legacy) based on:
337 * - non-bonded kernel selector environment variables
338 * - GPU SM version TODO ???
340 static int pick_nbnxn_kernel_version()
342 bool bLegacyKernel, bDefaultKernel, bCUDA40, bCUDA32;
346 /* legacy kernel (former k2), kept for now for backward compatibility,
347 faster than the default with CUDA 3.2/4.0 (TODO: on Kepler?). */
348 bLegacyKernel = (getenv("GMX_CUDA_NB_LEGACY") != NULL);
349 /* default kernel (former k3). */
350 bDefaultKernel = (getenv("GMX_CUDA_NB_DEFAULT") != NULL);
352 if ((unsigned)(bLegacyKernel + bDefaultKernel) > 1)
354 gmx_fatal(FARGS, "Multiple CUDA non-bonded kernels requested; to manually pick a kernel set only one \n"
355 "of the following environment variables: \n"
356 "GMX_CUDA_NB_DEFAULT, GMX_CUDA_NB_LEGACY");
359 bCUDA32 = bCUDA40 = false;
360 #if CUDA_VERSION == 3200
362 sprintf(sbuf, "3.2");
363 #elif CUDA_VERSION == 4000
365 sprintf(sbuf, "4.0");
368 /* default is default ;) */
369 kver = eNbnxnCuKDefault;
371 if (bCUDA32 || bCUDA40)
373 /* use legacy kernel unless something else is forced by an env. var */
377 "\nNOTE: CUDA %s compilation detected; with this compiler version the legacy\n"
378 " non-bonded kernels perform best. However, the default kernels were\n"
379 " selected by the GMX_CUDA_NB_DEFAULT environment variable.\n"
380 " For best performance upgrade your CUDA toolkit.",
385 kver = eNbnxnCuKLegacy;
390 /* issue not if the non-default kernel is forced by an env. var */
394 "\nNOTE: Legacy non-bonded CUDA kernels were selected by the GMX_CUDA_NB_LEGACY\n"
395 " env. var. Consider using using the default kernels which should be faster!\n");
397 kver = eNbnxnCuKLegacy;
404 void nbnxn_cuda_init(FILE *fplog,
405 nbnxn_cuda_ptr_t *p_cu_nb,
406 gmx_gpu_info_t *gpu_info, int my_gpu_index,
407 gmx_bool bLocalAndNonlocal)
412 bool bStreamSync, bNoStreamSync, bTMPIAtomics, bX86, bOldDriver;
417 if (p_cu_nb == NULL) return;
421 snew(nb->nbparam, 1);
422 snew(nb->plist[eintLocal], 1);
423 if (bLocalAndNonlocal)
425 snew(nb->plist[eintNonlocal], 1);
428 nb->bUseTwoStreams = bLocalAndNonlocal;
431 snew(nb->timings, 1);
434 pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
435 pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
436 pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
438 init_plist(nb->plist[eintLocal]);
440 /* local/non-local GPU streams */
441 stat = cudaStreamCreate(&nb->stream[eintLocal]);
442 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
443 if (nb->bUseTwoStreams)
445 init_plist(nb->plist[eintNonlocal]);
446 stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
447 CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
450 /* init events for sychronization (timing disabled for performance reasons!) */
451 stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
452 CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
453 stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
454 CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
456 /* set device info, just point it to the right GPU among the detected ones */
457 nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, my_gpu_index)];
459 /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
460 * (which increases with shorter time/step) caused by a known CUDA driver bug.
461 * To work around the issue we'll use an (admittedly fragile) memory polling
462 * waiting to preserve performance. This requires support for atomic
463 * operations and only works on x86/x86_64.
464 * With polling wait event-timing also needs to be disabled.
466 * The overhead is greatly reduced in API v5.0 drivers and the improvement
467 $ is independent of runtime version. Hence, with API v5.0 drivers and later
468 * we won't switch to polling.
470 * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
471 * ranks so we will also disable it in that case.
474 bStreamSync = getenv("GMX_CUDA_STREAMSYNC") != NULL;
475 bNoStreamSync = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
480 bTMPIAtomics = false;
483 #if defined(i386) || defined(__x86_64__)
489 if (bStreamSync && bNoStreamSync)
491 gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
494 stat = cudaDriverGetVersion(&cuda_drv_ver);
495 CU_RET_ERR(stat, "cudaDriverGetVersion failed");
496 bOldDriver = (cuda_drv_ver < 5000);
498 if (nb->dev_info->prop.ECCEnabled == 1)
502 nb->bUseStreamSync = true;
504 /* only warn if polling should be used */
505 if (bOldDriver && !gpu_info->bDevShare)
508 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
509 " cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
514 /* Can/should turn of cudaStreamSynchronize wait only if
515 * - we're on x86/x86_64
516 * - atomics are available
517 * - GPUs are not being shared
518 * - and driver is old. */
520 (bX86 && bTMPIAtomics && !gpu_info->bDevShare && bOldDriver) ?
523 if (nb->bUseStreamSync)
526 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
527 " cause performance loss. Switching to the alternative polling GPU waiting.\n"
528 " If you encounter issues, switch back to standard GPU waiting by setting\n"
529 " the GMX_CUDA_STREAMSYNC environment variable.\n");
533 /* Tell the user that the ECC+old driver combination can be bad */
535 "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0. A bug in this\n"
536 " driver can cause performance loss.\n"
537 " However, the polling waiting workaround can not be used because\n%s\n"
538 " Consider updating the driver or turning ECC off.",
539 (!bX86 || !bTMPIAtomics) ?
540 " atomic operations are not supported by the platform/CPU+compiler." :
541 " GPU(s) are being oversubscribed.");
542 md_print_warn(fplog, sbuf);
550 nb->bUseStreamSync = false;
553 "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
557 /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
558 nb->bUseStreamSync = true;
562 /* CUDA timing disabled as event timers don't work:
563 - with multiple streams = domain-decomposition;
564 - with the polling waiting hack (without cudaStreamSynchronize);
565 - when turned off by GMX_DISABLE_CUDA_TIMING.
567 nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
568 (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
572 init_timers(nb->timers, nb->bUseTwoStreams);
573 init_timings(nb->timings);
576 /* set the kernel type for the current GPU */
577 nb->kernel_ver = pick_nbnxn_kernel_version();
578 /* pick L1 cache configuration */
579 nbnxn_cuda_set_cacheconfig(nb->dev_info);
585 fprintf(debug, "Initialized CUDA data structures.\n");
589 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t cu_nb,
590 const interaction_const_t *ic,
591 const nonbonded_verlet_t *nbv)
593 init_atomdata_first(cu_nb->atdat, nbv->grp[0].nbat->ntype);
594 init_nbparam(cu_nb->nbparam, ic, nbv);
596 /* clear energy and shift force outputs */
597 nbnxn_cuda_clear_e_fshift(cu_nb);
600 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t cu_nb,
601 const nbnxn_pairlist_t *h_plist,
606 bool bDoTime = cu_nb->bDoTime;
607 cudaStream_t stream = cu_nb->stream[iloc];
608 cu_plist_t *d_plist = cu_nb->plist[iloc];
610 if (d_plist->na_c < 0)
612 d_plist->na_c = h_plist->na_ci;
616 if (d_plist->na_c != h_plist->na_ci)
618 sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
619 d_plist->na_c, h_plist->na_ci);
626 stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
627 CU_RET_ERR(stat, "cudaEventRecord failed");
630 cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
631 &d_plist->nsci, &d_plist->sci_nalloc,
635 cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
636 &d_plist->ncj4, &d_plist->cj4_nalloc,
640 cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
641 &d_plist->nexcl, &d_plist->excl_nalloc,
647 stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
648 CU_RET_ERR(stat, "cudaEventRecord failed");
651 /* need to prune the pair list during the next step */
652 d_plist->bDoPrune = true;
655 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t cu_nb,
656 const nbnxn_atomdata_t *nbatom)
658 cu_atomdata_t *adat = cu_nb->atdat;
659 cudaStream_t ls = cu_nb->stream[eintLocal];
661 /* only if we have a dynamic box */
662 if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
664 cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
665 SHIFTS * sizeof(*adat->shift_vec), ls);
666 adat->bShiftVecUploaded = true;
670 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
671 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
674 cu_atomdata_t *adat = cu_nb->atdat;
675 cudaStream_t ls = cu_nb->stream[eintLocal];
677 stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
678 CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
681 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
682 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
685 cu_atomdata_t *adat = cu_nb->atdat;
686 cudaStream_t ls = cu_nb->stream[eintLocal];
688 stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
689 CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
690 stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
691 CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
692 stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
693 CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
696 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
698 nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
699 /* clear shift force array and energies if the outputs were
700 used in the current step */
701 if (flags & GMX_FORCE_VIRIAL)
703 nbnxn_cuda_clear_e_fshift(cu_nb);
707 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t cu_nb,
708 const nbnxn_atomdata_t *nbat)
713 bool bDoTime = cu_nb->bDoTime;
714 cu_timers_t *timers = cu_nb->timers;
715 cu_atomdata_t *d_atdat = cu_nb->atdat;
716 cudaStream_t ls = cu_nb->stream[eintLocal];
718 natoms = nbat->natoms;
723 /* time async copy */
724 stat = cudaEventRecord(timers->start_atdat, ls);
725 CU_RET_ERR(stat, "cudaEventRecord failed");
728 /* need to reallocate if we have to copy more atoms than the amount of space
729 available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
730 if (natoms > d_atdat->nalloc)
732 nalloc = over_alloc_small(natoms);
734 /* free up first if the arrays have already been initialized */
735 if (d_atdat->nalloc != -1)
737 cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
738 cu_free_buffered(d_atdat->xq);
739 cu_free_buffered(d_atdat->atom_types);
742 stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
743 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
744 stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
745 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
747 stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
748 CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
750 d_atdat->nalloc = nalloc;
754 d_atdat->natoms = natoms;
755 d_atdat->natoms_local = nbat->natoms_local;
757 /* need to clear GPU f output if realloc happened */
760 nbnxn_cuda_clear_f(cu_nb, nalloc);
763 cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
764 natoms*sizeof(*d_atdat->atom_types), ls);
768 stat = cudaEventRecord(timers->stop_atdat, ls);
769 CU_RET_ERR(stat, "cudaEventRecord failed");
773 void nbnxn_cuda_free(FILE *fplog, nbnxn_cuda_ptr_t cu_nb)
776 cu_atomdata_t *atdat;
777 cu_nbparam_t *nbparam;
778 cu_plist_t *plist, *plist_nl;
781 if (cu_nb == NULL) return;
783 atdat = cu_nb->atdat;
784 nbparam = cu_nb->nbparam;
785 plist = cu_nb->plist[eintLocal];
786 plist_nl = cu_nb->plist[eintNonlocal];
787 timers = cu_nb->timers;
789 if (nbparam->eeltype == eelCuEWALD || nbparam->eeltype == eelCuEWALD_TWIN)
791 stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
792 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
793 cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
796 stat = cudaEventDestroy(cu_nb->nonlocal_done);
797 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
798 stat = cudaEventDestroy(cu_nb->misc_ops_done);
799 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
803 stat = cudaEventDestroy(timers->start_atdat);
804 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
805 stat = cudaEventDestroy(timers->stop_atdat);
806 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
808 /* The non-local counters/stream (second in the array) are needed only with DD. */
809 for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
811 stat = cudaEventDestroy(timers->start_nb_k[i]);
812 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
813 stat = cudaEventDestroy(timers->stop_nb_k[i]);
814 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
816 stat = cudaEventDestroy(timers->start_pl_h2d[i]);
817 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
818 stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
819 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
821 stat = cudaStreamDestroy(cu_nb->stream[i]);
822 CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
824 stat = cudaEventDestroy(timers->start_nb_h2d[i]);
825 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
826 stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
827 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
829 stat = cudaEventDestroy(timers->start_nb_d2h[i]);
830 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
831 stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
832 CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
836 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
837 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab failed");
838 cu_free_buffered(nbparam->nbfp);
840 stat = cudaFree(atdat->shift_vec);
841 CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
842 stat = cudaFree(atdat->fshift);
843 CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
845 stat = cudaFree(atdat->e_lj);
846 CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
847 stat = cudaFree(atdat->e_el);
848 CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
850 cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
851 cu_free_buffered(atdat->xq);
852 cu_free_buffered(atdat->atom_types, &atdat->ntypes);
854 cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
855 cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
856 cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
857 if (cu_nb->bUseTwoStreams)
859 cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
860 cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
861 cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
867 if (cu_nb->bUseTwoStreams)
872 sfree(cu_nb->timings);
877 fprintf(debug, "Cleaned up CUDA data structures.\n");
881 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
884 cudaStream_t stream = cu_nb->stream[iloc];
886 stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
887 CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
890 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
892 return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
895 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
899 init_timings(cu_nb->timings);
903 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
905 return cu_nb != NULL ?
906 gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;