Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_data_mgmt.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,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 <assert.h>
40 #include <stdarg.h>
41 #include <stdlib.h>
42 #include <stdio.h>
43
44 #include <cuda.h>
45
46 #include "tables.h"
47 #include "typedefs.h"
48 #include "types/enums.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 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/smalloc.h"
65
66 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
67
68 /* This is a heuristically determined parameter for the Fermi architecture for
69  * the minimum size of ci lists by multiplying this constant with the # of
70  * multiprocessors on the current device.
71  */
72 static unsigned int gpu_min_ci_balanced_factor = 40;
73
74 /* Functions from nbnxn_cuda.cu */
75 extern void nbnxn_cuda_set_cacheconfig(cuda_dev_info_t *devinfo);
76 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref();
77 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref();
78 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref();
79
80 /* We should actually be using md_print_warn in md_logging.c,
81  * but we can't include mpi.h in CUDA code.
82  */
83 static void md_print_warn(FILE       *fplog,
84                           const char *fmt, ...)
85 {
86     va_list ap;
87
88     if (fplog != NULL)
89     {
90         /* We should only print to stderr on the master node,
91          * in most cases fplog is only set on the master node, so this works.
92          */
93         va_start(ap, fmt);
94         fprintf(stderr, "\n");
95         vfprintf(stderr, fmt, ap);
96         fprintf(stderr, "\n");
97         va_end(ap);
98
99         va_start(ap, fmt);
100         fprintf(fplog, "\n");
101         vfprintf(fplog, fmt, ap);
102         fprintf(fplog, "\n");
103         va_end(ap);
104     }
105 }
106
107
108 /* Fw. decl. */
109 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb);
110
111
112 /*! Tabulates the Ewald Coulomb force and initializes the size/scale
113     and the table GPU array. If called with an already allocated table,
114     it just re-uploads the table.
115  */
116 static void init_ewald_coulomb_force_table(cu_nbparam_t          *nbp,
117                                            const cuda_dev_info_t *dev_info)
118 {
119     float       *ftmp, *coul_tab;
120     int          tabsize;
121     double       tabscale;
122     cudaError_t  stat;
123
124     tabsize     = GPU_EWALD_COULOMB_FORCE_TABLE_SIZE;
125     /* Subtract 2 iso 1 to avoid access out of range due to rounding */
126     tabscale    = (tabsize - 2) / sqrt(nbp->rcoulomb_sq);
127
128     pmalloc((void**)&ftmp, tabsize*sizeof(*ftmp));
129
130     table_spline3_fill_ewald_lr(ftmp, NULL, NULL, tabsize,
131                                 1/tabscale, nbp->ewald_beta, v_q_ewald_lr);
132
133     /* If the table pointer == NULL the table is generated the first time =>
134        the array pointer will be saved to nbparam and the texture is bound.
135      */
136     coul_tab = nbp->coulomb_tab;
137     if (coul_tab == NULL)
138     {
139         stat = cudaMalloc((void **)&coul_tab, tabsize*sizeof(*coul_tab));
140         CU_RET_ERR(stat, "cudaMalloc failed on coul_tab");
141
142         nbp->coulomb_tab = coul_tab;
143
144 #ifdef TEXOBJ_SUPPORTED
145         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
146         if (dev_info->prop.major >= 3)
147         {
148             cudaResourceDesc rd;
149             memset(&rd, 0, sizeof(rd));
150             rd.resType                  = cudaResourceTypeLinear;
151             rd.res.linear.devPtr        = nbp->coulomb_tab;
152             rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
153             rd.res.linear.desc.x        = 32;
154             rd.res.linear.sizeInBytes   = tabsize*sizeof(*coul_tab);
155
156             cudaTextureDesc td;
157             memset(&td, 0, sizeof(td));
158             td.readMode                 = cudaReadModeElementType;
159             stat = cudaCreateTextureObject(&nbp->coulomb_tab_texobj, &rd, &td, NULL);
160             CU_RET_ERR(stat, "cudaCreateTextureObject on coulomb_tab_texobj failed");
161         }
162         else
163 #endif
164         {
165             GMX_UNUSED_VALUE(dev_info);
166             cudaChannelFormatDesc cd   = cudaCreateChannelDesc<float>();
167             stat = cudaBindTexture(NULL, &nbnxn_cuda_get_coulomb_tab_texref(),
168                                    coul_tab, &cd, tabsize*sizeof(*coul_tab));
169             CU_RET_ERR(stat, "cudaBindTexture on coulomb_tab_texref failed");
170         }
171     }
172
173     cu_copy_H2D(coul_tab, ftmp, tabsize*sizeof(*coul_tab));
174
175     nbp->coulomb_tab_size     = tabsize;
176     nbp->coulomb_tab_scale    = tabscale;
177
178     pfree(ftmp);
179 }
180
181
182 /*! Initializes the atomdata structure first time, it only gets filled at
183     pair-search. */
184 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
185 {
186     cudaError_t stat;
187
188     ad->ntypes  = ntypes;
189     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
190     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
191     ad->bShiftVecUploaded = false;
192
193     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
194     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
195
196     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
197     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
198     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
199     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
200
201     /* initialize to NULL poiters to data that is not allocated here and will
202        need reallocation in nbnxn_cuda_init_atomdata */
203     ad->xq = NULL;
204     ad->f  = NULL;
205
206     /* size -1 indicates that the respective array hasn't been initialized yet */
207     ad->natoms = -1;
208     ad->nalloc = -1;
209 }
210
211 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
212     earlier GPUs, single or twin cut-off. */
213 static int pick_ewald_kernel_type(bool                   bTwinCut,
214                                   const cuda_dev_info_t *dev_info)
215 {
216     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
217     int  kernel_type;
218
219     /* Benchmarking/development environment variables to force the use of
220        analytical or tabulated Ewald kernel. */
221     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
222     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
223
224     if (bForceAnalyticalEwald && bForceTabulatedEwald)
225     {
226         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
227                    "requested through environment variables.");
228     }
229
230     /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
231     if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
232     {
233         bUseAnalyticalEwald = true;
234
235         if (debug)
236         {
237             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
238         }
239     }
240     else
241     {
242         bUseAnalyticalEwald = false;
243
244         if (debug)
245         {
246             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
247         }
248     }
249
250     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
251        forces it (use it for debugging/benchmarking only). */
252     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
253     {
254         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
255     }
256     else
257     {
258         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
259     }
260
261     return kernel_type;
262 }
263
264 /*! Copies all parameters related to the cut-off from ic to nbp */
265 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
266                                   const interaction_const_t *ic)
267 {
268     nbp->ewald_beta       = ic->ewaldcoeff_q;
269     nbp->sh_ewald         = ic->sh_ewald;
270     nbp->epsfac           = ic->epsfac;
271     nbp->two_k_rf         = 2.0 * ic->k_rf;
272     nbp->c_rf             = ic->c_rf;
273     nbp->rvdw_sq          = ic->rvdw * ic->rvdw;
274     nbp->rcoulomb_sq      = ic->rcoulomb * ic->rcoulomb;
275     nbp->rlist_sq         = ic->rlist * ic->rlist;
276
277     nbp->sh_lj_ewald      = ic->sh_lj_ewald;
278     nbp->ewaldcoeff_lj    = ic->ewaldcoeff_lj;
279
280     nbp->rvdw_switch      = ic->rvdw_switch;
281     nbp->dispersion_shift = ic->dispersion_shift;
282     nbp->repulsion_shift  = ic->repulsion_shift;
283     nbp->vdw_switch       = ic->vdw_switch;
284 }
285
286 /*! Initializes the nonbonded parameter data structure. */
287 static void init_nbparam(cu_nbparam_t              *nbp,
288                          const interaction_const_t *ic,
289                          const nbnxn_atomdata_t    *nbat,
290                          const cuda_dev_info_t     *dev_info)
291 {
292     cudaError_t stat;
293     int         ntypes, nnbfp, nnbfp_comb;
294
295     ntypes  = nbat->ntype;
296
297     set_cutoff_parameters(nbp, ic);
298
299     if (ic->vdwtype == evdwCUT)
300     {
301         switch (ic->vdw_modifier)
302         {
303             case eintmodNONE:
304             case eintmodPOTSHIFT:
305                 nbp->vdwtype = evdwCuCUT;
306                 break;
307             case eintmodFORCESWITCH:
308                 nbp->vdwtype = evdwCuFSWITCH;
309                 break;
310             case eintmodPOTSWITCH:
311                 nbp->vdwtype = evdwCuPSWITCH;
312                 break;
313             default:
314                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
315                 break;
316         }
317     }
318     else if (ic->vdwtype == evdwPME)
319     {
320         if (ic->ljpme_comb_rule == ljcrGEOM)
321         {
322             assert(nbat->comb_rule == ljcrGEOM);
323             nbp->vdwtype = evdwCuEWALDGEOM;
324         }
325         else
326         {
327             assert(nbat->comb_rule == ljcrLB);
328             nbp->vdwtype = evdwCuEWALDLB;
329         }
330     }
331     else
332     {
333         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
334     }
335
336     if (ic->eeltype == eelCUT)
337     {
338         nbp->eeltype = eelCuCUT;
339     }
340     else if (EEL_RF(ic->eeltype))
341     {
342         nbp->eeltype = eelCuRF;
343     }
344     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
345     {
346         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
347         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
348     }
349     else
350     {
351         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
352         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
353     }
354
355     /* generate table for PME */
356     nbp->coulomb_tab = NULL;
357     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
358     {
359         init_ewald_coulomb_force_table(nbp, dev_info);
360     }
361
362     nnbfp      = 2*ntypes*ntypes;
363     nnbfp_comb = 2*ntypes;
364
365     stat  = cudaMalloc((void **)&nbp->nbfp, nnbfp*sizeof(*nbp->nbfp));
366     CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp");
367     cu_copy_H2D(nbp->nbfp, nbat->nbfp, nnbfp*sizeof(*nbp->nbfp));
368
369
370     if (ic->vdwtype == evdwPME)
371     {
372         stat  = cudaMalloc((void **)&nbp->nbfp_comb, nnbfp_comb*sizeof(*nbp->nbfp_comb));
373         CU_RET_ERR(stat, "cudaMalloc failed on nbp->nbfp_comb");
374         cu_copy_H2D(nbp->nbfp_comb, nbat->nbfp_comb, nnbfp_comb*sizeof(*nbp->nbfp_comb));
375     }
376
377 #ifdef TEXOBJ_SUPPORTED
378     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
379     if (dev_info->prop.major >= 3)
380     {
381         cudaResourceDesc rd;
382         cudaTextureDesc  td;
383
384         memset(&rd, 0, sizeof(rd));
385         rd.resType                  = cudaResourceTypeLinear;
386         rd.res.linear.devPtr        = nbp->nbfp;
387         rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
388         rd.res.linear.desc.x        = 32;
389         rd.res.linear.sizeInBytes   = nnbfp*sizeof(*nbp->nbfp);
390
391         memset(&td, 0, sizeof(td));
392         td.readMode                 = cudaReadModeElementType;
393         stat = cudaCreateTextureObject(&nbp->nbfp_texobj, &rd, &td, NULL);
394         CU_RET_ERR(stat, "cudaCreateTextureObject on nbfp_texobj failed");
395
396         if (ic->vdwtype == evdwPME)
397         {
398             memset(&rd, 0, sizeof(rd));
399             rd.resType                  = cudaResourceTypeLinear;
400             rd.res.linear.devPtr        = nbp->nbfp_comb;
401             rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
402             rd.res.linear.desc.x        = 32;
403             rd.res.linear.sizeInBytes   = nnbfp_comb*sizeof(*nbp->nbfp_comb);
404
405             memset(&td, 0, sizeof(td));
406             td.readMode = cudaReadModeElementType;
407             stat        = cudaCreateTextureObject(&nbp->nbfp_comb_texobj, &rd, &td, NULL);
408             CU_RET_ERR(stat, "cudaCreateTextureObject on nbfp_comb_texobj failed");
409         }
410     }
411     else
412 #endif
413     {
414         cudaChannelFormatDesc cd = cudaCreateChannelDesc<float>();
415         stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_texref(),
416                                nbp->nbfp, &cd, nnbfp*sizeof(*nbp->nbfp));
417         CU_RET_ERR(stat, "cudaBindTexture on nbfp_texref failed");
418
419         if (ic->vdwtype == evdwPME)
420         {
421             stat = cudaBindTexture(NULL, &nbnxn_cuda_get_nbfp_comb_texref(),
422                                    nbp->nbfp_comb, &cd, nnbfp_comb*sizeof(*nbp->nbfp_comb));
423             CU_RET_ERR(stat, "cudaBindTexture on nbfp_comb_texref failed");
424         }
425     }
426 }
427
428 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
429  *  electrostatic type switching to twin cut-off (or back) if needed. */
430 void nbnxn_cuda_pme_loadbal_update_param(nbnxn_cuda_ptr_t           cu_nb,
431                                          const interaction_const_t *ic)
432 {
433     cu_nbparam_t *nbp = cu_nb->nbparam;
434
435     set_cutoff_parameters(nbp, ic);
436
437     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
438                                                  cu_nb->dev_info);
439
440     init_ewald_coulomb_force_table(cu_nb->nbparam, cu_nb->dev_info);
441 }
442
443 /*! Initializes the pair list data structure. */
444 static void init_plist(cu_plist_t *pl)
445 {
446     /* initialize to NULL pointers to data that is not allocated here and will
447        need reallocation in nbnxn_cuda_init_pairlist */
448     pl->sci     = NULL;
449     pl->cj4     = NULL;
450     pl->excl    = NULL;
451
452     /* size -1 indicates that the respective array hasn't been initialized yet */
453     pl->na_c        = -1;
454     pl->nsci        = -1;
455     pl->sci_nalloc  = -1;
456     pl->ncj4        = -1;
457     pl->cj4_nalloc  = -1;
458     pl->nexcl       = -1;
459     pl->excl_nalloc = -1;
460     pl->bDoPrune    = false;
461 }
462
463 /*! Initializes the timer data structure. */
464 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
465 {
466     cudaError_t stat;
467     int         eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
468
469     stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
470     CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
471     stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
472     CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
473
474     /* The non-local counters/stream (second in the array) are needed only with DD. */
475     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
476     {
477         stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
478         CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
479         stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
480         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
481
482
483         stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
484         CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
485         stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
486         CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
487
488         stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
489         CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
490         stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
491         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
492
493         stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
494         CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
495         stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
496         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
497     }
498 }
499
500 /*! Initializes the timings data structure. */
501 static void init_timings(wallclock_gpu_t *t)
502 {
503     int i, j;
504
505     t->nb_h2d_t = 0.0;
506     t->nb_d2h_t = 0.0;
507     t->nb_c     = 0;
508     t->pl_h2d_t = 0.0;
509     t->pl_h2d_c = 0;
510     for (i = 0; i < 2; i++)
511     {
512         for (j = 0; j < 2; j++)
513         {
514             t->ktime[i][j].t = 0.0;
515             t->ktime[i][j].c = 0;
516         }
517     }
518 }
519
520 void nbnxn_cuda_init(FILE                 *fplog,
521                      nbnxn_cuda_ptr_t     *p_cu_nb,
522                      const gmx_gpu_info_t *gpu_info,
523                      const gmx_gpu_opt_t  *gpu_opt,
524                      int                   my_gpu_index,
525                      gmx_bool              bLocalAndNonlocal)
526 {
527     cudaError_t       stat;
528     nbnxn_cuda_ptr_t  nb;
529     char              sbuf[STRLEN];
530     bool              bStreamSync, bNoStreamSync, bTMPIAtomics, bX86, bOldDriver;
531     int               cuda_drv_ver;
532
533     assert(gpu_info);
534
535     if (p_cu_nb == NULL)
536     {
537         return;
538     }
539
540     snew(nb, 1);
541     snew(nb->atdat, 1);
542     snew(nb->nbparam, 1);
543     snew(nb->plist[eintLocal], 1);
544     if (bLocalAndNonlocal)
545     {
546         snew(nb->plist[eintNonlocal], 1);
547     }
548
549     nb->bUseTwoStreams = bLocalAndNonlocal;
550
551     snew(nb->timers, 1);
552     snew(nb->timings, 1);
553
554     /* init nbst */
555     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
556     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
557     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
558
559     init_plist(nb->plist[eintLocal]);
560
561     /* set device info, just point it to the right GPU among the detected ones */
562     nb->dev_info = &gpu_info->cuda_dev[get_gpu_device_id(gpu_info, gpu_opt, my_gpu_index)];
563
564     /* local/non-local GPU streams */
565     stat = cudaStreamCreate(&nb->stream[eintLocal]);
566     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
567     if (nb->bUseTwoStreams)
568     {
569         init_plist(nb->plist[eintNonlocal]);
570
571         /* CUDA stream priority available in the CUDA RT 5.5 API.
572          * Note that the device we're running on does not have to support
573          * priorities, because we are querying the priority range which in this
574          * case will be a single value.
575          */
576 #if CUDA_VERSION >= 5500
577         {
578             int highest_priority;
579             stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
580             CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
581
582             stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
583                                                 cudaStreamDefault,
584                                                 highest_priority);
585             CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
586         }
587 #else
588         stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
589         CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
590 #endif
591     }
592
593     /* init events for sychronization (timing disabled for performance reasons!) */
594     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
595     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
596     stat = cudaEventCreateWithFlags(&nb->misc_ops_done, cudaEventDisableTiming);
597     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_one failed");
598
599     /* On GPUs with ECC enabled, cudaStreamSynchronize shows a large overhead
600      * (which increases with shorter time/step) caused by a known CUDA driver bug.
601      * To work around the issue we'll use an (admittedly fragile) memory polling
602      * waiting to preserve performance. This requires support for atomic
603      * operations and only works on x86/x86_64.
604      * With polling wait event-timing also needs to be disabled.
605      *
606      * The overhead is greatly reduced in API v5.0 drivers and the improvement
607      * is independent of runtime version. Hence, with API v5.0 drivers and later
608      * we won't switch to polling.
609      *
610      * NOTE: Unfortunately, this is known to fail when GPUs are shared by (t)MPI,
611      * ranks so we will also disable it in that case.
612      */
613
614     bStreamSync    = getenv("GMX_CUDA_STREAMSYNC") != NULL;
615     bNoStreamSync  = getenv("GMX_NO_CUDA_STREAMSYNC") != NULL;
616
617 #ifdef TMPI_ATOMICS
618     bTMPIAtomics = true;
619 #else
620     bTMPIAtomics = false;
621 #endif
622
623 #ifdef GMX_TARGET_X86
624     bX86 = true;
625 #else
626     bX86 = false;
627 #endif
628
629     if (bStreamSync && bNoStreamSync)
630     {
631         gmx_fatal(FARGS, "Conflicting environment variables: both GMX_CUDA_STREAMSYNC and GMX_NO_CUDA_STREAMSYNC defined");
632     }
633
634     stat = cudaDriverGetVersion(&cuda_drv_ver);
635     CU_RET_ERR(stat, "cudaDriverGetVersion failed");
636
637     bOldDriver = (cuda_drv_ver < 5000);
638
639     if ((nb->dev_info->prop.ECCEnabled == 1) && bOldDriver)
640     {
641         /* Polling wait should be used instead of cudaStreamSynchronize only if:
642          *   - ECC is ON & driver is old (checked above),
643          *   - we're on x86/x86_64,
644          *   - atomics are available, and
645          *   - GPUs are not being shared.
646          */
647         bool bShouldUsePollSync = (bX86 && bTMPIAtomics &&
648                                    (gmx_count_gpu_dev_shared(gpu_opt) < 1));
649
650         if (bStreamSync)
651         {
652             nb->bUseStreamSync = true;
653
654             /* only warn if polling should be used */
655             if (bShouldUsePollSync)
656             {
657                 md_print_warn(fplog,
658                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, but\n"
659                               "      cudaStreamSynchronize waiting is forced by the GMX_CUDA_STREAMSYNC env. var.\n");
660             }
661         }
662         else
663         {
664             nb->bUseStreamSync = !bShouldUsePollSync;
665
666             if (bShouldUsePollSync)
667             {
668                 md_print_warn(fplog,
669                               "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0, known to\n"
670                               "      cause performance loss. Switching to the alternative polling GPU wait.\n"
671                               "      If you encounter issues, switch back to standard GPU waiting by setting\n"
672                               "      the GMX_CUDA_STREAMSYNC environment variable.\n");
673             }
674             else
675             {
676                 /* Tell the user that the ECC+old driver combination can be bad */
677                 sprintf(sbuf,
678                         "NOTE: Using a GPU with ECC enabled and CUDA driver API version <5.0.\n"
679                         "      A known bug in this driver version can cause performance loss.\n"
680                         "      However, the polling wait workaround can not be used because\n%s\n"
681                         "      Consider updating the driver or turning ECC off.",
682                         (bX86 && bTMPIAtomics) ?
683                         "      GPU(s) are being oversubscribed." :
684                         "      atomic operations are not supported by the platform/CPU+compiler.");
685                 md_print_warn(fplog, sbuf);
686             }
687         }
688     }
689     else
690     {
691         if (bNoStreamSync)
692         {
693             nb->bUseStreamSync = false;
694
695             md_print_warn(fplog,
696                           "NOTE: Polling wait for GPU synchronization requested by GMX_NO_CUDA_STREAMSYNC\n");
697         }
698         else
699         {
700             /* no/off ECC, cudaStreamSynchronize not turned off by env. var. */
701             nb->bUseStreamSync = true;
702         }
703     }
704
705     /* CUDA timing disabled as event timers don't work:
706        - with multiple streams = domain-decomposition;
707        - with the polling waiting hack (without cudaStreamSynchronize);
708        - when turned off by GMX_DISABLE_CUDA_TIMING.
709      */
710     nb->bDoTime = (!nb->bUseTwoStreams && nb->bUseStreamSync &&
711                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL));
712
713     if (nb->bDoTime)
714     {
715         init_timers(nb->timers, nb->bUseTwoStreams);
716         init_timings(nb->timings);
717     }
718
719     /* set the kernel type for the current GPU */
720     /* pick L1 cache configuration */
721     nbnxn_cuda_set_cacheconfig(nb->dev_info);
722
723     *p_cu_nb = nb;
724
725     if (debug)
726     {
727         fprintf(debug, "Initialized CUDA data structures.\n");
728     }
729 }
730
731 void nbnxn_cuda_init_const(nbnxn_cuda_ptr_t                cu_nb,
732                            const interaction_const_t      *ic,
733                            const nonbonded_verlet_group_t *nbv_group)
734 {
735     init_atomdata_first(cu_nb->atdat, nbv_group[0].nbat->ntype);
736     init_nbparam(cu_nb->nbparam, ic, nbv_group[0].nbat, cu_nb->dev_info);
737
738     /* clear energy and shift force outputs */
739     nbnxn_cuda_clear_e_fshift(cu_nb);
740 }
741
742 void nbnxn_cuda_init_pairlist(nbnxn_cuda_ptr_t        cu_nb,
743                               const nbnxn_pairlist_t *h_plist,
744                               int                     iloc)
745 {
746     char          sbuf[STRLEN];
747     cudaError_t   stat;
748     bool          bDoTime    = cu_nb->bDoTime;
749     cudaStream_t  stream     = cu_nb->stream[iloc];
750     cu_plist_t   *d_plist    = cu_nb->plist[iloc];
751
752     if (d_plist->na_c < 0)
753     {
754         d_plist->na_c = h_plist->na_ci;
755     }
756     else
757     {
758         if (d_plist->na_c != h_plist->na_ci)
759         {
760             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
761                     d_plist->na_c, h_plist->na_ci);
762             gmx_incons(sbuf);
763         }
764     }
765
766     if (bDoTime)
767     {
768         stat = cudaEventRecord(cu_nb->timers->start_pl_h2d[iloc], stream);
769         CU_RET_ERR(stat, "cudaEventRecord failed");
770     }
771
772     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
773                         &d_plist->nsci, &d_plist->sci_nalloc,
774                         h_plist->nsci,
775                         stream, true);
776
777     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
778                         &d_plist->ncj4, &d_plist->cj4_nalloc,
779                         h_plist->ncj4,
780                         stream, true);
781
782     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
783                         &d_plist->nexcl, &d_plist->excl_nalloc,
784                         h_plist->nexcl,
785                         stream, true);
786
787     if (bDoTime)
788     {
789         stat = cudaEventRecord(cu_nb->timers->stop_pl_h2d[iloc], stream);
790         CU_RET_ERR(stat, "cudaEventRecord failed");
791     }
792
793     /* need to prune the pair list during the next step */
794     d_plist->bDoPrune = true;
795 }
796
797 void nbnxn_cuda_upload_shiftvec(nbnxn_cuda_ptr_t        cu_nb,
798                                 const nbnxn_atomdata_t *nbatom)
799 {
800     cu_atomdata_t *adat  = cu_nb->atdat;
801     cudaStream_t   ls    = cu_nb->stream[eintLocal];
802
803     /* only if we have a dynamic box */
804     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
805     {
806         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
807                           SHIFTS * sizeof(*adat->shift_vec), ls);
808         adat->bShiftVecUploaded = true;
809     }
810 }
811
812 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
813 static void nbnxn_cuda_clear_f(nbnxn_cuda_ptr_t cu_nb, int natoms_clear)
814 {
815     cudaError_t    stat;
816     cu_atomdata_t *adat  = cu_nb->atdat;
817     cudaStream_t   ls    = cu_nb->stream[eintLocal];
818
819     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
820     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
821 }
822
823 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
824 static void nbnxn_cuda_clear_e_fshift(nbnxn_cuda_ptr_t cu_nb)
825 {
826     cudaError_t    stat;
827     cu_atomdata_t *adat  = cu_nb->atdat;
828     cudaStream_t   ls    = cu_nb->stream[eintLocal];
829
830     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
831     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
832     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
833     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
834     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
835     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
836 }
837
838 void nbnxn_cuda_clear_outputs(nbnxn_cuda_ptr_t cu_nb, int flags)
839 {
840     nbnxn_cuda_clear_f(cu_nb, cu_nb->atdat->natoms);
841     /* clear shift force array and energies if the outputs were
842        used in the current step */
843     if (flags & GMX_FORCE_VIRIAL)
844     {
845         nbnxn_cuda_clear_e_fshift(cu_nb);
846     }
847 }
848
849 void nbnxn_cuda_init_atomdata(nbnxn_cuda_ptr_t        cu_nb,
850                               const nbnxn_atomdata_t *nbat)
851 {
852     cudaError_t    stat;
853     int            nalloc, natoms;
854     bool           realloced;
855     bool           bDoTime   = cu_nb->bDoTime;
856     cu_timers_t   *timers    = cu_nb->timers;
857     cu_atomdata_t *d_atdat   = cu_nb->atdat;
858     cudaStream_t   ls        = cu_nb->stream[eintLocal];
859
860     natoms    = nbat->natoms;
861     realloced = false;
862
863     if (bDoTime)
864     {
865         /* time async copy */
866         stat = cudaEventRecord(timers->start_atdat, ls);
867         CU_RET_ERR(stat, "cudaEventRecord failed");
868     }
869
870     /* need to reallocate if we have to copy more atoms than the amount of space
871        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
872     if (natoms > d_atdat->nalloc)
873     {
874         nalloc = over_alloc_small(natoms);
875
876         /* free up first if the arrays have already been initialized */
877         if (d_atdat->nalloc != -1)
878         {
879             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
880             cu_free_buffered(d_atdat->xq);
881             cu_free_buffered(d_atdat->atom_types);
882         }
883
884         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
885         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
886         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
887         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
888
889         stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
890         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
891
892         d_atdat->nalloc = nalloc;
893         realloced       = true;
894     }
895
896     d_atdat->natoms       = natoms;
897     d_atdat->natoms_local = nbat->natoms_local;
898
899     /* need to clear GPU f output if realloc happened */
900     if (realloced)
901     {
902         nbnxn_cuda_clear_f(cu_nb, nalloc);
903     }
904
905     cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
906                       natoms*sizeof(*d_atdat->atom_types), ls);
907
908     if (bDoTime)
909     {
910         stat = cudaEventRecord(timers->stop_atdat, ls);
911         CU_RET_ERR(stat, "cudaEventRecord failed");
912     }
913 }
914
915 void nbnxn_cuda_free(nbnxn_cuda_ptr_t cu_nb)
916 {
917     cudaError_t      stat;
918     cu_atomdata_t   *atdat;
919     cu_nbparam_t    *nbparam;
920     cu_plist_t      *plist, *plist_nl;
921     cu_timers_t     *timers;
922
923     if (cu_nb == NULL)
924     {
925         return;
926     }
927
928     atdat       = cu_nb->atdat;
929     nbparam     = cu_nb->nbparam;
930     plist       = cu_nb->plist[eintLocal];
931     plist_nl    = cu_nb->plist[eintNonlocal];
932     timers      = cu_nb->timers;
933
934     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
935     {
936
937 #ifdef TEXOBJ_SUPPORTED
938         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
939         if (cu_nb->dev_info->prop.major >= 3)
940         {
941             stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
942             CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
943         }
944         else
945 #endif
946         {
947             stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
948             CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
949         }
950         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
951     }
952
953     stat = cudaEventDestroy(cu_nb->nonlocal_done);
954     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
955     stat = cudaEventDestroy(cu_nb->misc_ops_done);
956     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_done");
957
958     if (cu_nb->bDoTime)
959     {
960         stat = cudaEventDestroy(timers->start_atdat);
961         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
962         stat = cudaEventDestroy(timers->stop_atdat);
963         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
964
965         /* The non-local counters/stream (second in the array) are needed only with DD. */
966         for (int i = 0; i <= (cu_nb->bUseTwoStreams ? 1 : 0); i++)
967         {
968             stat = cudaEventDestroy(timers->start_nb_k[i]);
969             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
970             stat = cudaEventDestroy(timers->stop_nb_k[i]);
971             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
972
973             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
974             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
975             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
976             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
977
978             stat = cudaStreamDestroy(cu_nb->stream[i]);
979             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
980
981             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
982             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
983             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
984             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
985
986             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
987             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
988             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
989             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
990         }
991     }
992
993 #ifdef TEXOBJ_SUPPORTED
994     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
995     if (cu_nb->dev_info->prop.major >= 3)
996     {
997         stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
998         CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
999     }
1000     else
1001 #endif
1002     {
1003         stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
1004         CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
1005     }
1006     cu_free_buffered(nbparam->nbfp);
1007
1008     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
1009     {
1010 #ifdef TEXOBJ_SUPPORTED
1011         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
1012         if (cu_nb->dev_info->prop.major >= 3)
1013         {
1014             stat = cudaDestroyTextureObject(nbparam->nbfp_comb_texobj);
1015             CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_comb_texobj failed");
1016         }
1017         else
1018 #endif
1019         {
1020             stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_comb_texref());
1021             CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_comb_texref failed");
1022         }
1023         cu_free_buffered(nbparam->nbfp_comb);
1024     }
1025
1026     stat = cudaFree(atdat->shift_vec);
1027     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
1028     stat = cudaFree(atdat->fshift);
1029     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
1030
1031     stat = cudaFree(atdat->e_lj);
1032     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
1033     stat = cudaFree(atdat->e_el);
1034     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
1035
1036     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
1037     cu_free_buffered(atdat->xq);
1038     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
1039
1040     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
1041     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
1042     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
1043     if (cu_nb->bUseTwoStreams)
1044     {
1045         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
1046         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
1047         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
1048     }
1049
1050     sfree(atdat);
1051     sfree(nbparam);
1052     sfree(plist);
1053     if (cu_nb->bUseTwoStreams)
1054     {
1055         sfree(plist_nl);
1056     }
1057     sfree(timers);
1058     sfree(cu_nb->timings);
1059     sfree(cu_nb);
1060
1061     if (debug)
1062     {
1063         fprintf(debug, "Cleaned up CUDA data structures.\n");
1064     }
1065 }
1066
1067 void cu_synchstream_atdat(nbnxn_cuda_ptr_t cu_nb, int iloc)
1068 {
1069     cudaError_t  stat;
1070     cudaStream_t stream = cu_nb->stream[iloc];
1071
1072     stat = cudaStreamWaitEvent(stream, cu_nb->timers->stop_atdat, 0);
1073     CU_RET_ERR(stat, "cudaStreamWaitEvent failed");
1074 }
1075
1076 wallclock_gpu_t * nbnxn_cuda_get_timings(nbnxn_cuda_ptr_t cu_nb)
1077 {
1078     return (cu_nb != NULL && cu_nb->bDoTime) ? cu_nb->timings : NULL;
1079 }
1080
1081 void nbnxn_cuda_reset_timings(nbnxn_cuda_ptr_t cu_nb)
1082 {
1083     if (cu_nb->bDoTime)
1084     {
1085         init_timings(cu_nb->timings);
1086     }
1087 }
1088
1089 int nbnxn_cuda_min_ci_balanced(nbnxn_cuda_ptr_t cu_nb)
1090 {
1091     return cu_nb != NULL ?
1092            gpu_min_ci_balanced_factor*cu_nb->dev_info->prop.multiProcessorCount : 0;
1093
1094 }
1095
1096 gmx_bool nbnxn_cuda_is_kernel_ewald_analytical(const nbnxn_cuda_ptr_t cu_nb)
1097 {
1098     return ((cu_nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1099             (cu_nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
1100 }