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