Fix malformed CUDA version macro check
[alexxy/gromacs.git] / src / mdlib / nbnxn_cuda / nbnxn_cuda_data_mgmt.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2012, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdlib.h>
43 #include <stdio.h>
44 #include <assert.h>
45
46 #include <cuda.h>
47
48 #include "gmx_fatal.h"
49 #include "smalloc.h"
50 #include "tables.h"
51 #include "typedefs.h"
52 #include "types/nb_verlet.h"
53 #include "types/interaction_const.h"
54 #include "types/force_flags.h"
55 #include "../nbnxn_consts.h"
56 #include "gmx_detect_hardware.h"
57
58 #include "nbnxn_cuda_types.h"
59 #include "../../gmxlib/cuda_tools/cudautils.cuh"
60 #include "nbnxn_cuda_data_mgmt.h"
61 #include "pmalloc_cuda.h"
62 #include "gpu_utils.h"
63
64 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
65
66 /* This is a heuristically determined parameter for the Fermi architecture for
67  * the minimum size of ci lists by multiplying this constant with the # of
68  * multiprocessors on the current device.
69  */
70 static unsigned int gpu_min_ci_balanced_factor = 40;
71
72 /* Functions from nbnxn_cuda.cu */
73 extern void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo);
74 extern const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_nbfp_texref();
75 extern const struct texture<float, 1, cudaReadModeElementType>& nbnxn_cuda_get_coulomb_tab_texref();
76
77 /* We should actually be using md_print_warn in md_logging.c,
78  * but we can't include mpi.h in CUDA code.
79  */
80 static void md_print_warn(FILE       *fplog,
81                           const char *fmt, ...)
82 {
83     va_list ap;
84
85     if (fplog != NULL)
86     {
87         /* We should only print to stderr on the master node,
88          * in most cases fplog is only set on the master node, so this works.
89          */
90         va_start(ap, fmt);
91         fprintf(stderr, "\n");
92         vfprintf(stderr, fmt, ap);
93         fprintf(stderr, "\n");
94         va_end(ap);
95
96         va_start(ap, fmt);
97         fprintf(fplog, "\n");
98         vfprintf(fplog, fmt, ap);
99         fprintf(fplog, "\n");
100         va_end(ap);
101     }
102 }
103
104
105 /* Fw. decl. */
106 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb);
107
108
109 /*! Tabulates the Ewald Coulomb force and initializes the size/scale
110     and the table GPU array. If called with an already allocated table,
111     it just re-uploads the table.
112  */
113 static void init_ewald_coulomb_force_table(cu_nbparam_t          *nbp,
114                                            const cuda_dev_info_t *dev_info)
115 {
116     float       *ftmp, *coul_tab;
117     int         tabsize;
118     double      tabscale;
119     cudaError_t stat;
120
121     tabsize     = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
122     /* Subtract 2 iso 1 to avoid access out of range due to rounding */
123     tabscale    = (tabsize - 2) / sqrt(nbp->rcoulomb_sq);
124
125     pmalloc((void**)&ftmp, tabsize*sizeof(*ftmp));
126
127     table_spline3_fill_ewald_lr(ftmp, NULL, NULL, tabsize,
128                                 1/tabscale, nbp->ewald_beta);
129
130     /* If the table pointer == NULL the table is generated the first time =>
131        the array pointer will be saved to nbparam and the texture is bound.
132      */
133     coul_tab = nbp->coulomb_tab;
134     if (coul_tab == NULL)
135     {
136         stat = cudaMalloc((void **)&coul_tab, tabsize*sizeof(*coul_tab));
137         CU_RET_ERR(stat, "cudaMalloc failed on coul_tab");
138
139         nbp->coulomb_tab = coul_tab;
140
141 #ifdef TEXOBJ_SUPPORTED
142         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
143         if (dev_info->prop.major >= 3)
144         {
145             cudaResourceDesc rd;
146             memset(&rd, 0, sizeof(rd));
147             rd.resType                  = cudaResourceTypeLinear;
148             rd.res.linear.devPtr        = nbp->coulomb_tab;
149             rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
150             rd.res.linear.desc.x        = 32;
151             rd.res.linear.sizeInBytes   = tabsize*sizeof(*coul_tab);
152
153             cudaTextureDesc td;
154             memset(&td, 0, sizeof(td));
155             td.readMode                 = cudaReadModeElementType;
156             stat = cudaCreateTextureObject(&nbp->coulomb_tab_texobj, &rd, &td, NULL);
157             CU_RET_ERR(stat, "cudaCreateTextureObject on coulomb_tab_texobj failed");
158         }
159         else
160 #endif
161         {
162             cudaChannelFormatDesc cd   = cudaCreateChannelDesc<float>();
163             stat = cudaBindTexture(NULL, &nbnxn_cuda_get_coulomb_tab_texref(),
164                                    coul_tab, &cd, tabsize*sizeof(*coul_tab));
165             CU_RET_ERR(stat, "cudaBindTexture on coulomb_tab_texref failed");
166         }
167     }
168
169     cu_copy_H2D(coul_tab, ftmp, tabsize*sizeof(*coul_tab));
170
171     nbp->coulomb_tab_size     = tabsize;
172     nbp->coulomb_tab_scale    = tabscale;
173
174     pfree(ftmp);
175 }
176
177
178 /*! Initializes the atomdata structure first time, it only gets filled at
179     pair-search. */
180 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
181 {
182     cudaError_t stat;
183
184     ad->ntypes  = ntypes;
185     stat = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
186     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
187     ad->bShiftVecUploaded = false;
188
189     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
190     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
191
192     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
193     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
194     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
195     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
196
197     /* initialize to NULL poiters to data that is not allocated here and will
198        need reallocation in nbnxn_cuda_init_atomdata */
199     ad->xq = NULL;
200     ad->f  = NULL;
201
202     /* size -1 indicates that the respective array hasn't been initialized yet */
203     ad->natoms = -1;
204     ad->nalloc = -1;
205 }
206
207 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
208     earlier GPUs, single or twin cut-off. */
209 static int pick_ewald_kernel_type(bool                   bTwinCut,
210                                   const cuda_dev_info_t *dev_info)
211 {
212     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
213     int  kernel_type;
214
215     /* Benchmarking/development environment variables to force the use of
216        analytical or tabulated Ewald kernel. */
217     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
218     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
219
220     if (bForceAnalyticalEwald && bForceTabulatedEwald)
221     {
222         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
223                    "requested through environment variables.");
224     }
225
226     /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
227     if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
228     {
229         bUseAnalyticalEwald = true;
230
231         if (debug)
232         {
233             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
234         }
235     }
236     else
237     {
238         bUseAnalyticalEwald = false;
239
240         if (debug)
241         {
242             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
243         }
244     }
245
246     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
247        forces it (use it for debugging/benchmarking only). */
248     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
249     {
250         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
251     }
252     else
253     {
254         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
255     }
256
257     return kernel_type;
258 }
259
260
261 /*! Initializes the nonbonded parameter data structure. */
262 static void init_nbparam(cu_nbparam_t *nbp,
263                          const interaction_const_t *ic,
264                          const nbnxn_atomdata_t *nbat,
265                          const cuda_dev_info_t *dev_info)
266 {
267     cudaError_t stat;
268     int         ntypes, nnbfp;
269
270     ntypes  = nbat->ntype;
271
272     nbp->ewald_beta = ic->ewaldcoeff;
273     nbp->sh_ewald   = ic->sh_ewald;
274     nbp->epsfac     = ic->epsfac;
275     nbp->two_k_rf   = 2.0 * ic->k_rf;
276     nbp->c_rf       = ic->c_rf;
277     nbp->rvdw_sq    = ic->rvdw * ic->rvdw;
278     nbp->rcoulomb_sq= ic->rcoulomb * ic->rcoulomb;
279     nbp->rlist_sq   = ic->rlist * ic->rlist;
280     nbp->sh_invrc6  = ic->sh_invrc6;
281
282     if (ic->eeltype == eelCUT)
283     {
284         nbp->eeltype = eelCuCUT;
285     }
286     else if (EEL_RF(ic->eeltype))
287     {
288         nbp->eeltype = eelCuRF;
289     }
290     else if ((EEL_PME(ic->eeltype) || ic->eeltype==eelEWALD))
291     {
292         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
293         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
294     }
295     else
296     {
297         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
298         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
299     }
300
301     /* generate table for PME */
302     nbp->coulomb_tab = NULL;
303     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
304     {
305         init_ewald_coulomb_force_table(nbp, dev_info);
306     }
307
308     nnbfp = 2*ntypes*ntypes;
309     stat = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
310     CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
311     cu_copy_H2D(nbp->nbfp, nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
312
313 #ifdef TEXOBJ_SUPPORTED
314         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
315         if (dev_info->prop.major >= 3)
316         {
317             cudaResourceDesc rd;
318             memset(&rd, 0, sizeof(rd));
319             rd.resType                  = cudaResourceTypeLinear;
320             rd.res.linear.devPtr        = nbp->nbfp;
321             rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
322             rd.res.linear.desc.x        = 32;
323             rd.res.linear.sizeInBytes   = nnbfp*sizeof(*nbp->nbfp);
324
325             cudaTextureDesc td;
326             memset(&td, 0, sizeof(td));
327             td.readMode                 = cudaReadModeElementType;
328             stat = cudaCreateTextureObject(&nbp->nbfp_texobj, &rd, &td, NULL);
329             CU_RET_ERR(stat, "cudaCreateTextureObject on nbfp_texobj failed");
330         }
331         else
332 #endif
333         {
334             cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
335             stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
336                                    nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
337             CU_RET_ERR(stat, "cudaBindTexture on nbfp_texref failed");
338         }
339 }
340
341 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
342  *  electrostatic type switching to twin cut-off (or back) if needed. */
343 void nbnxn_cuda_pme_loadbal_update_param(nbnxn_cuda_ptr_t cu_nb,
344                                          const interaction_const_t *ic)
345 {
346     cu_nbparam_t *nbp = cu_nb->nbparam;
347
348     nbp->rlist_sq       = ic->rlist * ic->rlist;
349     nbp->rcoulomb_sq    = ic->rcoulomb * ic->rcoulomb;
350     nbp->ewald_beta     = ic->ewaldcoeff;
351
352     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
353                                                  cu_nb->dev_info);
354
355     init_ewald_coulomb_force_table(cu_nb->nbparam, cu_nb->dev_info);
356 }
357
358 /*! Initializes the pair list data structure. */
359 static void init_plist(cu_plist_t *pl)
360 {
361     /* initialize to NULL pointers to data that is not allocated here and will
362        need reallocation in nbnxn_cuda_init_pairlist */
363     pl->sci     = NULL;
364     pl->cj4     = NULL;
365     pl->excl    = NULL;
366
367     /* size -1 indicates that the respective array hasn't been initialized yet */
368     pl->na_c        = -1;
369     pl->nsci        = -1;
370     pl->sci_nalloc  = -1;
371     pl->ncj4        = -1;
372     pl->cj4_nalloc  = -1;
373     pl->nexcl       = -1;
374     pl->excl_nalloc = -1;
375     pl->bDoPrune    = false;
376 }
377
378 /*! Initializes the timer data structure. */
379 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
380 {
381     cudaError_t stat;
382     int eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync: cudaEventDefault );
383
384     stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
385     CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
386     stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
387     CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
388
389     /* The non-local counters/stream (second in the array) are needed only with DD. */
390     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
391     {
392         stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
393         CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
394         stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
395         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
396
397
398         stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
399         CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
400         stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
401         CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
402
403         stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
404         CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
405         stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
406         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
407
408         stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
409         CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
410         stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
411         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
412     }
413 }
414
415 /*! Initializes the timings data structure. */
416 static void init_timings(wallclock_gpu_t *t)
417 {
418     int i, j;
419
420     t->nb_h2d_t = 0.0;
421     t->nb_d2h_t = 0.0;
422     t->nb_c    = 0;
423     t->pl_h2d_t = 0.0;
424     t->pl_h2d_c = 0;
425     for (i = 0; i < 2; i++)
426     {
427         for(j = 0; j < 2; j++)
428         {
429             t->ktime[i][j].t = 0.0;
430             t->ktime[i][j].c = 0;
431         }
432     }
433 }
434
435 /* Decide which kernel version to use (default or legacy) based on:
436  *  - CUDA version used for compilation
437  *  - non-bonded kernel selector environment variables
438  *  - GPU architecture version
439  */
440 static int pick_nbnxn_kernel_version(FILE            *fplog,
441                                      cuda_dev_info_t *devinfo)
442 {
443     bool bForceLegacyKernel, bForceDefaultKernel, bCUDA40, bCUDA32;
444     char sbuf[STRLEN];
445     int  kver;
446
447     /* Legacy kernel (former k2), kept for backward compatibility as it is
448        faster than the default with CUDA 3.2/4.0 on Fermi (not on Kepler). */
449     bForceLegacyKernel  = (getenv("GMX_CUDA_NB_LEGACY") != NULL);
450     /* default kernel (former k3). */
451     bForceDefaultKernel = (getenv("GMX_CUDA_NB_DEFAULT") != NULL);
452
453     if ((unsigned)(bForceLegacyKernel + bForceDefaultKernel) > 1)
454     {
455         gmx_fatal(FARGS, "Multiple CUDA non-bonded kernels requested; to manually pick a kernel set only one \n"
456                   "of the following environment variables: \n"
457                   "GMX_CUDA_NB_DEFAULT, GMX_CUDA_NB_LEGACY");
458     }
459
460     bCUDA32 = bCUDA40 = false;
461 #if CUDA_VERSION == 3020
462     bCUDA32 = true;
463     sprintf(sbuf, "3.2");
464 #elif CUDA_VERSION == 4000
465     bCUDA40 = true;
466     sprintf(sbuf, "4.0");
467 #endif
468
469     /* default is default ;) */
470     kver = eNbnxnCuKDefault;
471
472     /* Consider switching to legacy kernels only on Fermi */
473     if (devinfo->prop.major < 3 && (bCUDA32 || bCUDA40))
474     {
475         /* use legacy kernel unless something else is forced by an env. var */
476         if (bForceDefaultKernel)
477         {
478             md_print_warn(fplog,
479                           "NOTE: CUDA %s compilation detected; with this compiler version the legacy\n"
480                           "      non-bonded kernels perform best. However, the default kernels were\n"
481                           "      selected by the GMX_CUDA_NB_DEFAULT environment variable.\n"
482                           "      For best performance upgrade your CUDA toolkit.\n",
483                           sbuf);
484         }
485         else
486         {
487             kver = eNbnxnCuKLegacy;
488         }
489     }
490     else
491     {
492         /* issue note if the non-default kernel is forced by an env. var */
493         if (bForceLegacyKernel)
494         {
495             md_print_warn(fplog,
496                     "NOTE: Legacy non-bonded CUDA kernels selected by the GMX_CUDA_NB_LEGACY\n"
497                     "      env. var. Consider using using the default kernels which should be faster!\n");
498
499             kver = eNbnxnCuKLegacy;
500         }
501     }
502
503     return kver;
504 }
505
506 void nbnxn_cuda_init(FILE *fplog,
507                      nbnxn_cuda_ptr_t *p_cu_nb,
508                      const gmx_gpu_info_t *gpu_info,
509                      const gmx_gpu_opt_t *gpu_opt,
510                      int my_gpu_index,
511                      gmx_bool bLocalAndNonlocal)
512 {
513     cudaError_t stat;
514     nbnxn_cuda_ptr_t  nb;
515     char sbuf[STRLEN];
516     bool bStreamSync, bNoStreamSync, bTMPIAtomics, bX86, bOldDriver;
517     int cuda_drv_ver;
518
519     assert(gpu_info);
520
521     if (p_cu_nb == NULL) return;
522
523     snew(nb, 1);
524     snew(nb->atdat, 1);
525     snew(nb->nbparam, 1);
526     snew(nb->plist[eintLocal], 1);
527     if (bLocalAndNonlocal)
528     {
529         snew(nb->plist[eintNonlocal], 1);
530     }
531
532     nb->bUseTwoStreams = bLocalAndNonlocal;
533
534     snew(nb->timers, 1);
535     snew(nb->timings, 1);
536
537     /* init nbst */
538     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
539     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
540     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
541
542     init_plist(nb->plist[eintLocal]);
543
544     /* set device info, just point it to the right GPU among the detected ones */
545     nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, gpu_opt, my_gpu_index)];
546
547     /* local/non-local GPU streams */
548     stat = cudaStreamCreate(&nb->stream[eintLocal]);
549     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
550     if (nb->bUseTwoStreams)
551     {
552         init_plist(nb->plist[eintNonlocal]);
553
554         /* CUDA stream priority available in the CUDA RT 5.5 API.
555          * Note that the device we're running on does not have to support
556          * priorities, because we are querying the priority range which in this
557          * case will be a single value.
558          */
559 #if CUDA_VERSION >= 5050
560         {
561             int highest_priority;
562             stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
563             CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
564
565             stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
566                                                 cudaStreamDefault,
567                                                 highest_priority);
568             CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
569         }
570 #else
571         stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
572         CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
573 #endif
574     }
575
576     /* init events for sychronization (timing disabled for performance reasons!) */
577     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
578     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
579     stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
580     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
581
582     /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
583      * (which increases with shorter time/step) caused by a known CUDA driver bug.
584      * To work around the issue we'll use an (admittedly fragile) memory polling
585      * waiting to preserve performance. This requires support for atomic
586      * operations and only works on x86/x86_64.
587      * With polling wait event-timing also needs to be disabled.
588      *
589      * The overhead is greatly reduced in API v5.0 drivers and the improvement
590      $ is independent of runtime version. Hence, with API v5.0 drivers and later
591      * we won't switch to polling.
592      *
593      * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
594      * ranks so we will also disable it in that case.
595      */
596
597     bStreamSync    = getenv("GMX_CUDA_STREAMSYNC") != NULL;
598     bNoStreamSync  = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
599
600 #ifdef TMPI_ATOMICS
601     bTMPIAtomics = true;
602 #else
603     bTMPIAtomics = false;
604 #endif
605
606 #ifdef GMX_TARGET_X86
607     bX86 = true;
608 #else
609     bX86 = false;
610 #endif
611
612     if (bStreamSync && bNoStreamSync)
613     {
614         gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
615     }
616
617     stat = cudaDriverGetVersion(&cuda_drv_ver);
618     CU_RET_ERR(stat, "cudaDriverGetVersion failed");
619
620     bOldDriver = (cuda_drv_ver < 5000);
621
622     if ((nb->dev_info->prop.ECCEnabled == 1) && bOldDriver)
623     {
624         /* Polling wait should be used instead of cudaStreamSynchronize only if:
625          *   - ECC is ON & driver is old (checked above),
626          *   - we're on x86/x86_64,
627          *   - atomics are available, and
628          *   - GPUs are not being shared.
629          */
630         bool bShouldUsePollSync = (bX86 && bTMPIAtomics &&
631                                    (gmx_count_gpu_dev_shared(gpu_opt) < 1));
632
633         if (bStreamSync)
634         {
635             nb->bUseStreamSync = true;
636
637             /* only warn if polling should be used */
638             if (bShouldUsePollSync)
639             {
640                 md_print_warn(fplog,
641                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
642                               "      cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
643             }
644         }
645         else
646         {
647             nb->bUseStreamSync = !bShouldUsePollSync;
648
649             if (bShouldUsePollSync)
650             {
651                 md_print_warn(fplog,
652                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
653                               "      cause performance loss. Switching to the alternative polling GPU wait.\n"
654                               "      If you encounter issues, switch back to standard GPU waiting by setting\n"
655                               "      the GMX_CUDA_STREAMSYNC environment variable.\n");
656             }
657             else
658             {
659                 /* Tell the user that the ECC+old driver combination can be bad */
660                 sprintf(sbuf,
661                         "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0.\n"
662                         "      A known bug in this driver version can cause performance loss.\n"
663                         "      However, the polling wait workaround can not be used because\n%s\n"
664                         "      Consider updating the driver or turning ECC off.",
665                         (bX86 && bTMPIAtomics) ?
666                             "      GPU(s) are being oversubscribed." :
667                             "      atomic operations are not supported by the platform/CPU+compiler.");
668                 md_print_warn(fplog, sbuf);
669             }
670         }
671     }
672     else
673     {
674         if (bNoStreamSync)
675         {
676             nb->bUseStreamSync = false;
677
678             md_print_warn(fplog,
679                           "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
680         }
681         else
682         {
683             /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
684             nb->bUseStreamSync = true;
685         }
686     }
687
688     /* CUDA timing disabled as event timers don't work:
689        - with multiple streams = domain-decomposition;
690        - with the polling waiting hack (without cudaStreamSynchronize);
691        - when turned off by GMX_DISABLE_CUDA_TIMING.
692      */
693     nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
694                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
695
696     if (nb->bDoTime)
697     {
698         init_timers(nb->timers, nb->bUseTwoStreams);
699         init_timings(nb->timings);
700     }
701
702     /* set the kernel type for the current GPU */
703     nb->kernel_ver = pick_nbnxn_kernel_version(fplog, nb->dev_info);
704     /* pick L1 cache configuration */
705     nbnxn_cuda_set_cacheconfig(nb->dev_info);
706
707     *p_cu_nb = nb;
708
709     if (debug)
710     {
711         fprintf(debug, "Initialized CUDA data structures.\n");
712     }
713 }
714
715 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t                cu_nb,
716                            const interaction_const_t      *ic,
717                            const nonbonded_verlet_group_t *nbv_group)
718 {
719     init_atomdata_first(cu_nb->atdat, nbv_group[0].nbat->ntype);
720     init_nbparam(cu_nb->nbparam, ic, nbv_group[0].nbat, cu_nb->dev_info);
721
722     /* clear energy and shift force outputs */
723     nbnxn_cuda_clear_e_fshift(cu_nb);
724 }
725
726 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t cu_nb,
727                               const nbnxn_pairlist_t *h_plist,
728                               int iloc)
729 {
730     char         sbuf[STRLEN];
731     cudaError_t  stat;
732     bool         bDoTime    = cu_nb->bDoTime;
733     cudaStream_t stream     = cu_nb->stream[iloc];
734     cu_plist_t   *d_plist   = cu_nb->plist[iloc];
735
736     if (d_plist->na_c < 0)
737     {
738         d_plist->na_c = h_plist->na_ci;
739     }
740     else
741     {
742         if (d_plist->na_c != h_plist->na_ci)
743         {
744             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
745                     d_plist->na_c, h_plist->na_ci);
746             gmx_incons(sbuf);
747         }
748     }
749
750     if (bDoTime)
751     {
752         stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
753         CU_RET_ERR(stat, "cudaEventRecord failed");
754     }
755
756     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
757                          &d_plist->nsci, &d_plist->sci_nalloc,
758                          h_plist->nsci,
759                          stream, true);
760
761     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
762                          &d_plist->ncj4, &d_plist->cj4_nalloc,
763                          h_plist->ncj4,
764                          stream, true);
765
766     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
767                          &d_plist->nexcl, &d_plist->excl_nalloc,
768                          h_plist->nexcl,
769                          stream, true);
770
771     if (bDoTime)
772     {
773         stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
774         CU_RET_ERR(stat, "cudaEventRecord failed");
775     }
776
777     /* need to prune the pair list during the next step */
778     d_plist->bDoPrune = true;
779 }
780
781 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t cu_nb,
782                                 const nbnxn_atomdata_t *nbatom)
783 {
784     cu_atomdata_t *adat = cu_nb->atdat;
785     cudaStream_t  ls    = cu_nb->stream[eintLocal];
786
787     /* only if we have a dynamic box */
788     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
789     {
790         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec, 
791                           SHIFTS * sizeof(*adat->shift_vec), ls);
792         adat->bShiftVecUploaded = true;
793     }
794 }
795
796 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
797 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
798 {
799     cudaError_t   stat;
800     cu_atomdata_t *adat = cu_nb->atdat;
801     cudaStream_t  ls    = cu_nb->stream[eintLocal];
802
803     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
804     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
805 }
806
807 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
808 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
809 {
810     cudaError_t   stat;
811     cu_atomdata_t *adat = cu_nb->atdat;
812     cudaStream_t  ls    = cu_nb->stream[eintLocal];
813
814     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
815     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
816     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
817     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
818     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
819     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
820 }
821
822 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
823 {
824     nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
825     /* clear shift force array and energies if the outputs were 
826        used in the current step */
827     if (flags & GMX_FORCE_VIRIAL)
828     {
829         nbnxn_cuda_clear_e_fshift(cu_nb);
830     }
831 }
832
833 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t cu_nb,
834                               const nbnxn_atomdata_t *nbat)
835 {
836     cudaError_t   stat;
837     int           nalloc, natoms;
838     bool          realloced;
839     bool          bDoTime   = cu_nb->bDoTime;
840     cu_timers_t   *timers   = cu_nb->timers;
841     cu_atomdata_t *d_atdat  = cu_nb->atdat;
842     cudaStream_t  ls        = cu_nb->stream[eintLocal];
843
844     natoms = nbat->natoms;
845     realloced = false;
846
847     if (bDoTime)
848     {
849         /* time async copy */
850         stat = cudaEventRecord(timers->start_atdat, ls);
851         CU_RET_ERR(stat, "cudaEventRecord failed");
852     }
853
854     /* need to reallocate if we have to copy more atoms than the amount of space
855        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
856     if (natoms > d_atdat->nalloc)
857     {
858         nalloc = over_alloc_small(natoms);
859
860         /* free up first if the arrays have already been initialized */
861         if (d_atdat->nalloc != -1)
862         {
863             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
864             cu_free_buffered(d_atdat->xq);
865             cu_free_buffered(d_atdat->atom_types);
866         }
867
868         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
869         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
870         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
871         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
872
873         stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
874         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
875
876         d_atdat->nalloc = nalloc;
877         realloced = true;
878     }
879
880     d_atdat->natoms = natoms;
881     d_atdat->natoms_local = nbat->natoms_local;
882
883     /* need to clear GPU f output if realloc happened */
884     if (realloced)
885     {
886         nbnxn_cuda_clear_f(cu_nb, nalloc);
887     }
888
889     cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
890                       natoms*sizeof(*d_atdat->atom_types), ls);
891
892     if (bDoTime)
893     {
894         stat = cudaEventRecord(timers->stop_atdat, ls);
895         CU_RET_ERR(stat, "cudaEventRecord failed");
896     }
897 }
898
899 void nbnxn_cuda_free(nbnxn_cuda_ptr_t cu_nb)
900 {
901     cudaError_t     stat;
902     cu_atomdata_t   *atdat;
903     cu_nbparam_t    *nbparam;
904     cu_plist_t      *plist, *plist_nl;
905     cu_timers_t     *timers;
906
907     if (cu_nb == NULL) return;
908
909     atdat       = cu_nb->atdat;
910     nbparam     = cu_nb->nbparam;
911     plist       = cu_nb->plist[eintLocal];
912     plist_nl    = cu_nb->plist[eintNonlocal];
913     timers      = cu_nb->timers;
914
915     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
916     {
917
918 #ifdef TEXOBJ_SUPPORTED
919         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
920         if (cu_nb->dev_info->prop.major >= 3)
921         {
922             stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
923             CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
924         }
925         else
926 #endif
927         {
928             stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
929             CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
930         }
931         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
932     }
933
934     stat = cudaEventDestroy(cu_nb->nonlocal_done);
935     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
936     stat = cudaEventDestroy(cu_nb->misc_ops_done);
937     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
938
939     if (cu_nb->bDoTime)
940     {
941         stat = cudaEventDestroy(timers->start_atdat);
942         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
943         stat = cudaEventDestroy(timers->stop_atdat);
944         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
945
946         /* The non-local counters/stream (second in the array) are needed only with DD. */
947         for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
948         {
949             stat = cudaEventDestroy(timers->start_nb_k[i]);
950             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
951             stat = cudaEventDestroy(timers->stop_nb_k[i]);
952             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
953
954             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
955             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
956             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
957             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
958
959             stat = cudaStreamDestroy(cu_nb->stream[i]);
960             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
961
962             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
963             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
964             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
965             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
966
967             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
968             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
969             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
970             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
971         }
972     }
973
974 #ifdef TEXOBJ_SUPPORTED
975     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
976     if (cu_nb->dev_info->prop.major >= 3)
977     {
978         stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
979         CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
980     }
981     else
982 #endif
983     {
984         stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
985         CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
986     }
987     cu_free_buffered(nbparam->nbfp);
988
989     stat = cudaFree(atdat->shift_vec);
990     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
991     stat = cudaFree(atdat->fshift);
992     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
993
994     stat = cudaFree(atdat->e_lj);
995     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
996     stat = cudaFree(atdat->e_el);
997     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
998
999     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
1000     cu_free_buffered(atdat->xq);
1001     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
1002
1003     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
1004     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
1005     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
1006     if (cu_nb->bUseTwoStreams)
1007     {
1008         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
1009         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
1010         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
1011     }
1012
1013     sfree(atdat);
1014     sfree(nbparam);
1015     sfree(plist);
1016     if (cu_nb->bUseTwoStreams)
1017     {
1018         sfree(plist_nl);
1019     }
1020     sfree(timers);
1021     sfree(cu_nb->timings);
1022     sfree(cu_nb);
1023
1024     if (debug)
1025     {
1026         fprintf(debug, "Cleaned up CUDA data structures.\n");
1027     }
1028 }
1029
1030 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
1031 {
1032     cudaError_t stat;
1033     cudaStream_t stream = cu_nb->stream[iloc];
1034
1035     stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
1036     CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
1037 }
1038
1039 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
1040 {
1041     return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
1042 }
1043
1044 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
1045 {
1046     if (cu_nb->bDoTime)
1047     {
1048         init_timings(cu_nb->timings);
1049     }
1050 }
1051
1052 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
1053 {
1054     return cu_nb != NULL ?
1055         gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;
1056
1057 }
1058
1059 gmx_bool nbnxn_cuda_is_kernel_ewald_analytical(const nbnxn_cuda_ptr_t cu_nb)
1060 {
1061     return ((cu_nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1062             (cu_nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
1063 }