52e255de560dc788b5f363beb5f303b8e9e93251
[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,2015,2016,2017, 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 /*! \file
36  *  \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
37  *
38  *  \author Szilard Pall <pall.szilard@gmail.com>
39  */
40 #include "gmxpre.h"
41
42 #include <assert.h>
43 #include <stdarg.h>
44 #include <stdio.h>
45 #include <stdlib.h>
46
47 #include "gromacs/gpu_utils/cudautils.cuh"
48 #include "gromacs/gpu_utils/gpu_utils.h"
49 #include "gromacs/gpu_utils/pmalloc_cuda.h"
50 #include "gromacs/hardware/gpu_hw_info.h"
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/mdlib/force_flags.h"
53 #include "gromacs/mdlib/nb_verlet.h"
54 #include "gromacs/mdlib/nbnxn_consts.h"
55 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/timing/gpu_timing.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/real.h"
64 #include "gromacs/utility/smalloc.h"
65
66 #include "nbnxn_cuda_types.h"
67
68 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
69
70 /* This is a heuristically determined parameter for the Fermi, Kepler
71  * and Maxwell architectures for the minimum size of ci lists by multiplying
72  * this constant with the # of multiprocessors on the current device.
73  * Since the maximum number of blocks per multiprocessor is 16, the ideal
74  * count for small systems is 32 or 48 blocks per multiprocessor. Because
75  * there is a bit of fluctuations in the generated block counts, we use
76  * a target of 44 instead of the ideal value of 48.
77  */
78 static unsigned int gpu_min_ci_balanced_factor = 44;
79
80 /* Functions from nbnxn_cuda.cu */
81 extern void nbnxn_cuda_set_cacheconfig(const gmx_device_info_t *devinfo);
82 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref();
83 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref();
84 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref();
85
86
87 /* Fw. decl. */
88 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
89
90 /* Fw. decl, */
91 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
92                                           const gmx_device_info_t *dev_info);
93
94
95 /*! \brief Return whether texture objects are used on this device.
96  *
97  * \param[in]   pointer to the GPU device info structure to inspect for texture objects support
98  * \return      true if texture objects are used on this device
99  */
100 static bool use_texobj(const gmx_device_info_t *dev_info)
101 {
102     assert(!c_disableCudaTextures);
103     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
104     return (dev_info->prop.major >= 3);
105 }
106
107 /*! \brief Return whether combination rules are used.
108  *
109  * \param[in]   pointer to nonbonded paramter struct
110  * \return      true if combination rules are used in this run, false otherwise
111  */
112 static inline bool useLjCombRule(const cu_nbparam_t  *nbparam)
113 {
114     return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
115             nbparam->vdwtype == evdwCuCUTCOMBLB);
116 }
117
118 /*! \brief Set up float texture object.
119  *
120  * Set up texture object for float data and bind it to the device memory
121  * \p devPtr points to.
122  *
123  * \param[out] texObj   texture object to initialize
124  * \param[in]  devPtr   pointer to device global memory to bind \p texObj to
125  * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
126  */
127 static void setup1DFloatTexture(cudaTextureObject_t &texObj,
128                                 void                *devPtr,
129                                 size_t               sizeInBytes)
130 {
131     assert(!c_disableCudaTextures);
132
133     cudaError_t      stat;
134     cudaResourceDesc rd;
135     cudaTextureDesc  td;
136
137     memset(&rd, 0, sizeof(rd));
138     rd.resType                  = cudaResourceTypeLinear;
139     rd.res.linear.devPtr        = devPtr;
140     rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
141     rd.res.linear.desc.x        = 32;
142     rd.res.linear.sizeInBytes   = sizeInBytes;
143
144     memset(&td, 0, sizeof(td));
145     td.readMode                 = cudaReadModeElementType;
146     stat = cudaCreateTextureObject(&texObj, &rd, &td, NULL);
147     CU_RET_ERR(stat, "cudaCreateTextureObject failed");
148 }
149
150 /*! \brief Set up float texture reference.
151  *
152  * Set up texture object for float data and bind it to the device memory
153  * \p devPtr points to.
154  *
155  * \param[out] texObj   texture reference to initialize
156  * \param[in]  devPtr   pointer to device global memory to bind \p texObj to
157  * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
158  */
159 static void setup1DFloatTexture(const struct texture<float, 1, cudaReadModeElementType> *texRef,
160                                 const void                                              *devPtr,
161                                 size_t                                                   sizeInBytes)
162 {
163     assert(!c_disableCudaTextures);
164
165     cudaError_t           stat;
166     cudaChannelFormatDesc cd;
167
168     cd   = cudaCreateChannelDesc<float>();
169     stat = cudaBindTexture(NULL, texRef, devPtr, &cd, sizeInBytes);
170     CU_RET_ERR(stat, "cudaBindTexture failed");
171 }
172
173
174 /*! \brief Initialized the Ewald Coulomb correction GPU table.
175
176     Tabulates the Ewald Coulomb force and initializes the size/scale
177     and the table GPU array. If called with an already allocated table,
178     it just re-uploads the table.
179  */
180 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
181                                            cu_nbparam_t              *nbp,
182                                            const gmx_device_info_t   *dev_info)
183 {
184     float       *coul_tab;
185     cudaError_t  stat;
186
187     if (nbp->coulomb_tab != NULL)
188     {
189         nbnxn_cuda_free_nbparam_table(nbp, dev_info);
190     }
191
192     /* initialize table data in nbp and crete/copy into in global mem */
193     stat = cudaMalloc((void **)&coul_tab, ic->tabq_size*sizeof(*coul_tab));
194     CU_RET_ERR(stat, "cudaMalloc failed on coulumb_tab");
195     cu_copy_H2D(coul_tab, ic->tabq_coul_F, ic->tabq_size*sizeof(*coul_tab));
196
197     nbp->coulomb_tab       = coul_tab;
198     nbp->coulomb_tab_size  = ic->tabq_size;
199     nbp->coulomb_tab_scale = ic->tabq_scale;
200
201     if (!c_disableCudaTextures)
202     {
203         if (use_texobj(dev_info))
204         {
205             setup1DFloatTexture(nbp->coulomb_tab_texobj, nbp->coulomb_tab,
206                                 nbp->coulomb_tab_size*sizeof(*nbp->coulomb_tab));
207         }
208         else
209         {
210             setup1DFloatTexture(&nbnxn_cuda_get_coulomb_tab_texref(), nbp->coulomb_tab,
211                                 nbp->coulomb_tab_size*sizeof(*nbp->coulomb_tab));
212         }
213     }
214 }
215
216
217 /*! Initializes the atomdata structure first time, it only gets filled at
218     pair-search. */
219 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
220 {
221     cudaError_t stat;
222
223     ad->ntypes  = ntypes;
224     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
225     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
226     ad->bShiftVecUploaded = false;
227
228     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
229     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
230
231     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
232     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
233     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
234     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
235
236     /* initialize to NULL poiters to data that is not allocated here and will
237        need reallocation in nbnxn_cuda_init_atomdata */
238     ad->xq = NULL;
239     ad->f  = NULL;
240
241     /* size -1 indicates that the respective array hasn't been initialized yet */
242     ad->natoms = -1;
243     ad->nalloc = -1;
244 }
245
246 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
247     earlier GPUs, single or twin cut-off. */
248 static int pick_ewald_kernel_type(bool                     bTwinCut,
249                                   const gmx_device_info_t *dev_info)
250 {
251     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
252     int  kernel_type;
253
254     /* Benchmarking/development environment variables to force the use of
255        analytical or tabulated Ewald kernel. */
256     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
257     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
258
259     if (bForceAnalyticalEwald && bForceTabulatedEwald)
260     {
261         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
262                    "requested through environment variables.");
263     }
264
265     /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
266     if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
267     {
268         bUseAnalyticalEwald = true;
269
270         if (debug)
271         {
272             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
273         }
274     }
275     else
276     {
277         bUseAnalyticalEwald = false;
278
279         if (debug)
280         {
281             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
282         }
283     }
284
285     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
286        forces it (use it for debugging/benchmarking only). */
287     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
288     {
289         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
290     }
291     else
292     {
293         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
294     }
295
296     return kernel_type;
297 }
298
299 /*! Copies all parameters related to the cut-off from ic to nbp */
300 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
301                                   const interaction_const_t *ic)
302 {
303     nbp->ewald_beta       = ic->ewaldcoeff_q;
304     nbp->sh_ewald         = ic->sh_ewald;
305     nbp->epsfac           = ic->epsfac;
306     nbp->two_k_rf         = 2.0 * ic->k_rf;
307     nbp->c_rf             = ic->c_rf;
308     nbp->rvdw_sq          = ic->rvdw * ic->rvdw;
309     nbp->rcoulomb_sq      = ic->rcoulomb * ic->rcoulomb;
310     nbp->rlist_sq         = ic->rlist * ic->rlist;
311
312     nbp->sh_lj_ewald      = ic->sh_lj_ewald;
313     nbp->ewaldcoeff_lj    = ic->ewaldcoeff_lj;
314
315     nbp->rvdw_switch      = ic->rvdw_switch;
316     nbp->dispersion_shift = ic->dispersion_shift;
317     nbp->repulsion_shift  = ic->repulsion_shift;
318     nbp->vdw_switch       = ic->vdw_switch;
319 }
320
321 /*! \brief Initialize LJ parameter lookup table.
322  *
323  * Initializes device memory and copies data from host an binds
324  * a texture to allocated device memory to be used for LJ parameter
325  * lookup.
326  *
327  * \param[out] devPtr    device pointer to the memory to be allocated
328  * \param[out] texObj    texture object to be initialized
329  * \param[out] texRef    texture reference to be initialized
330  * \param[in]  hostPtr   pointer to the host memory to be uploaded to the device
331  * \param[in]  numElem   number of elements in the hostPtr
332  * \param[in]  devInfo   pointer to the info struct of the device in use
333  */
334 static void initParamLookupTable(float                    * &devPtr,
335                                  cudaTextureObject_t       &texObj,
336                                  const struct texture<float, 1, cudaReadModeElementType> *texRef,
337                                  const float               *hostPtr,
338                                  int                        numElem,
339                                  const gmx_device_info_t   *devInfo)
340 {
341     cudaError_t stat;
342
343     size_t      sizeInBytes = numElem*sizeof(*devPtr);
344
345     stat  = cudaMalloc((void **)&devPtr, sizeInBytes);
346     CU_RET_ERR(stat, "cudaMalloc failed in initParamLookupTable");
347     cu_copy_H2D(devPtr, (void *)hostPtr, sizeInBytes);
348
349     if (!c_disableCudaTextures)
350     {
351         if (use_texobj(devInfo))
352         {
353             setup1DFloatTexture(texObj, devPtr, sizeInBytes);
354         }
355         else
356         {
357             setup1DFloatTexture(texRef, devPtr, sizeInBytes);
358         }
359     }
360 }
361
362 /*! Initializes the nonbonded parameter data structure. */
363 static void init_nbparam(cu_nbparam_t              *nbp,
364                          const interaction_const_t *ic,
365                          const nbnxn_atomdata_t    *nbat,
366                          const gmx_device_info_t   *dev_info)
367 {
368     int         ntypes;
369
370     ntypes  = nbat->ntype;
371
372     set_cutoff_parameters(nbp, ic);
373
374     /* The kernel code supports LJ combination rules (geometric and LB) for
375      * all kernel types, but we only generate useful combination rule kernels.
376      * We currently only use LJ combination rule (geometric and LB) kernels
377      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
378      * with PME and 20% with RF, the other kernels speed up about half as much.
379      * For LJ force-switch the geometric rule would give 7% speed-up, but this
380      * combination is rarely used. LJ force-switch with LB rule is more common,
381      * but gives only 1% speed-up.
382      */
383     if (ic->vdwtype == evdwCUT)
384     {
385         switch (ic->vdw_modifier)
386         {
387             case eintmodNONE:
388             case eintmodPOTSHIFT:
389                 switch (nbat->comb_rule)
390                 {
391                     case ljcrNONE:
392                         nbp->vdwtype = evdwCuCUT;
393                         break;
394                     case ljcrGEOM:
395                         nbp->vdwtype = evdwCuCUTCOMBGEOM;
396                         break;
397                     case ljcrLB:
398                         nbp->vdwtype = evdwCuCUTCOMBLB;
399                         break;
400                     default:
401                         gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
402                         break;
403                 }
404                 break;
405             case eintmodFORCESWITCH:
406                 nbp->vdwtype = evdwCuFSWITCH;
407                 break;
408             case eintmodPOTSWITCH:
409                 nbp->vdwtype = evdwCuPSWITCH;
410                 break;
411             default:
412                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
413                 break;
414         }
415     }
416     else if (ic->vdwtype == evdwPME)
417     {
418         if (ic->ljpme_comb_rule == ljcrGEOM)
419         {
420             assert(nbat->comb_rule == ljcrGEOM);
421             nbp->vdwtype = evdwCuEWALDGEOM;
422         }
423         else
424         {
425             assert(nbat->comb_rule == ljcrLB);
426             nbp->vdwtype = evdwCuEWALDLB;
427         }
428     }
429     else
430     {
431         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
432     }
433
434     if (ic->eeltype == eelCUT)
435     {
436         nbp->eeltype = eelCuCUT;
437     }
438     else if (EEL_RF(ic->eeltype))
439     {
440         nbp->eeltype = eelCuRF;
441     }
442     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
443     {
444         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
445         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
446     }
447     else
448     {
449         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
450         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
451     }
452
453     /* generate table for PME */
454     nbp->coulomb_tab = NULL;
455     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
456     {
457         init_ewald_coulomb_force_table(ic, nbp, dev_info);
458     }
459
460     /* set up LJ parameter lookup table */
461     if (!useLjCombRule(nbp))
462     {
463         initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
464                              &nbnxn_cuda_get_nbfp_texref(),
465                              nbat->nbfp, 2*ntypes*ntypes, dev_info);
466     }
467
468     /* set up LJ-PME parameter lookup table */
469     if (ic->vdwtype == evdwPME)
470     {
471         initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
472                              &nbnxn_cuda_get_nbfp_comb_texref(),
473                              nbat->nbfp_comb, 2*ntypes, dev_info);
474     }
475 }
476
477 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
478  *  electrostatic type switching to twin cut-off (or back) if needed. */
479 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t    *nbv,
480                                         const interaction_const_t   *ic)
481 {
482     if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
483     {
484         return;
485     }
486     gmx_nbnxn_cuda_t *nb    = nbv->gpu_nbv;
487     cu_nbparam_t     *nbp   = nb->nbparam;
488
489     set_cutoff_parameters(nbp, ic);
490
491     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
492                                                  nb->dev_info);
493
494     init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
495 }
496
497 /*! Initializes the pair list data structure. */
498 static void init_plist(cu_plist_t *pl)
499 {
500     /* initialize to NULL pointers to data that is not allocated here and will
501        need reallocation in nbnxn_gpu_init_pairlist */
502     pl->sci     = NULL;
503     pl->cj4     = NULL;
504     pl->excl    = NULL;
505
506     /* size -1 indicates that the respective array hasn't been initialized yet */
507     pl->na_c        = -1;
508     pl->nsci        = -1;
509     pl->sci_nalloc  = -1;
510     pl->ncj4        = -1;
511     pl->cj4_nalloc  = -1;
512     pl->nexcl       = -1;
513     pl->excl_nalloc = -1;
514     pl->bDoPrune    = false;
515 }
516
517 /*! Initializes the timer data structure. */
518 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
519 {
520     cudaError_t stat;
521     int         eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
522
523     stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
524     CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
525     stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
526     CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
527
528     /* The non-local counters/stream (second in the array) are needed only with DD. */
529     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
530     {
531         stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
532         CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
533         stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
534         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
535
536
537         stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
538         CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
539         stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
540         CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
541
542         stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
543         CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
544         stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
545         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
546
547         stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
548         CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
549         stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
550         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
551     }
552 }
553
554 /*! Initializes the timings data structure. */
555 static void init_timings(gmx_wallclock_gpu_t *t)
556 {
557     int i, j;
558
559     t->nb_h2d_t = 0.0;
560     t->nb_d2h_t = 0.0;
561     t->nb_c     = 0;
562     t->pl_h2d_t = 0.0;
563     t->pl_h2d_c = 0;
564     for (i = 0; i < 2; i++)
565     {
566         for (j = 0; j < 2; j++)
567         {
568             t->ktime[i][j].t = 0.0;
569             t->ktime[i][j].c = 0;
570         }
571     }
572 }
573
574 /*! Initializes simulation constant data. */
575 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t               *nb,
576                                   const interaction_const_t      *ic,
577                                   const nonbonded_verlet_group_t *nbv_group)
578 {
579     init_atomdata_first(nb->atdat, nbv_group[0].nbat->ntype);
580     init_nbparam(nb->nbparam, ic, nbv_group[0].nbat, nb->dev_info);
581
582     /* clear energy and shift force outputs */
583     nbnxn_cuda_clear_e_fshift(nb);
584 }
585
586 void nbnxn_gpu_init(gmx_nbnxn_cuda_t         **p_nb,
587                     const gmx_device_info_t   *deviceInfo,
588                     const interaction_const_t *ic,
589                     nonbonded_verlet_group_t  *nbv_grp,
590                     int                        /*rank*/,
591                     gmx_bool                   bLocalAndNonlocal)
592 {
593     cudaError_t       stat;
594     gmx_nbnxn_cuda_t *nb;
595
596     if (p_nb == NULL)
597     {
598         return;
599     }
600
601     snew(nb, 1);
602     snew(nb->atdat, 1);
603     snew(nb->nbparam, 1);
604     snew(nb->plist[eintLocal], 1);
605     if (bLocalAndNonlocal)
606     {
607         snew(nb->plist[eintNonlocal], 1);
608     }
609
610     nb->bUseTwoStreams = bLocalAndNonlocal;
611
612     snew(nb->timers, 1);
613     snew(nb->timings, 1);
614
615     /* init nbst */
616     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
617     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
618     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
619
620     init_plist(nb->plist[eintLocal]);
621
622     /* set device info, just point it to the right GPU among the detected ones */
623     nb->dev_info = deviceInfo;
624
625     /* local/non-local GPU streams */
626     stat = cudaStreamCreate(&nb->stream[eintLocal]);
627     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
628     if (nb->bUseTwoStreams)
629     {
630         init_plist(nb->plist[eintNonlocal]);
631
632         /* Note that the device we're running on does not have to support
633          * priorities, because we are querying the priority range which in this
634          * case will be a single value.
635          */
636         int highest_priority;
637         stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
638         CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
639
640         stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
641                                             cudaStreamDefault,
642                                             highest_priority);
643         CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
644     }
645
646     /* init events for sychronization (timing disabled for performance reasons!) */
647     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
648     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
649     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
650     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
651
652     /* CUDA timing disabled as event timers don't work:
653        - with multiple streams = domain-decomposition;
654        - when turned off by GMX_DISABLE_CUDA_TIMING/GMX_DISABLE_GPU_TIMING.
655      */
656     nb->bDoTime = (!nb->bUseTwoStreams &&
657                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL) &&
658                    (getenv("GMX_DISABLE_GPU_TIMING") == NULL));
659
660     if (nb->bDoTime)
661     {
662         init_timers(nb->timers, nb->bUseTwoStreams);
663         init_timings(nb->timings);
664     }
665
666     /* set the kernel type for the current GPU */
667     /* pick L1 cache configuration */
668     nbnxn_cuda_set_cacheconfig(nb->dev_info);
669
670     nbnxn_cuda_init_const(nb, ic, nbv_grp);
671
672     *p_nb = nb;
673
674     if (debug)
675     {
676         fprintf(debug, "Initialized CUDA data structures.\n");
677     }
678 }
679
680 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t       *nb,
681                              const nbnxn_pairlist_t *h_plist,
682                              int                     iloc)
683 {
684     char          sbuf[STRLEN];
685     cudaError_t   stat;
686     bool          bDoTime    = nb->bDoTime;
687     cudaStream_t  stream     = nb->stream[iloc];
688     cu_plist_t   *d_plist    = nb->plist[iloc];
689
690     if (d_plist->na_c < 0)
691     {
692         d_plist->na_c = h_plist->na_ci;
693     }
694     else
695     {
696         if (d_plist->na_c != h_plist->na_ci)
697         {
698             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
699                     d_plist->na_c, h_plist->na_ci);
700             gmx_incons(sbuf);
701         }
702     }
703
704     if (bDoTime)
705     {
706         stat = cudaEventRecord(nb->timers->start_pl_h2d[iloc], stream);
707         CU_RET_ERR(stat, "cudaEventRecord failed");
708     }
709
710     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
711                         &d_plist->nsci, &d_plist->sci_nalloc,
712                         h_plist->nsci,
713                         stream, true);
714
715     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
716                         &d_plist->ncj4, &d_plist->cj4_nalloc,
717                         h_plist->ncj4,
718                         stream, true);
719
720     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
721                         &d_plist->nexcl, &d_plist->excl_nalloc,
722                         h_plist->nexcl,
723                         stream, true);
724
725     if (bDoTime)
726     {
727         stat = cudaEventRecord(nb->timers->stop_pl_h2d[iloc], stream);
728         CU_RET_ERR(stat, "cudaEventRecord failed");
729     }
730
731     /* need to prune the pair list during the next step */
732     d_plist->bDoPrune = true;
733 }
734
735 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t       *nb,
736                                const nbnxn_atomdata_t *nbatom)
737 {
738     cu_atomdata_t *adat  = nb->atdat;
739     cudaStream_t   ls    = nb->stream[eintLocal];
740
741     /* only if we have a dynamic box */
742     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
743     {
744         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
745                           SHIFTS * sizeof(*adat->shift_vec), ls);
746         adat->bShiftVecUploaded = true;
747     }
748 }
749
750 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
751 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
752 {
753     cudaError_t    stat;
754     cu_atomdata_t *adat  = nb->atdat;
755     cudaStream_t   ls    = nb->stream[eintLocal];
756
757     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
758     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
759 }
760
761 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
762 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
763 {
764     cudaError_t    stat;
765     cu_atomdata_t *adat  = nb->atdat;
766     cudaStream_t   ls    = nb->stream[eintLocal];
767
768     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
769     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
770     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
771     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
772     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
773     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
774 }
775
776 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
777 {
778     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
779     /* clear shift force array and energies if the outputs were
780        used in the current step */
781     if (flags & GMX_FORCE_VIRIAL)
782     {
783         nbnxn_cuda_clear_e_fshift(nb);
784     }
785 }
786
787 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t              *nb,
788                              const struct nbnxn_atomdata_t *nbat)
789 {
790     cudaError_t    stat;
791     int            nalloc, natoms;
792     bool           realloced;
793     bool           bDoTime   = nb->bDoTime;
794     cu_timers_t   *timers    = nb->timers;
795     cu_atomdata_t *d_atdat   = nb->atdat;
796     cudaStream_t   ls        = nb->stream[eintLocal];
797
798     natoms    = nbat->natoms;
799     realloced = false;
800
801     if (bDoTime)
802     {
803         /* time async copy */
804         stat = cudaEventRecord(timers->start_atdat, ls);
805         CU_RET_ERR(stat, "cudaEventRecord failed");
806     }
807
808     /* need to reallocate if we have to copy more atoms than the amount of space
809        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
810     if (natoms > d_atdat->nalloc)
811     {
812         nalloc = over_alloc_small(natoms);
813
814         /* free up first if the arrays have already been initialized */
815         if (d_atdat->nalloc != -1)
816         {
817             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
818             cu_free_buffered(d_atdat->xq);
819             cu_free_buffered(d_atdat->atom_types);
820             cu_free_buffered(d_atdat->lj_comb);
821         }
822
823         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
824         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
825         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
826         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
827         if (useLjCombRule(nb->nbparam))
828         {
829             stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
830             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
831         }
832         else
833         {
834             stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
835             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
836         }
837
838         d_atdat->nalloc = nalloc;
839         realloced       = true;
840     }
841
842     d_atdat->natoms       = natoms;
843     d_atdat->natoms_local = nbat->natoms_local;
844
845     /* need to clear GPU f output if realloc happened */
846     if (realloced)
847     {
848         nbnxn_cuda_clear_f(nb, nalloc);
849     }
850
851     if (useLjCombRule(nb->nbparam))
852     {
853         cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
854                           natoms*sizeof(*d_atdat->lj_comb), ls);
855     }
856     else
857     {
858         cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
859                           natoms*sizeof(*d_atdat->atom_types), ls);
860     }
861
862     if (bDoTime)
863     {
864         stat = cudaEventRecord(timers->stop_atdat, ls);
865         CU_RET_ERR(stat, "cudaEventRecord failed");
866     }
867 }
868
869 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
870                                           const gmx_device_info_t *dev_info)
871 {
872     cudaError_t stat;
873
874     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
875     {
876         if (!c_disableCudaTextures)
877         {
878             /* Only device CC >= 3.0 (Kepler and later) support texture objects */
879             if (use_texobj(dev_info))
880             {
881                 stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
882                 CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
883             }
884             else
885             {
886                 GMX_UNUSED_VALUE(dev_info);
887                 stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
888                 CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
889             }
890         }
891         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
892     }
893 }
894
895 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
896 {
897     cudaError_t      stat;
898     cu_atomdata_t   *atdat;
899     cu_nbparam_t    *nbparam;
900     cu_plist_t      *plist, *plist_nl;
901     cu_timers_t     *timers;
902
903     if (nb == NULL)
904     {
905         return;
906     }
907
908     atdat       = nb->atdat;
909     nbparam     = nb->nbparam;
910     plist       = nb->plist[eintLocal];
911     plist_nl    = nb->plist[eintNonlocal];
912     timers      = nb->timers;
913
914     nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
915
916     stat = cudaEventDestroy(nb->nonlocal_done);
917     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
918     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
919     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
920
921     if (nb->bDoTime)
922     {
923         stat = cudaEventDestroy(timers->start_atdat);
924         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
925         stat = cudaEventDestroy(timers->stop_atdat);
926         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
927
928         /* The non-local counters/stream (second in the array) are needed only with DD. */
929         for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
930         {
931             stat = cudaEventDestroy(timers->start_nb_k[i]);
932             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
933             stat = cudaEventDestroy(timers->stop_nb_k[i]);
934             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
935
936             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
937             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
938             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
939             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
940
941             stat = cudaStreamDestroy(nb->stream[i]);
942             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
943
944             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
945             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
946             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
947             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
948
949             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
950             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
951             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
952             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
953         }
954     }
955
956     if (!useLjCombRule(nb->nbparam))
957     {
958         if (!c_disableCudaTextures)
959         {
960             /* Only device CC >= 3.0 (Kepler and later) support texture objects */
961             if (use_texobj(nb->dev_info))
962             {
963                 stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
964                 CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
965             }
966             else
967             {
968                 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
969                 CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
970             }
971         }
972         cu_free_buffered(nbparam->nbfp);
973     }
974
975     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
976     {
977         if (!c_disableCudaTextures)
978         {
979             /* Only device CC >= 3.0 (Kepler and later) support texture objects */
980             if (use_texobj(nb->dev_info))
981             {
982                 stat = cudaDestroyTextureObject(nbparam->nbfp_comb_texobj);
983                 CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_comb_texobj failed");
984             }
985             else
986             {
987                 stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_comb_texref());
988                 CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_comb_texref failed");
989             }
990         }
991         cu_free_buffered(nbparam->nbfp_comb);
992     }
993
994     stat = cudaFree(atdat->shift_vec);
995     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
996     stat = cudaFree(atdat->fshift);
997     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
998
999     stat = cudaFree(atdat->e_lj);
1000     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
1001     stat = cudaFree(atdat->e_el);
1002     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
1003
1004     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
1005     cu_free_buffered(atdat->xq);
1006     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
1007     cu_free_buffered(atdat->lj_comb);
1008
1009     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
1010     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
1011     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
1012     if (nb->bUseTwoStreams)
1013     {
1014         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
1015         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
1016         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
1017     }
1018
1019     sfree(atdat);
1020     sfree(nbparam);
1021     sfree(plist);
1022     if (nb->bUseTwoStreams)
1023     {
1024         sfree(plist_nl);
1025     }
1026     sfree(timers);
1027     sfree(nb->timings);
1028     sfree(nb);
1029
1030     if (debug)
1031     {
1032         fprintf(debug, "Cleaned up CUDA data structures.\n");
1033     }
1034 }
1035
1036 gmx_wallclock_gpu_t * nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
1037 {
1038     return (nb != NULL && nb->bDoTime) ? nb->timings : NULL;
1039 }
1040
1041 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
1042 {
1043     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
1044     {
1045         init_timings(nbv->gpu_nbv->timings);
1046     }
1047 }
1048
1049 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
1050 {
1051     return nb != NULL ?
1052            gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
1053
1054 }
1055
1056 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
1057 {
1058     return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1059             (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
1060 }