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