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