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