ea9b3b8fe4e616059cfe05c793f7c6a8167d9f84
[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,2018, 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_common_utils.h"
56 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
57 #include "gromacs/mdtypes/interaction_const.h"
58 #include "gromacs/mdtypes/md_enums.h"
59 #include "gromacs/pbcutil/ishift.h"
60 #include "gromacs/timing/gpu_timing.h"
61 #include "gromacs/utility/basedefinitions.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/fatalerror.h"
64 #include "gromacs/utility/real.h"
65 #include "gromacs/utility/smalloc.h"
66
67 #include "nbnxn_cuda.h"
68 #include "nbnxn_cuda_types.h"
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 /* Fw. decl. */
81 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
82
83 /* Fw. decl, */
84 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
85                                           const gmx_device_info_t *dev_info);
86
87 /*! \brief Return whether combination rules are used.
88  *
89  * \param[in]   pointer to nonbonded paramter struct
90  * \return      true if combination rules are used in this run, false otherwise
91  */
92 static inline bool useLjCombRule(const cu_nbparam_t  *nbparam)
93 {
94     return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
95             nbparam->vdwtype == evdwCuCUTCOMBLB);
96 }
97
98 /*! \brief Initialized the Ewald Coulomb correction GPU table.
99
100     Tabulates the Ewald Coulomb force and initializes the size/scale
101     and the table GPU array. If called with an already allocated table,
102     it just re-uploads the table.
103  */
104 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
105                                            cu_nbparam_t              *nbp,
106                                            const gmx_device_info_t   *dev_info)
107 {
108     if (nbp->coulomb_tab != NULL)
109     {
110         nbnxn_cuda_free_nbparam_table(nbp, dev_info);
111     }
112
113     nbp->coulomb_tab_scale = ic->tabq_scale;
114     initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
115                          ic->tabq_coul_F, ic->tabq_size, dev_info);
116 }
117
118
119 /*! Initializes the atomdata structure first time, it only gets filled at
120     pair-search. */
121 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
122 {
123     cudaError_t stat;
124
125     ad->ntypes  = ntypes;
126     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
127     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
128     ad->bShiftVecUploaded = false;
129
130     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
131     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
132
133     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
134     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
135     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
136     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
137
138     /* initialize to NULL poiters to data that is not allocated here and will
139        need reallocation in nbnxn_cuda_init_atomdata */
140     ad->xq = NULL;
141     ad->f  = NULL;
142
143     /* size -1 indicates that the respective array hasn't been initialized yet */
144     ad->natoms = -1;
145     ad->nalloc = -1;
146 }
147
148 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
149     earlier GPUs, single or twin cut-off. */
150 static int pick_ewald_kernel_type(bool                     bTwinCut,
151                                   const gmx_device_info_t *dev_info)
152 {
153     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
154     int  kernel_type;
155
156     /* Benchmarking/development environment variables to force the use of
157        analytical or tabulated Ewald kernel. */
158     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
159     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
160
161     if (bForceAnalyticalEwald && bForceTabulatedEwald)
162     {
163         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
164                    "requested through environment variables.");
165     }
166
167     /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
168     if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
169     {
170         bUseAnalyticalEwald = true;
171
172         if (debug)
173         {
174             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
175         }
176     }
177     else
178     {
179         bUseAnalyticalEwald = false;
180
181         if (debug)
182         {
183             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
184         }
185     }
186
187     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
188        forces it (use it for debugging/benchmarking only). */
189     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
190     {
191         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
192     }
193     else
194     {
195         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
196     }
197
198     return kernel_type;
199 }
200
201 /*! Copies all parameters related to the cut-off from ic to nbp */
202 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
203                                   const interaction_const_t *ic,
204                                   const NbnxnListParameters *listParams)
205 {
206     nbp->ewald_beta        = ic->ewaldcoeff_q;
207     nbp->sh_ewald          = ic->sh_ewald;
208     nbp->epsfac            = ic->epsfac;
209     nbp->two_k_rf          = 2.0 * ic->k_rf;
210     nbp->c_rf              = ic->c_rf;
211     nbp->rvdw_sq           = ic->rvdw * ic->rvdw;
212     nbp->rcoulomb_sq       = ic->rcoulomb * ic->rcoulomb;
213     nbp->rlistOuter_sq     = listParams->rlistOuter * listParams->rlistOuter;
214     nbp->rlistInner_sq     = listParams->rlistInner * listParams->rlistInner;
215     nbp->useDynamicPruning = listParams->useDynamicPruning;
216
217     nbp->sh_lj_ewald       = ic->sh_lj_ewald;
218     nbp->ewaldcoeff_lj     = ic->ewaldcoeff_lj;
219
220     nbp->rvdw_switch       = ic->rvdw_switch;
221     nbp->dispersion_shift  = ic->dispersion_shift;
222     nbp->repulsion_shift   = ic->repulsion_shift;
223     nbp->vdw_switch        = ic->vdw_switch;
224 }
225
226 /*! Initializes the nonbonded parameter data structure. */
227 static void init_nbparam(cu_nbparam_t              *nbp,
228                          const interaction_const_t *ic,
229                          const NbnxnListParameters *listParams,
230                          const nbnxn_atomdata_t    *nbat,
231                          const gmx_device_info_t   *dev_info)
232 {
233     int         ntypes;
234
235     ntypes  = nbat->ntype;
236
237     set_cutoff_parameters(nbp, ic, listParams);
238
239     /* The kernel code supports LJ combination rules (geometric and LB) for
240      * all kernel types, but we only generate useful combination rule kernels.
241      * We currently only use LJ combination rule (geometric and LB) kernels
242      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
243      * with PME and 20% with RF, the other kernels speed up about half as much.
244      * For LJ force-switch the geometric rule would give 7% speed-up, but this
245      * combination is rarely used. LJ force-switch with LB rule is more common,
246      * but gives only 1% speed-up.
247      */
248     if (ic->vdwtype == evdwCUT)
249     {
250         switch (ic->vdw_modifier)
251         {
252             case eintmodNONE:
253             case eintmodPOTSHIFT:
254                 switch (nbat->comb_rule)
255                 {
256                     case ljcrNONE:
257                         nbp->vdwtype = evdwCuCUT;
258                         break;
259                     case ljcrGEOM:
260                         nbp->vdwtype = evdwCuCUTCOMBGEOM;
261                         break;
262                     case ljcrLB:
263                         nbp->vdwtype = evdwCuCUTCOMBLB;
264                         break;
265                     default:
266                         gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
267                         break;
268                 }
269                 break;
270             case eintmodFORCESWITCH:
271                 nbp->vdwtype = evdwCuFSWITCH;
272                 break;
273             case eintmodPOTSWITCH:
274                 nbp->vdwtype = evdwCuPSWITCH;
275                 break;
276             default:
277                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
278                 break;
279         }
280     }
281     else if (ic->vdwtype == evdwPME)
282     {
283         if (ic->ljpme_comb_rule == ljcrGEOM)
284         {
285             assert(nbat->comb_rule == ljcrGEOM);
286             nbp->vdwtype = evdwCuEWALDGEOM;
287         }
288         else
289         {
290             assert(nbat->comb_rule == ljcrLB);
291             nbp->vdwtype = evdwCuEWALDLB;
292         }
293     }
294     else
295     {
296         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
297     }
298
299     if (ic->eeltype == eelCUT)
300     {
301         nbp->eeltype = eelCuCUT;
302     }
303     else if (EEL_RF(ic->eeltype))
304     {
305         nbp->eeltype = eelCuRF;
306     }
307     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
308     {
309         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
310         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
311     }
312     else
313     {
314         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
315         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
316     }
317
318     /* generate table for PME */
319     nbp->coulomb_tab = NULL;
320     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
321     {
322         init_ewald_coulomb_force_table(ic, nbp, dev_info);
323     }
324
325     /* set up LJ parameter lookup table */
326     if (!useLjCombRule(nbp))
327     {
328         initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
329                              nbat->nbfp, 2*ntypes*ntypes, dev_info);
330     }
331
332     /* set up LJ-PME parameter lookup table */
333     if (ic->vdwtype == evdwPME)
334     {
335         initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
336                              nbat->nbfp_comb, 2*ntypes, dev_info);
337     }
338 }
339
340 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
341  *  electrostatic type switching to twin cut-off (or back) if needed. */
342 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t    *nbv,
343                                         const interaction_const_t   *ic,
344                                         const NbnxnListParameters   *listParams)
345 {
346     if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
347     {
348         return;
349     }
350     gmx_nbnxn_cuda_t *nb    = nbv->gpu_nbv;
351     cu_nbparam_t     *nbp   = nb->nbparam;
352
353     set_cutoff_parameters(nbp, ic, listParams);
354
355     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
356                                                  nb->dev_info);
357
358     init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
359 }
360
361 /*! Initializes the pair list data structure. */
362 static void init_plist(cu_plist_t *pl)
363 {
364     /* initialize to NULL pointers to data that is not allocated here and will
365        need reallocation in nbnxn_gpu_init_pairlist */
366     pl->sci      = NULL;
367     pl->cj4      = NULL;
368     pl->imask    = NULL;
369     pl->excl     = NULL;
370
371     /* size -1 indicates that the respective array hasn't been initialized yet */
372     pl->na_c           = -1;
373     pl->nsci           = -1;
374     pl->sci_nalloc     = -1;
375     pl->ncj4           = -1;
376     pl->cj4_nalloc     = -1;
377     pl->nimask         = -1;
378     pl->imask_nalloc   = -1;
379     pl->nexcl          = -1;
380     pl->excl_nalloc    = -1;
381     pl->haveFreshList  = false;
382 }
383
384 /*! Initializes the timer data structure. */
385 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
386 {
387     /* The non-local counters/stream (second in the array) are needed only with DD. */
388     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
389     {
390         t->didPairlistH2D[i]  = false;
391         t->didPrune[i]        = false;
392         t->didRollingPrune[i] = false;
393     }
394 }
395
396 /*! Initializes the timings data structure. */
397 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
398 {
399     int i, j;
400
401     t->nb_h2d_t = 0.0;
402     t->nb_d2h_t = 0.0;
403     t->nb_c     = 0;
404     t->pl_h2d_t = 0.0;
405     t->pl_h2d_c = 0;
406     for (i = 0; i < 2; i++)
407     {
408         for (j = 0; j < 2; j++)
409         {
410             t->ktime[i][j].t = 0.0;
411             t->ktime[i][j].c = 0;
412         }
413     }
414     t->pruneTime.c        = 0;
415     t->pruneTime.t        = 0.0;
416     t->dynamicPruneTime.c = 0;
417     t->dynamicPruneTime.t = 0.0;
418 }
419
420 /*! Initializes simulation constant data. */
421 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t               *nb,
422                                   const interaction_const_t      *ic,
423                                   const NbnxnListParameters      *listParams,
424                                   const nbnxn_atomdata_t         *nbat)
425 {
426     init_atomdata_first(nb->atdat, nbat->ntype);
427     init_nbparam(nb->nbparam, ic, listParams, nbat, nb->dev_info);
428
429     /* clear energy and shift force outputs */
430     nbnxn_cuda_clear_e_fshift(nb);
431 }
432
433 void nbnxn_gpu_init(gmx_nbnxn_cuda_t         **p_nb,
434                     const gmx_device_info_t   *deviceInfo,
435                     const interaction_const_t *ic,
436                     const NbnxnListParameters *listParams,
437                     const nbnxn_atomdata_t    *nbat,
438                     int                        /*rank*/,
439                     gmx_bool                   bLocalAndNonlocal)
440 {
441     cudaError_t       stat;
442     gmx_nbnxn_cuda_t *nb;
443
444     if (p_nb == NULL)
445     {
446         return;
447     }
448
449     snew(nb, 1);
450     snew(nb->atdat, 1);
451     snew(nb->nbparam, 1);
452     snew(nb->plist[eintLocal], 1);
453     if (bLocalAndNonlocal)
454     {
455         snew(nb->plist[eintNonlocal], 1);
456     }
457
458     nb->bUseTwoStreams = bLocalAndNonlocal;
459
460     nb->timers = new cu_timers_t();
461     snew(nb->timings, 1);
462
463     /* init nbst */
464     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
465     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
466     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
467
468     init_plist(nb->plist[eintLocal]);
469
470     /* set device info, just point it to the right GPU among the detected ones */
471     nb->dev_info = deviceInfo;
472
473     /* local/non-local GPU streams */
474     stat = cudaStreamCreate(&nb->stream[eintLocal]);
475     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
476     if (nb->bUseTwoStreams)
477     {
478         init_plist(nb->plist[eintNonlocal]);
479
480         /* Note that the device we're running on does not have to support
481          * priorities, because we are querying the priority range which in this
482          * case will be a single value.
483          */
484         int highest_priority;
485         stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
486         CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
487
488         stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
489                                             cudaStreamDefault,
490                                             highest_priority);
491         CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
492     }
493
494     /* init events for sychronization (timing disabled for performance reasons!) */
495     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
496     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
497     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
498     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
499
500     /* WARNING: CUDA timings are incorrect with multiple streams.
501      *          This is the main reason why they are disabled by default.
502      */
503     // TODO: Consider turning on by default when we can detect nr of streams.
504     nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != NULL);
505
506     if (nb->bDoTime)
507     {
508         init_timers(nb->timers, nb->bUseTwoStreams);
509         init_timings(nb->timings);
510     }
511
512     /* set the kernel type for the current GPU */
513     /* pick L1 cache configuration */
514     nbnxn_cuda_set_cacheconfig(nb->dev_info);
515
516     nbnxn_cuda_init_const(nb, ic, listParams, nbat);
517
518     *p_nb = nb;
519
520     if (debug)
521     {
522         fprintf(debug, "Initialized CUDA data structures.\n");
523     }
524 }
525
526 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t       *nb,
527                              const nbnxn_pairlist_t *h_plist,
528                              int                     iloc)
529 {
530     if (canSkipWork(nb, iloc))
531     {
532         return;
533     }
534
535     char          sbuf[STRLEN];
536     bool          bDoTime    = nb->bDoTime;
537     cudaStream_t  stream     = nb->stream[iloc];
538     cu_plist_t   *d_plist    = nb->plist[iloc];
539
540     if (d_plist->na_c < 0)
541     {
542         d_plist->na_c = h_plist->na_ci;
543     }
544     else
545     {
546         if (d_plist->na_c != h_plist->na_ci)
547         {
548             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
549                     d_plist->na_c, h_plist->na_ci);
550             gmx_incons(sbuf);
551         }
552     }
553
554     if (bDoTime)
555     {
556         nb->timers->pl_h2d[iloc].openTimingRegion(stream);
557         nb->timers->didPairlistH2D[iloc] = true;
558     }
559
560     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
561                         &d_plist->nsci, &d_plist->sci_nalloc,
562                         h_plist->nsci,
563                         stream, true);
564
565     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
566                         &d_plist->ncj4, &d_plist->cj4_nalloc,
567                         h_plist->ncj4,
568                         stream, true);
569
570     /* this call only allocates space on the device (no data is transferred) */
571     cu_realloc_buffered((void **)&d_plist->imask, NULL, sizeof(*d_plist->imask),
572                         &d_plist->nimask, &d_plist->imask_nalloc,
573                         h_plist->ncj4*c_nbnxnGpuClusterpairSplit,
574                         stream, true);
575
576     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
577                         &d_plist->nexcl, &d_plist->excl_nalloc,
578                         h_plist->nexcl,
579                         stream, true);
580
581     if (bDoTime)
582     {
583         nb->timers->pl_h2d[iloc].closeTimingRegion(stream);
584     }
585
586     /* the next use of thist list we be the first one, so we need to prune */
587     d_plist->haveFreshList = true;
588 }
589
590 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t       *nb,
591                                const nbnxn_atomdata_t *nbatom)
592 {
593     cu_atomdata_t *adat  = nb->atdat;
594     cudaStream_t   ls    = nb->stream[eintLocal];
595
596     /* only if we have a dynamic box */
597     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
598     {
599         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
600                           SHIFTS * sizeof(*adat->shift_vec), ls);
601         adat->bShiftVecUploaded = true;
602     }
603 }
604
605 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
606 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
607 {
608     cudaError_t    stat;
609     cu_atomdata_t *adat  = nb->atdat;
610     cudaStream_t   ls    = nb->stream[eintLocal];
611
612     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
613     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
614 }
615
616 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
617 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
618 {
619     cudaError_t    stat;
620     cu_atomdata_t *adat  = nb->atdat;
621     cudaStream_t   ls    = nb->stream[eintLocal];
622
623     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
624     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
625     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
626     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
627     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
628     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
629 }
630
631 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
632 {
633     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
634     /* clear shift force array and energies if the outputs were
635        used in the current step */
636     if (flags & GMX_FORCE_VIRIAL)
637     {
638         nbnxn_cuda_clear_e_fshift(nb);
639     }
640 }
641
642 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t              *nb,
643                              const nbnxn_atomdata_t        *nbat)
644 {
645     cudaError_t    stat;
646     int            nalloc, natoms;
647     bool           realloced;
648     bool           bDoTime   = nb->bDoTime;
649     cu_timers_t   *timers    = nb->timers;
650     cu_atomdata_t *d_atdat   = nb->atdat;
651     cudaStream_t   ls        = nb->stream[eintLocal];
652
653     natoms    = nbat->natoms;
654     realloced = false;
655
656     if (bDoTime)
657     {
658         /* time async copy */
659         timers->atdat.openTimingRegion(ls);
660     }
661
662     /* need to reallocate if we have to copy more atoms than the amount of space
663        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
664     if (natoms > d_atdat->nalloc)
665     {
666         nalloc = over_alloc_small(natoms);
667
668         /* free up first if the arrays have already been initialized */
669         if (d_atdat->nalloc != -1)
670         {
671             freeDeviceBuffer(&d_atdat->f);
672             freeDeviceBuffer(&d_atdat->xq);
673             freeDeviceBuffer(&d_atdat->atom_types);
674             freeDeviceBuffer(&d_atdat->lj_comb);
675         }
676
677         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
678         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
679         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
680         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
681         if (useLjCombRule(nb->nbparam))
682         {
683             stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
684             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
685         }
686         else
687         {
688             stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
689             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
690         }
691
692         d_atdat->nalloc = nalloc;
693         realloced       = true;
694     }
695
696     d_atdat->natoms       = natoms;
697     d_atdat->natoms_local = nbat->natoms_local;
698
699     /* need to clear GPU f output if realloc happened */
700     if (realloced)
701     {
702         nbnxn_cuda_clear_f(nb, nalloc);
703     }
704
705     if (useLjCombRule(nb->nbparam))
706     {
707         cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
708                           natoms*sizeof(*d_atdat->lj_comb), ls);
709     }
710     else
711     {
712         cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
713                           natoms*sizeof(*d_atdat->atom_types), ls);
714     }
715
716     if (bDoTime)
717     {
718         timers->atdat.closeTimingRegion(ls);
719     }
720 }
721
722 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
723                                           const gmx_device_info_t *dev_info)
724 {
725     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
726     {
727         destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj,
728                                 dev_info);
729     }
730 }
731
732 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
733 {
734     cudaError_t      stat;
735     cu_atomdata_t   *atdat;
736     cu_nbparam_t    *nbparam;
737
738     if (nb == NULL)
739     {
740         return;
741     }
742
743     atdat       = nb->atdat;
744     nbparam     = nb->nbparam;
745
746     nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
747
748     stat = cudaEventDestroy(nb->nonlocal_done);
749     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
750     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
751     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
752
753     delete nb->timers;
754     if (nb->bDoTime)
755     {
756         /* The non-local counters/stream (second in the array) are needed only with DD. */
757         for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
758         {
759             stat = cudaStreamDestroy(nb->stream[i]);
760             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
761         }
762     }
763
764     if (!useLjCombRule(nb->nbparam))
765     {
766         destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj,
767                                 nb->dev_info);
768
769     }
770
771     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
772     {
773         destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj,
774                                 nb->dev_info);
775     }
776
777     stat = cudaFree(atdat->shift_vec);
778     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
779     stat = cudaFree(atdat->fshift);
780     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
781
782     stat = cudaFree(atdat->e_lj);
783     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
784     stat = cudaFree(atdat->e_el);
785     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
786
787     freeDeviceBuffer(&atdat->f);
788     freeDeviceBuffer(&atdat->xq);
789     freeDeviceBuffer(&atdat->atom_types);
790     freeDeviceBuffer(&atdat->lj_comb);
791
792     /* Free plist */
793     auto *plist = nb->plist[eintLocal];
794     freeDeviceBuffer(&plist->sci);
795     freeDeviceBuffer(&plist->cj4);
796     freeDeviceBuffer(&plist->imask);
797     freeDeviceBuffer(&plist->excl);
798     sfree(plist);
799     if (nb->bUseTwoStreams)
800     {
801         auto *plist_nl = nb->plist[eintNonlocal];
802         freeDeviceBuffer(&plist_nl->sci);
803         freeDeviceBuffer(&plist_nl->cj4);
804         freeDeviceBuffer(&plist_nl->imask);
805         freeDeviceBuffer(&plist_nl->excl);
806         sfree(plist_nl);
807     }
808
809     sfree(atdat);
810     sfree(nbparam);
811     sfree(nb->timings);
812     sfree(nb);
813
814     if (debug)
815     {
816         fprintf(debug, "Cleaned up CUDA data structures.\n");
817     }
818 }
819
820 //! This function is documented in the header file
821 gmx_wallclock_gpu_nbnxn_t *nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
822 {
823     return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
824 }
825
826 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
827 {
828     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
829     {
830         init_timings(nbv->gpu_nbv->timings);
831     }
832 }
833
834 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
835 {
836     return nb != NULL ?
837            gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
838
839 }
840
841 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
842 {
843     return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
844             (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
845 }