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