Implemented nbnxn LJ switch functions
[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,2014, 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     /* TODO: implemented LJ force- and potential-switch CUDA kernels */
283     if (!(ic->vdw_modifier == eintmodNONE ||
284           ic->vdw_modifier == eintmodPOTSHIFT))
285     {
286         gmx_fatal(FARGS, "The CUDA kernels do not yet support switched LJ interactions");
287     }
288
289     if (ic->eeltype == eelCUT)
290     {
291         nbp->eeltype = eelCuCUT;
292     }
293     else if (EEL_RF(ic->eeltype))
294     {
295         nbp->eeltype = eelCuRF;
296     }
297     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
298     {
299         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
300         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
301     }
302     else
303     {
304         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
305         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
306     }
307
308     /* generate table for PME */
309     nbp->coulomb_tab = NULL;
310     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
311     {
312         init_ewald_coulomb_force_table(nbp, dev_info);
313     }
314
315     nnbfp = 2*ntypes*ntypes;
316     stat  = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
317     CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
318     cu_copy_H2D(nbp->nbfp, nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
319
320 #ifdef TEXOBJ_SUPPORTED
321     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
322     if (dev_info->prop.major >= 3)
323     {
324         cudaResourceDesc rd;
325         memset(&rd, 0, sizeof(rd));
326         rd.resType                  = cudaResourceTypeLinear;
327         rd.res.linear.devPtr        = nbp->nbfp;
328         rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
329         rd.res.linear.desc.x        = 32;
330         rd.res.linear.sizeInBytes   = nnbfp*sizeof(*nbp->nbfp);
331
332         cudaTextureDesc td;
333         memset(&td, 0, sizeof(td));
334         td.readMode                 = cudaReadModeElementType;
335         stat = cudaCreateTextureObject(&nbp->nbfp_texobj, &rd, &td, NULL);
336         CU_RET_ERR(stat, "cudaCreateTextureObject on nbfp_texobj failed");
337     }
338     else
339 #endif
340     {
341         cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
342         stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
343                                nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
344         CU_RET_ERR(stat, "cudaBindTexture on nbfp_texref failed");
345     }
346 }
347
348 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
349  *  electrostatic type switching to twin cut-off (or back) if needed. */
350 void nbnxn_cuda_pme_loadbal_update_param(nbnxn_cuda_ptr_t           cu_nb,
351                                          const interaction_const_t *ic)
352 {
353     cu_nbparam_t *nbp = cu_nb->nbparam;
354
355     nbp->rlist_sq       = ic->rlist * ic->rlist;
356     nbp->rcoulomb_sq    = ic->rcoulomb * ic->rcoulomb;
357     nbp->ewald_beta     = ic->ewaldcoeff_q;
358
359     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
360                                                  cu_nb->dev_info);
361
362     init_ewald_coulomb_force_table(cu_nb->nbparam, cu_nb->dev_info);
363 }
364
365 /*! Initializes the pair list data structure. */
366 static void init_plist(cu_plist_t *pl)
367 {
368     /* initialize to NULL pointers to data that is not allocated here and will
369        need reallocation in nbnxn_cuda_init_pairlist */
370     pl->sci     = NULL;
371     pl->cj4     = NULL;
372     pl->excl    = NULL;
373
374     /* size -1 indicates that the respective array hasn't been initialized yet */
375     pl->na_c        = -1;
376     pl->nsci        = -1;
377     pl->sci_nalloc  = -1;
378     pl->ncj4        = -1;
379     pl->cj4_nalloc  = -1;
380     pl->nexcl       = -1;
381     pl->excl_nalloc = -1;
382     pl->bDoPrune    = false;
383 }
384
385 /*! Initializes the timer data structure. */
386 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
387 {
388     cudaError_t stat;
389     int         eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
390
391     stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
392     CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
393     stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
394     CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
395
396     /* The non-local counters/stream (second in the array) are needed only with DD. */
397     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
398     {
399         stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
400         CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
401         stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
402         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
403
404
405         stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
406         CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
407         stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
408         CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
409
410         stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
411         CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
412         stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
413         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
414
415         stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
416         CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
417         stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
418         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
419     }
420 }
421
422 /*! Initializes the timings data structure. */
423 static void init_timings(wallclock_gpu_t *t)
424 {
425     int i, j;
426
427     t->nb_h2d_t = 0.0;
428     t->nb_d2h_t = 0.0;
429     t->nb_c     = 0;
430     t->pl_h2d_t = 0.0;
431     t->pl_h2d_c = 0;
432     for (i = 0; i < 2; i++)
433     {
434         for (j = 0; j < 2; j++)
435         {
436             t->ktime[i][j].t = 0.0;
437             t->ktime[i][j].c = 0;
438         }
439     }
440 }
441
442 void nbnxn_cuda_init(FILE                 *fplog,
443                      nbnxn_cuda_ptr_t     *p_cu_nb,
444                      const gmx_gpu_info_t *gpu_info,
445                      const gmx_gpu_opt_t  *gpu_opt,
446                      int                   my_gpu_index,
447                      gmx_bool              bLocalAndNonlocal)
448 {
449     cudaError_t       stat;
450     nbnxn_cuda_ptr_t  nb;
451     char              sbuf[STRLEN];
452     bool              bStreamSync, bNoStreamSync, bTMPIAtomics, bX86, bOldDriver;
453     int               cuda_drv_ver;
454
455     assert(gpu_info);
456
457     if (p_cu_nb == NULL)
458     {
459         return;
460     }
461
462     snew(nb, 1);
463     snew(nb->atdat, 1);
464     snew(nb->nbparam, 1);
465     snew(nb->plist[eintLocal], 1);
466     if (bLocalAndNonlocal)
467     {
468         snew(nb->plist[eintNonlocal], 1);
469     }
470
471     nb->bUseTwoStreams = bLocalAndNonlocal;
472
473     snew(nb->timers, 1);
474     snew(nb->timings, 1);
475
476     /* init nbst */
477     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
478     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
479     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
480
481     init_plist(nb->plist[eintLocal]);
482
483     /* set device info, just point it to the right GPU among the detected ones */
484     nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, gpu_opt, my_gpu_index)];
485
486     /* local/non-local GPU streams */
487     stat = cudaStreamCreate(&nb->stream[eintLocal]);
488     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
489     if (nb->bUseTwoStreams)
490     {
491         init_plist(nb->plist[eintNonlocal]);
492
493         /* CUDA stream priority available in the CUDA RT 5.5 API.
494          * Note that the device we're running on does not have to support
495          * priorities, because we are querying the priority range which in this
496          * case will be a single value.
497          */
498 #if CUDA_VERSION >= 5500
499         {
500             int highest_priority;
501             stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
502             CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
503
504             stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
505                                                 cudaStreamDefault,
506                                                 highest_priority);
507             CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
508         }
509 #else
510         stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
511         CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
512 #endif
513     }
514
515     /* init events for sychronization (timing disabled for performance reasons!) */
516     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
517     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
518     stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
519     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
520
521     /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
522      * (which increases with shorter time/step) caused by a known CUDA driver bug.
523      * To work around the issue we'll use an (admittedly fragile) memory polling
524      * waiting to preserve performance. This requires support for atomic
525      * operations and only works on x86/x86_64.
526      * With polling wait event-timing also needs to be disabled.
527      *
528      * The overhead is greatly reduced in API v5.0 drivers and the improvement
529      * is independent of runtime version. Hence, with API v5.0 drivers and later
530      * we won't switch to polling.
531      *
532      * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
533      * ranks so we will also disable it in that case.
534      */
535
536     bStreamSync    = getenv("GMX_CUDA_STREAMSYNC") != NULL;
537     bNoStreamSync  = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
538
539 #ifdef TMPI_ATOMICS
540     bTMPIAtomics = true;
541 #else
542     bTMPIAtomics = false;
543 #endif
544
545 #ifdef GMX_TARGET_X86
546     bX86 = true;
547 #else
548     bX86 = false;
549 #endif
550
551     if (bStreamSync && bNoStreamSync)
552     {
553         gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
554     }
555
556     stat = cudaDriverGetVersion(&cuda_drv_ver);
557     CU_RET_ERR(stat, "cudaDriverGetVersion failed");
558
559     bOldDriver = (cuda_drv_ver < 5000);
560
561     if ((nb->dev_info->prop.ECCEnabled == 1) && bOldDriver)
562     {
563         /* Polling wait should be used instead of cudaStreamSynchronize only if:
564          *   - ECC is ON & driver is old (checked above),
565          *   - we're on x86/x86_64,
566          *   - atomics are available, and
567          *   - GPUs are not being shared.
568          */
569         bool bShouldUsePollSync = (bX86 && bTMPIAtomics &&
570                                    (gmx_count_gpu_dev_shared(gpu_opt) < 1));
571
572         if (bStreamSync)
573         {
574             nb->bUseStreamSync = true;
575
576             /* only warn if polling should be used */
577             if (bShouldUsePollSync)
578             {
579                 md_print_warn(fplog,
580                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
581                               "      cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
582             }
583         }
584         else
585         {
586             nb->bUseStreamSync = !bShouldUsePollSync;
587
588             if (bShouldUsePollSync)
589             {
590                 md_print_warn(fplog,
591                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
592                               "      cause performance loss. Switching to the alternative polling GPU wait.\n"
593                               "      If you encounter issues, switch back to standard GPU waiting by setting\n"
594                               "      the GMX_CUDA_STREAMSYNC environment variable.\n");
595             }
596             else
597             {
598                 /* Tell the user that the ECC+old driver combination can be bad */
599                 sprintf(sbuf,
600                         "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0.\n"
601                         "      A known bug in this driver version can cause performance loss.\n"
602                         "      However, the polling wait workaround can not be used because\n%s\n"
603                         "      Consider updating the driver or turning ECC off.",
604                         (bX86 && bTMPIAtomics) ?
605                         "      GPU(s) are being oversubscribed." :
606                         "      atomic operations are not supported by the platform/CPU+compiler.");
607                 md_print_warn(fplog, sbuf);
608             }
609         }
610     }
611     else
612     {
613         if (bNoStreamSync)
614         {
615             nb->bUseStreamSync = false;
616
617             md_print_warn(fplog,
618                           "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
619         }
620         else
621         {
622             /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
623             nb->bUseStreamSync = true;
624         }
625     }
626
627     /* CUDA timing disabled as event timers don't work:
628        - with multiple streams = domain-decomposition;
629        - with the polling waiting hack (without cudaStreamSynchronize);
630        - when turned off by GMX_DISABLE_CUDA_TIMING.
631      */
632     nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
633                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
634
635     if (nb->bDoTime)
636     {
637         init_timers(nb->timers, nb->bUseTwoStreams);
638         init_timings(nb->timings);
639     }
640
641     /* set the kernel type for the current GPU */
642     /* pick L1 cache configuration */
643     nbnxn_cuda_set_cacheconfig(nb->dev_info);
644
645     *p_cu_nb = nb;
646
647     if (debug)
648     {
649         fprintf(debug, "Initialized CUDA data structures.\n");
650     }
651 }
652
653 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t                cu_nb,
654                            const interaction_const_t      *ic,
655                            const nonbonded_verlet_group_t *nbv_group)
656 {
657     init_atomdata_first(cu_nb->atdat, nbv_group[0].nbat->ntype);
658     init_nbparam(cu_nb->nbparam, ic, nbv_group[0].nbat, cu_nb->dev_info);
659
660     /* clear energy and shift force outputs */
661     nbnxn_cuda_clear_e_fshift(cu_nb);
662 }
663
664 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t        cu_nb,
665                               const nbnxn_pairlist_t *h_plist,
666                               int                     iloc)
667 {
668     char          sbuf[STRLEN];
669     cudaError_t   stat;
670     bool          bDoTime    = cu_nb->bDoTime;
671     cudaStream_t  stream     = cu_nb->stream[iloc];
672     cu_plist_t   *d_plist    = cu_nb->plist[iloc];
673
674     if (d_plist->na_c < 0)
675     {
676         d_plist->na_c = h_plist->na_ci;
677     }
678     else
679     {
680         if (d_plist->na_c != h_plist->na_ci)
681         {
682             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
683                     d_plist->na_c, h_plist->na_ci);
684             gmx_incons(sbuf);
685         }
686     }
687
688     if (bDoTime)
689     {
690         stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
691         CU_RET_ERR(stat, "cudaEventRecord failed");
692     }
693
694     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
695                         &d_plist->nsci, &d_plist->sci_nalloc,
696                         h_plist->nsci,
697                         stream, true);
698
699     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
700                         &d_plist->ncj4, &d_plist->cj4_nalloc,
701                         h_plist->ncj4,
702                         stream, true);
703
704     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
705                         &d_plist->nexcl, &d_plist->excl_nalloc,
706                         h_plist->nexcl,
707                         stream, true);
708
709     if (bDoTime)
710     {
711         stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
712         CU_RET_ERR(stat, "cudaEventRecord failed");
713     }
714
715     /* need to prune the pair list during the next step */
716     d_plist->bDoPrune = true;
717 }
718
719 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t        cu_nb,
720                                 const nbnxn_atomdata_t *nbatom)
721 {
722     cu_atomdata_t *adat  = cu_nb->atdat;
723     cudaStream_t   ls    = cu_nb->stream[eintLocal];
724
725     /* only if we have a dynamic box */
726     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
727     {
728         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
729                           SHIFTS * sizeof(*adat->shift_vec), ls);
730         adat->bShiftVecUploaded = true;
731     }
732 }
733
734 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
735 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
736 {
737     cudaError_t    stat;
738     cu_atomdata_t *adat  = cu_nb->atdat;
739     cudaStream_t   ls    = cu_nb->stream[eintLocal];
740
741     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
742     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
743 }
744
745 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
746 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
747 {
748     cudaError_t    stat;
749     cu_atomdata_t *adat  = cu_nb->atdat;
750     cudaStream_t   ls    = cu_nb->stream[eintLocal];
751
752     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
753     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
754     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
755     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
756     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
757     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
758 }
759
760 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
761 {
762     nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
763     /* clear shift force array and energies if the outputs were
764        used in the current step */
765     if (flags & GMX_FORCE_VIRIAL)
766     {
767         nbnxn_cuda_clear_e_fshift(cu_nb);
768     }
769 }
770
771 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t        cu_nb,
772                               const nbnxn_atomdata_t *nbat)
773 {
774     cudaError_t    stat;
775     int            nalloc, natoms;
776     bool           realloced;
777     bool           bDoTime   = cu_nb->bDoTime;
778     cu_timers_t   *timers    = cu_nb->timers;
779     cu_atomdata_t *d_atdat   = cu_nb->atdat;
780     cudaStream_t   ls        = cu_nb->stream[eintLocal];
781
782     natoms    = nbat->natoms;
783     realloced = false;
784
785     if (bDoTime)
786     {
787         /* time async copy */
788         stat = cudaEventRecord(timers->start_atdat, ls);
789         CU_RET_ERR(stat, "cudaEventRecord failed");
790     }
791
792     /* need to reallocate if we have to copy more atoms than the amount of space
793        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
794     if (natoms > d_atdat->nalloc)
795     {
796         nalloc = over_alloc_small(natoms);
797
798         /* free up first if the arrays have already been initialized */
799         if (d_atdat->nalloc != -1)
800         {
801             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
802             cu_free_buffered(d_atdat->xq);
803             cu_free_buffered(d_atdat->atom_types);
804         }
805
806         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
807         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
808         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
809         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
810
811         stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
812         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
813
814         d_atdat->nalloc = nalloc;
815         realloced       = true;
816     }
817
818     d_atdat->natoms       = natoms;
819     d_atdat->natoms_local = nbat->natoms_local;
820
821     /* need to clear GPU f output if realloc happened */
822     if (realloced)
823     {
824         nbnxn_cuda_clear_f(cu_nb, nalloc);
825     }
826
827     cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
828                       natoms*sizeof(*d_atdat->atom_types), ls);
829
830     if (bDoTime)
831     {
832         stat = cudaEventRecord(timers->stop_atdat, ls);
833         CU_RET_ERR(stat, "cudaEventRecord failed");
834     }
835 }
836
837 void nbnxn_cuda_free(nbnxn_cuda_ptr_t cu_nb)
838 {
839     cudaError_t      stat;
840     cu_atomdata_t   *atdat;
841     cu_nbparam_t    *nbparam;
842     cu_plist_t      *plist, *plist_nl;
843     cu_timers_t     *timers;
844
845     if (cu_nb == NULL)
846     {
847         return;
848     }
849
850     atdat       = cu_nb->atdat;
851     nbparam     = cu_nb->nbparam;
852     plist       = cu_nb->plist[eintLocal];
853     plist_nl    = cu_nb->plist[eintNonlocal];
854     timers      = cu_nb->timers;
855
856     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
857     {
858
859 #ifdef TEXOBJ_SUPPORTED
860         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
861         if (cu_nb->dev_info->prop.major >= 3)
862         {
863             stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
864             CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
865         }
866         else
867 #endif
868         {
869             stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
870             CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
871         }
872         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
873     }
874
875     stat = cudaEventDestroy(cu_nb->nonlocal_done);
876     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
877     stat = cudaEventDestroy(cu_nb->misc_ops_done);
878     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
879
880     if (cu_nb->bDoTime)
881     {
882         stat = cudaEventDestroy(timers->start_atdat);
883         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
884         stat = cudaEventDestroy(timers->stop_atdat);
885         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
886
887         /* The non-local counters/stream (second in the array) are needed only with DD. */
888         for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
889         {
890             stat = cudaEventDestroy(timers->start_nb_k[i]);
891             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
892             stat = cudaEventDestroy(timers->stop_nb_k[i]);
893             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
894
895             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
896             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
897             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
898             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
899
900             stat = cudaStreamDestroy(cu_nb->stream[i]);
901             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
902
903             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
904             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
905             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
906             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
907
908             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
909             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
910             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
911             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
912         }
913     }
914
915 #ifdef TEXOBJ_SUPPORTED
916     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
917     if (cu_nb->dev_info->prop.major >= 3)
918     {
919         stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
920         CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
921     }
922     else
923 #endif
924     {
925         stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
926         CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
927     }
928     cu_free_buffered(nbparam->nbfp);
929
930     stat = cudaFree(atdat->shift_vec);
931     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
932     stat = cudaFree(atdat->fshift);
933     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
934
935     stat = cudaFree(atdat->e_lj);
936     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
937     stat = cudaFree(atdat->e_el);
938     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
939
940     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
941     cu_free_buffered(atdat->xq);
942     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
943
944     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
945     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
946     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
947     if (cu_nb->bUseTwoStreams)
948     {
949         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
950         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
951         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
952     }
953
954     sfree(atdat);
955     sfree(nbparam);
956     sfree(plist);
957     if (cu_nb->bUseTwoStreams)
958     {
959         sfree(plist_nl);
960     }
961     sfree(timers);
962     sfree(cu_nb->timings);
963     sfree(cu_nb);
964
965     if (debug)
966     {
967         fprintf(debug, "Cleaned up CUDA data structures.\n");
968     }
969 }
970
971 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
972 {
973     cudaError_t  stat;
974     cudaStream_t stream = cu_nb->stream[iloc];
975
976     stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
977     CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
978 }
979
980 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
981 {
982     return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
983 }
984
985 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
986 {
987     if (cu_nb->bDoTime)
988     {
989         init_timings(cu_nb->timings);
990     }
991 }
992
993 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
994 {
995     return cu_nb != NULL ?
996            gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;
997
998 }
999
1000 gmx_bool nbnxn_cuda_is_kernel_ewald_analytical(const nbnxn_cuda_ptr_t cu_nb)
1001 {
1002     return ((cu_nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1003             (cu_nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
1004 }