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