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