Uncrustify all files
[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)
451     {
452         return;
453     }
454
455     snew(nb, 1);
456     snew(nb->atdat, 1);
457     snew(nb->nbparam, 1);
458     snew(nb->plist[eintLocal], 1);
459     if (bLocalAndNonlocal)
460     {
461         snew(nb->plist[eintNonlocal], 1);
462     }
463
464     nb->bUseTwoStreams = bLocalAndNonlocal;
465
466     snew(nb->timers, 1);
467     snew(nb->timings, 1);
468
469     /* init nbst */
470     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
471     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
472     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
473
474     init_plist(nb->plist[eintLocal]);
475
476     /* set device info, just point it to the right GPU among the detected ones */
477     nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, gpu_opt, my_gpu_index)];
478
479     /* local/non-local GPU streams */
480     stat = cudaStreamCreate(&nb->stream[eintLocal]);
481     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
482     if (nb->bUseTwoStreams)
483     {
484         init_plist(nb->plist[eintNonlocal]);
485
486         /* CUDA stream priority available in the CUDA RT 5.5 API.
487          * Note that the device we're running on does not have to support
488          * priorities, because we are querying the priority range which in this
489          * case will be a single value.
490          */
491 #if CUDA_VERSION >= 5500
492         {
493             int highest_priority;
494             stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
495             CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
496
497             stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
498                                                 cudaStreamDefault,
499                                                 highest_priority);
500             CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
501         }
502 #else
503         stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
504         CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
505 #endif
506     }
507
508     /* init events for sychronization (timing disabled for performance reasons!) */
509     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
510     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
511     stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
512     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
513
514     /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
515      * (which increases with shorter time/step) caused by a known CUDA driver bug.
516      * To work around the issue we'll use an (admittedly fragile) memory polling
517      * waiting to preserve performance. This requires support for atomic
518      * operations and only works on x86/x86_64.
519      * With polling wait event-timing also needs to be disabled.
520      *
521      * The overhead is greatly reduced in API v5.0 drivers and the improvement
522      * is independent of runtime version. Hence, with API v5.0 drivers and later
523      * we won't switch to polling.
524      *
525      * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
526      * ranks so we will also disable it in that case.
527      */
528
529     bStreamSync    = getenv("GMX_CUDA_STREAMSYNC") != NULL;
530     bNoStreamSync  = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
531
532 #ifdef TMPI_ATOMICS
533     bTMPIAtomics = true;
534 #else
535     bTMPIAtomics = false;
536 #endif
537
538 #ifdef GMX_TARGET_X86
539     bX86 = true;
540 #else
541     bX86 = false;
542 #endif
543
544     if (bStreamSync && bNoStreamSync)
545     {
546         gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
547     }
548
549     stat = cudaDriverGetVersion(&cuda_drv_ver);
550     CU_RET_ERR(stat, "cudaDriverGetVersion failed");
551
552     bOldDriver = (cuda_drv_ver < 5000);
553
554     if ((nb->dev_info->prop.ECCEnabled == 1) && bOldDriver)
555     {
556         /* Polling wait should be used instead of cudaStreamSynchronize only if:
557          *   - ECC is ON & driver is old (checked above),
558          *   - we're on x86/x86_64,
559          *   - atomics are available, and
560          *   - GPUs are not being shared.
561          */
562         bool bShouldUsePollSync = (bX86 && bTMPIAtomics &&
563                                    (gmx_count_gpu_dev_shared(gpu_opt) < 1));
564
565         if (bStreamSync)
566         {
567             nb->bUseStreamSync = true;
568
569             /* only warn if polling should be used */
570             if (bShouldUsePollSync)
571             {
572                 md_print_warn(fplog,
573                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
574                               "      cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
575             }
576         }
577         else
578         {
579             nb->bUseStreamSync = !bShouldUsePollSync;
580
581             if (bShouldUsePollSync)
582             {
583                 md_print_warn(fplog,
584                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
585                               "      cause performance loss. Switching to the alternative polling GPU wait.\n"
586                               "      If you encounter issues, switch back to standard GPU waiting by setting\n"
587                               "      the GMX_CUDA_STREAMSYNC environment variable.\n");
588             }
589             else
590             {
591                 /* Tell the user that the ECC+old driver combination can be bad */
592                 sprintf(sbuf,
593                         "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0.\n"
594                         "      A known bug in this driver version can cause performance loss.\n"
595                         "      However, the polling wait workaround can not be used because\n%s\n"
596                         "      Consider updating the driver or turning ECC off.",
597                         (bX86 && bTMPIAtomics) ?
598                         "      GPU(s) are being oversubscribed." :
599                         "      atomic operations are not supported by the platform/CPU+compiler.");
600                 md_print_warn(fplog, sbuf);
601             }
602         }
603     }
604     else
605     {
606         if (bNoStreamSync)
607         {
608             nb->bUseStreamSync = false;
609
610             md_print_warn(fplog,
611                           "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
612         }
613         else
614         {
615             /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
616             nb->bUseStreamSync = true;
617         }
618     }
619
620     /* CUDA timing disabled as event timers don't work:
621        - with multiple streams = domain-decomposition;
622        - with the polling waiting hack (without cudaStreamSynchronize);
623        - when turned off by GMX_DISABLE_CUDA_TIMING.
624      */
625     nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
626                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
627
628     if (nb->bDoTime)
629     {
630         init_timers(nb->timers, nb->bUseTwoStreams);
631         init_timings(nb->timings);
632     }
633
634     /* set the kernel type for the current GPU */
635     /* pick L1 cache configuration */
636     nbnxn_cuda_set_cacheconfig(nb->dev_info);
637
638     *p_cu_nb = nb;
639
640     if (debug)
641     {
642         fprintf(debug, "Initialized CUDA data structures.\n");
643     }
644 }
645
646 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t                cu_nb,
647                            const interaction_const_t      *ic,
648                            const nonbonded_verlet_group_t *nbv_group)
649 {
650     init_atomdata_first(cu_nb->atdat, nbv_group[0].nbat->ntype);
651     init_nbparam(cu_nb->nbparam, ic, nbv_group[0].nbat, cu_nb->dev_info);
652
653     /* clear energy and shift force outputs */
654     nbnxn_cuda_clear_e_fshift(cu_nb);
655 }
656
657 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t        cu_nb,
658                               const nbnxn_pairlist_t *h_plist,
659                               int                     iloc)
660 {
661     char          sbuf[STRLEN];
662     cudaError_t   stat;
663     bool          bDoTime    = cu_nb->bDoTime;
664     cudaStream_t  stream     = cu_nb->stream[iloc];
665     cu_plist_t   *d_plist    = cu_nb->plist[iloc];
666
667     if (d_plist->na_c < 0)
668     {
669         d_plist->na_c = h_plist->na_ci;
670     }
671     else
672     {
673         if (d_plist->na_c != h_plist->na_ci)
674         {
675             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
676                     d_plist->na_c, h_plist->na_ci);
677             gmx_incons(sbuf);
678         }
679     }
680
681     if (bDoTime)
682     {
683         stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
684         CU_RET_ERR(stat, "cudaEventRecord failed");
685     }
686
687     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
688                         &d_plist->nsci, &d_plist->sci_nalloc,
689                         h_plist->nsci,
690                         stream, true);
691
692     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
693                         &d_plist->ncj4, &d_plist->cj4_nalloc,
694                         h_plist->ncj4,
695                         stream, true);
696
697     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
698                         &d_plist->nexcl, &d_plist->excl_nalloc,
699                         h_plist->nexcl,
700                         stream, true);
701
702     if (bDoTime)
703     {
704         stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
705         CU_RET_ERR(stat, "cudaEventRecord failed");
706     }
707
708     /* need to prune the pair list during the next step */
709     d_plist->bDoPrune = true;
710 }
711
712 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t        cu_nb,
713                                 const nbnxn_atomdata_t *nbatom)
714 {
715     cu_atomdata_t *adat  = cu_nb->atdat;
716     cudaStream_t   ls    = cu_nb->stream[eintLocal];
717
718     /* only if we have a dynamic box */
719     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
720     {
721         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
722                           SHIFTS * sizeof(*adat->shift_vec), ls);
723         adat->bShiftVecUploaded = true;
724     }
725 }
726
727 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
728 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
729 {
730     cudaError_t    stat;
731     cu_atomdata_t *adat  = cu_nb->atdat;
732     cudaStream_t   ls    = cu_nb->stream[eintLocal];
733
734     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
735     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
736 }
737
738 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
739 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
740 {
741     cudaError_t    stat;
742     cu_atomdata_t *adat  = cu_nb->atdat;
743     cudaStream_t   ls    = cu_nb->stream[eintLocal];
744
745     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
746     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
747     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
748     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
749     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
750     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
751 }
752
753 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
754 {
755     nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
756     /* clear shift force array and energies if the outputs were
757        used in the current step */
758     if (flags & GMX_FORCE_VIRIAL)
759     {
760         nbnxn_cuda_clear_e_fshift(cu_nb);
761     }
762 }
763
764 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t        cu_nb,
765                               const nbnxn_atomdata_t *nbat)
766 {
767     cudaError_t    stat;
768     int            nalloc, natoms;
769     bool           realloced;
770     bool           bDoTime   = cu_nb->bDoTime;
771     cu_timers_t   *timers    = cu_nb->timers;
772     cu_atomdata_t *d_atdat   = cu_nb->atdat;
773     cudaStream_t   ls        = cu_nb->stream[eintLocal];
774
775     natoms    = nbat->natoms;
776     realloced = false;
777
778     if (bDoTime)
779     {
780         /* time async copy */
781         stat = cudaEventRecord(timers->start_atdat, ls);
782         CU_RET_ERR(stat, "cudaEventRecord failed");
783     }
784
785     /* need to reallocate if we have to copy more atoms than the amount of space
786        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
787     if (natoms > d_atdat->nalloc)
788     {
789         nalloc = over_alloc_small(natoms);
790
791         /* free up first if the arrays have already been initialized */
792         if (d_atdat->nalloc != -1)
793         {
794             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
795             cu_free_buffered(d_atdat->xq);
796             cu_free_buffered(d_atdat->atom_types);
797         }
798
799         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
800         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
801         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
802         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
803
804         stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
805         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
806
807         d_atdat->nalloc = nalloc;
808         realloced       = true;
809     }
810
811     d_atdat->natoms       = natoms;
812     d_atdat->natoms_local = nbat->natoms_local;
813
814     /* need to clear GPU f output if realloc happened */
815     if (realloced)
816     {
817         nbnxn_cuda_clear_f(cu_nb, nalloc);
818     }
819
820     cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
821                       natoms*sizeof(*d_atdat->atom_types), ls);
822
823     if (bDoTime)
824     {
825         stat = cudaEventRecord(timers->stop_atdat, ls);
826         CU_RET_ERR(stat, "cudaEventRecord failed");
827     }
828 }
829
830 void nbnxn_cuda_free(nbnxn_cuda_ptr_t cu_nb)
831 {
832     cudaError_t      stat;
833     cu_atomdata_t   *atdat;
834     cu_nbparam_t    *nbparam;
835     cu_plist_t      *plist, *plist_nl;
836     cu_timers_t     *timers;
837
838     if (cu_nb == NULL)
839     {
840         return;
841     }
842
843     atdat       = cu_nb->atdat;
844     nbparam     = cu_nb->nbparam;
845     plist       = cu_nb->plist[eintLocal];
846     plist_nl    = cu_nb->plist[eintNonlocal];
847     timers      = cu_nb->timers;
848
849     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
850     {
851
852 #ifdef TEXOBJ_SUPPORTED
853         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
854         if (cu_nb->dev_info->prop.major >= 3)
855         {
856             stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
857             CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
858         }
859         else
860 #endif
861         {
862             stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
863             CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
864         }
865         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
866     }
867
868     stat = cudaEventDestroy(cu_nb->nonlocal_done);
869     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
870     stat = cudaEventDestroy(cu_nb->misc_ops_done);
871     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
872
873     if (cu_nb->bDoTime)
874     {
875         stat = cudaEventDestroy(timers->start_atdat);
876         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
877         stat = cudaEventDestroy(timers->stop_atdat);
878         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
879
880         /* The non-local counters/stream (second in the array) are needed only with DD. */
881         for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
882         {
883             stat = cudaEventDestroy(timers->start_nb_k[i]);
884             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
885             stat = cudaEventDestroy(timers->stop_nb_k[i]);
886             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
887
888             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
889             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
890             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
891             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
892
893             stat = cudaStreamDestroy(cu_nb->stream[i]);
894             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
895
896             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
897             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
898             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
899             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
900
901             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
902             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
903             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
904             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
905         }
906     }
907
908 #ifdef TEXOBJ_SUPPORTED
909     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
910     if (cu_nb->dev_info->prop.major >= 3)
911     {
912         stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
913         CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
914     }
915     else
916 #endif
917     {
918         stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
919         CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
920     }
921     cu_free_buffered(nbparam->nbfp);
922
923     stat = cudaFree(atdat->shift_vec);
924     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
925     stat = cudaFree(atdat->fshift);
926     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
927
928     stat = cudaFree(atdat->e_lj);
929     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
930     stat = cudaFree(atdat->e_el);
931     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
932
933     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
934     cu_free_buffered(atdat->xq);
935     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
936
937     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
938     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
939     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
940     if (cu_nb->bUseTwoStreams)
941     {
942         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
943         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
944         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
945     }
946
947     sfree(atdat);
948     sfree(nbparam);
949     sfree(plist);
950     if (cu_nb->bUseTwoStreams)
951     {
952         sfree(plist_nl);
953     }
954     sfree(timers);
955     sfree(cu_nb->timings);
956     sfree(cu_nb);
957
958     if (debug)
959     {
960         fprintf(debug, "Cleaned up CUDA data structures.\n");
961     }
962 }
963
964 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
965 {
966     cudaError_t  stat;
967     cudaStream_t stream = cu_nb->stream[iloc];
968
969     stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
970     CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
971 }
972
973 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
974 {
975     return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
976 }
977
978 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
979 {
980     if (cu_nb->bDoTime)
981     {
982         init_timings(cu_nb->timings);
983     }
984 }
985
986 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
987 {
988     return cu_nb != NULL ?
989            gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;
990
991 }
992
993 gmx_bool nbnxn_cuda_is_kernel_ewald_analytical(const nbnxn_cuda_ptr_t cu_nb)
994 {
995     return ((cu_nb->nbparam->eeltype == eelCuEWALD_ANA) ||
996             (cu_nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
997 }