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