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