StatePropagatorDataGpu object to manage GPU forces, positions and velocities buffers
[alexxy/gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_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,2019, 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 // TODO We would like to move this down, but the way gmx_nbnxn_gpu_t
48 //      is currently declared means this has to be before gpu_types.h
49 #include "nbnxm_cuda_types.h"
50
51 // TODO Remove this comment when the above order issue is resolved
52 #include "gromacs/gpu_utils/cudautils.cuh"
53 #include "gromacs/gpu_utils/gpu_utils.h"
54 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
55 #include "gromacs/gpu_utils/pmalloc_cuda.h"
56 #include "gromacs/hardware/gpu_hw_info.h"
57 #include "gromacs/math/vectypes.h"
58 #include "gromacs/mdlib/force_flags.h"
59 #include "gromacs/mdtypes/interaction_const.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/nbnxm/atomdata.h"
62 #include "gromacs/nbnxm/gpu_data_mgmt.h"
63 #include "gromacs/nbnxm/gridset.h"
64 #include "gromacs/nbnxm/nbnxm.h"
65 #include "gromacs/nbnxm/nbnxm_gpu.h"
66 #include "gromacs/nbnxm/pairlistsets.h"
67 #include "gromacs/pbcutil/ishift.h"
68 #include "gromacs/timing/gpu_timing.h"
69 #include "gromacs/utility/basedefinitions.h"
70 #include "gromacs/utility/cstringutil.h"
71 #include "gromacs/utility/fatalerror.h"
72 #include "gromacs/utility/real.h"
73 #include "gromacs/utility/smalloc.h"
74
75 #include "nbnxm_cuda.h"
76
77 namespace Nbnxm
78 {
79
80 /* This is a heuristically determined parameter for the Kepler
81  * and Maxwell architectures for the minimum size of ci lists by multiplying
82  * this constant with the # of multiprocessors on the current device.
83  * Since the maximum number of blocks per multiprocessor is 16, the ideal
84  * count for small systems is 32 or 48 blocks per multiprocessor. Because
85  * there is a bit of fluctuations in the generated block counts, we use
86  * a target of 44 instead of the ideal value of 48.
87  */
88 static unsigned int gpu_min_ci_balanced_factor = 44;
89
90 /* Fw. decl. */
91 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
92
93 /* Fw. decl, */
94 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam);
95
96 /*! \brief Return whether combination rules are used.
97  *
98  * \param[in]   pointer to nonbonded paramter struct
99  * \return      true if combination rules are used in this run, false otherwise
100  */
101 static inline bool useLjCombRule(const cu_nbparam_t  *nbparam)
102 {
103     return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
104             nbparam->vdwtype == evdwCuCUTCOMBLB);
105 }
106
107 /*! \brief Initialized the Ewald Coulomb correction GPU table.
108
109     Tabulates the Ewald Coulomb force and initializes the size/scale
110     and the table GPU array. If called with an already allocated table,
111     it just re-uploads the table.
112  */
113 static void init_ewald_coulomb_force_table(const EwaldCorrectionTables &tables,
114                                            cu_nbparam_t                *nbp)
115 {
116     if (nbp->coulomb_tab != nullptr)
117     {
118         nbnxn_cuda_free_nbparam_table(nbp);
119     }
120
121     nbp->coulomb_tab_scale = tables.scale;
122     initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
123                          tables.tableF.data(), tables.tableF.size());
124 }
125
126
127 /*! Initializes the atomdata structure first time, it only gets filled at
128     pair-search. */
129 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
130 {
131     cudaError_t stat;
132
133     ad->ntypes  = ntypes;
134     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
135     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
136     ad->bShiftVecUploaded = false;
137
138     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
139     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
140
141     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
142     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
143     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
144     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
145
146     /* initialize to nullptr poiters to data that is not allocated here and will
147        need reallocation in nbnxn_cuda_init_atomdata */
148     ad->xq = nullptr;
149     ad->f  = nullptr;
150
151     /* size -1 indicates that the respective array hasn't been initialized yet */
152     ad->natoms = -1;
153     ad->nalloc = -1;
154 }
155
156 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
157     earlier GPUs, single or twin cut-off. */
158 static int pick_ewald_kernel_type(const interaction_const_t &ic)
159 {
160     bool bTwinCut = (ic.rcoulomb != ic.rvdw);
161     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
162     int  kernel_type;
163
164     /* Benchmarking/development environment variables to force the use of
165        analytical or tabulated Ewald kernel. */
166     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
167     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
168
169     if (bForceAnalyticalEwald && bForceTabulatedEwald)
170     {
171         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
172                    "requested through environment variables.");
173     }
174
175     /* By default use analytical Ewald. */
176     bUseAnalyticalEwald = true;
177     if (bForceAnalyticalEwald)
178     {
179         if (debug)
180         {
181             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
182         }
183     }
184     else if (bForceTabulatedEwald)
185     {
186         bUseAnalyticalEwald = false;
187
188         if (debug)
189         {
190             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
191         }
192     }
193
194     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
195        forces it (use it for debugging/benchmarking only). */
196     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
197     {
198         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
199     }
200     else
201     {
202         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
203     }
204
205     return kernel_type;
206 }
207
208 /*! Copies all parameters related to the cut-off from ic to nbp */
209 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
210                                   const interaction_const_t *ic,
211                                   const PairlistParams      &listParams)
212 {
213     nbp->ewald_beta        = ic->ewaldcoeff_q;
214     nbp->sh_ewald          = ic->sh_ewald;
215     nbp->epsfac            = ic->epsfac;
216     nbp->two_k_rf          = 2.0 * ic->k_rf;
217     nbp->c_rf              = ic->c_rf;
218     nbp->rvdw_sq           = ic->rvdw * ic->rvdw;
219     nbp->rcoulomb_sq       = ic->rcoulomb * ic->rcoulomb;
220     nbp->rlistOuter_sq     = listParams.rlistOuter * listParams.rlistOuter;
221     nbp->rlistInner_sq     = listParams.rlistInner * listParams.rlistInner;
222     nbp->useDynamicPruning = listParams.useDynamicPruning;
223
224     nbp->sh_lj_ewald       = ic->sh_lj_ewald;
225     nbp->ewaldcoeff_lj     = ic->ewaldcoeff_lj;
226
227     nbp->rvdw_switch       = ic->rvdw_switch;
228     nbp->dispersion_shift  = ic->dispersion_shift;
229     nbp->repulsion_shift   = ic->repulsion_shift;
230     nbp->vdw_switch        = ic->vdw_switch;
231 }
232
233 /*! Initializes the nonbonded parameter data structure. */
234 static void init_nbparam(cu_nbparam_t                   *nbp,
235                          const interaction_const_t      *ic,
236                          const PairlistParams           &listParams,
237                          const nbnxn_atomdata_t::Params &nbatParams)
238 {
239     int         ntypes;
240
241     ntypes  = nbatParams.numTypes;
242
243     set_cutoff_parameters(nbp, ic, listParams);
244
245     /* The kernel code supports LJ combination rules (geometric and LB) for
246      * all kernel types, but we only generate useful combination rule kernels.
247      * We currently only use LJ combination rule (geometric and LB) kernels
248      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
249      * with PME and 20% with RF, the other kernels speed up about half as much.
250      * For LJ force-switch the geometric rule would give 7% speed-up, but this
251      * combination is rarely used. LJ force-switch with LB rule is more common,
252      * but gives only 1% speed-up.
253      */
254     if (ic->vdwtype == evdwCUT)
255     {
256         switch (ic->vdw_modifier)
257         {
258             case eintmodNONE:
259             case eintmodPOTSHIFT:
260                 switch (nbatParams.comb_rule)
261                 {
262                     case ljcrNONE:
263                         nbp->vdwtype = evdwCuCUT;
264                         break;
265                     case ljcrGEOM:
266                         nbp->vdwtype = evdwCuCUTCOMBGEOM;
267                         break;
268                     case ljcrLB:
269                         nbp->vdwtype = evdwCuCUTCOMBLB;
270                         break;
271                     default:
272                         gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
273                 }
274                 break;
275             case eintmodFORCESWITCH:
276                 nbp->vdwtype = evdwCuFSWITCH;
277                 break;
278             case eintmodPOTSWITCH:
279                 nbp->vdwtype = evdwCuPSWITCH;
280                 break;
281             default:
282                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
283         }
284     }
285     else if (ic->vdwtype == evdwPME)
286     {
287         if (ic->ljpme_comb_rule == ljcrGEOM)
288         {
289             assert(nbatParams.comb_rule == ljcrGEOM);
290             nbp->vdwtype = evdwCuEWALDGEOM;
291         }
292         else
293         {
294             assert(nbatParams.comb_rule == ljcrLB);
295             nbp->vdwtype = evdwCuEWALDLB;
296         }
297     }
298     else
299     {
300         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
301     }
302
303     if (ic->eeltype == eelCUT)
304     {
305         nbp->eeltype = eelCuCUT;
306     }
307     else if (EEL_RF(ic->eeltype))
308     {
309         nbp->eeltype = eelCuRF;
310     }
311     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
312     {
313         nbp->eeltype = pick_ewald_kernel_type(*ic);
314     }
315     else
316     {
317         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
318         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
319     }
320
321     /* generate table for PME */
322     nbp->coulomb_tab = nullptr;
323     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
324     {
325         GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
326         init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
327     }
328
329     /* set up LJ parameter lookup table */
330     if (!useLjCombRule(nbp))
331     {
332         initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
333                              nbatParams.nbfp.data(), 2*ntypes*ntypes);
334     }
335
336     /* set up LJ-PME parameter lookup table */
337     if (ic->vdwtype == evdwPME)
338     {
339         initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
340                              nbatParams.nbfp_comb.data(), 2*ntypes);
341     }
342 }
343
344 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
345  *  electrostatic type switching to twin cut-off (or back) if needed. */
346 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t    *nbv,
347                                   const interaction_const_t   *ic)
348 {
349     if (!nbv || !nbv->useGpu())
350     {
351         return;
352     }
353     cu_nbparam_t *nbp   = nbv->gpu_nbv->nbparam;
354
355     set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
356
357     nbp->eeltype        = pick_ewald_kernel_type(*ic);
358
359     GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
360     init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
361 }
362
363 /*! Initializes the pair list data structure. */
364 static void init_plist(cu_plist_t *pl)
365 {
366     /* initialize to nullptr pointers to data that is not allocated here and will
367        need reallocation in nbnxn_gpu_init_pairlist */
368     pl->sci      = nullptr;
369     pl->cj4      = nullptr;
370     pl->imask    = nullptr;
371     pl->excl     = nullptr;
372
373     /* size -1 indicates that the respective array hasn't been initialized yet */
374     pl->na_c           = -1;
375     pl->nsci           = -1;
376     pl->sci_nalloc     = -1;
377     pl->ncj4           = -1;
378     pl->cj4_nalloc     = -1;
379     pl->nimask         = -1;
380     pl->imask_nalloc   = -1;
381     pl->nexcl          = -1;
382     pl->excl_nalloc    = -1;
383     pl->haveFreshList  = false;
384 }
385
386 /*! Initializes the timings data structure. */
387 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
388 {
389     int i, j;
390
391     t->nb_h2d_t = 0.0;
392     t->nb_d2h_t = 0.0;
393     t->nb_c     = 0;
394     t->pl_h2d_t = 0.0;
395     t->pl_h2d_c = 0;
396     for (i = 0; i < 2; i++)
397     {
398         for (j = 0; j < 2; j++)
399         {
400             t->ktime[i][j].t = 0.0;
401             t->ktime[i][j].c = 0;
402         }
403     }
404     t->pruneTime.c        = 0;
405     t->pruneTime.t        = 0.0;
406     t->dynamicPruneTime.c = 0;
407     t->dynamicPruneTime.t = 0.0;
408 }
409
410 /*! Initializes simulation constant data. */
411 static void cuda_init_const(gmx_nbnxn_cuda_t               *nb,
412                             const interaction_const_t      *ic,
413                             const PairlistParams           &listParams,
414                             const nbnxn_atomdata_t::Params &nbatParams)
415 {
416     init_atomdata_first(nb->atdat, nbatParams.numTypes);
417     init_nbparam(nb->nbparam, ic, listParams, nbatParams);
418
419     /* clear energy and shift force outputs */
420     nbnxn_cuda_clear_e_fshift(nb);
421 }
422
423 gmx_nbnxn_cuda_t *
424 gpu_init(const gmx_device_info_t   *deviceInfo,
425          const interaction_const_t *ic,
426          const PairlistParams      &listParams,
427          const nbnxn_atomdata_t    *nbat,
428          int                        /*rank*/,
429          gmx_bool                   bLocalAndNonlocal)
430 {
431     cudaError_t       stat;
432
433     gmx_nbnxn_cuda_t *nb;
434     snew(nb, 1);
435     snew(nb->atdat, 1);
436     snew(nb->nbparam, 1);
437     snew(nb->plist[InteractionLocality::Local], 1);
438     if (bLocalAndNonlocal)
439     {
440         snew(nb->plist[InteractionLocality::NonLocal], 1);
441     }
442
443     nb->bUseTwoStreams = bLocalAndNonlocal;
444
445     nb->timers = new cu_timers_t();
446     snew(nb->timings, 1);
447
448     /* init nbst */
449     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
450     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
451     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
452
453     init_plist(nb->plist[InteractionLocality::Local]);
454
455     /* set device info, just point it to the right GPU among the detected ones */
456     nb->dev_info = deviceInfo;
457
458     /* local/non-local GPU streams */
459     stat = cudaStreamCreate(&nb->stream[InteractionLocality::Local]);
460     CU_RET_ERR(stat, "cudaStreamCreate on stream[InterationLocality::Local] failed");
461     if (nb->bUseTwoStreams)
462     {
463         init_plist(nb->plist[InteractionLocality::NonLocal]);
464
465         /* Note that the device we're running on does not have to support
466          * priorities, because we are querying the priority range which in this
467          * case will be a single value.
468          */
469         int highest_priority;
470         stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
471         CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
472
473         stat = cudaStreamCreateWithPriority(&nb->stream[InteractionLocality::NonLocal],
474                                             cudaStreamDefault,
475                                             highest_priority);
476         CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[InteractionLocality::NonLocal] failed");
477     }
478
479     /* init events for sychronization (timing disabled for performance reasons!) */
480     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
481     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
482     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
483     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
484
485     nb->xAvailableOnDevice   = new GpuEventSynchronizer();
486     nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
487
488     /* WARNING: CUDA timings are incorrect with multiple streams.
489      *          This is the main reason why they are disabled by default.
490      */
491     // TODO: Consider turning on by default when we can detect nr of streams.
492     nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
493
494     if (nb->bDoTime)
495     {
496         init_timings(nb->timings);
497     }
498
499     /* set the kernel type for the current GPU */
500     /* pick L1 cache configuration */
501     cuda_set_cacheconfig();
502
503     cuda_init_const(nb, ic, listParams, nbat->params());
504
505     nb->atomIndicesSize          = 0;
506     nb->atomIndicesSize_alloc    = 0;
507     nb->ncxy_na                  = 0;
508     nb->ncxy_na_alloc            = 0;
509     nb->ncxy_ind                 = 0;
510     nb->ncxy_ind_alloc           = 0;
511     nb->ncell                    = 0;
512     nb->ncell_alloc              = 0;
513
514     if (debug)
515     {
516         fprintf(debug, "Initialized CUDA data structures.\n");
517     }
518
519     return nb;
520 }
521
522 void gpu_init_pairlist(gmx_nbnxn_cuda_t          *nb,
523                        const NbnxnPairlistGpu    *h_plist,
524                        const InteractionLocality  iloc)
525 {
526     char          sbuf[STRLEN];
527     bool          bDoTime    =  (nb->bDoTime && !h_plist->sci.empty());
528     cudaStream_t  stream     = nb->stream[iloc];
529     cu_plist_t   *d_plist    = nb->plist[iloc];
530
531     if (d_plist->na_c < 0)
532     {
533         d_plist->na_c = h_plist->na_ci;
534     }
535     else
536     {
537         if (d_plist->na_c != h_plist->na_ci)
538         {
539             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
540                     d_plist->na_c, h_plist->na_ci);
541             gmx_incons(sbuf);
542         }
543     }
544
545     gpu_timers_t::Interaction &iTimers = nb->timers->interaction[iloc];
546
547     if (bDoTime)
548     {
549         iTimers.pl_h2d.openTimingRegion(stream);
550         iTimers.didPairlistH2D = true;
551     }
552
553     Context context = nullptr;
554
555     reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(),
556                            &d_plist->nsci, &d_plist->sci_nalloc, context);
557     copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(),
558                        stream, GpuApiCallBehavior::Async,
559                        bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
560
561     reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(),
562                            &d_plist->ncj4, &d_plist->cj4_nalloc, context);
563     copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(),
564                        stream, GpuApiCallBehavior::Async,
565                        bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
566
567     reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size()*c_nbnxnGpuClusterpairSplit,
568                            &d_plist->nimask, &d_plist->imask_nalloc, context);
569
570     reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(),
571                            &d_plist->nexcl, &d_plist->excl_nalloc, context);
572     copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(),
573                        stream, GpuApiCallBehavior::Async,
574                        bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
575
576     if (bDoTime)
577     {
578         iTimers.pl_h2d.closeTimingRegion(stream);
579     }
580
581     /* the next use of thist list we be the first one, so we need to prune */
582     d_plist->haveFreshList = true;
583 }
584
585 void gpu_upload_shiftvec(gmx_nbnxn_cuda_t       *nb,
586                          const nbnxn_atomdata_t *nbatom)
587 {
588     cu_atomdata_t *adat  = nb->atdat;
589     cudaStream_t   ls    = nb->stream[InteractionLocality::Local];
590
591     /* only if we have a dynamic box */
592     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
593     {
594         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(),
595                           SHIFTS * sizeof(*adat->shift_vec), ls);
596         adat->bShiftVecUploaded = true;
597     }
598 }
599
600 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
601 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
602 {
603     cudaError_t    stat;
604     cu_atomdata_t *adat  = nb->atdat;
605     cudaStream_t   ls    = nb->stream[InteractionLocality::Local];
606
607     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
608     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
609 }
610
611 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
612 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
613 {
614     cudaError_t    stat;
615     cu_atomdata_t *adat  = nb->atdat;
616     cudaStream_t   ls    = nb->stream[InteractionLocality::Local];
617
618     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
619     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
620     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
621     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
622     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
623     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
624 }
625
626 void gpu_clear_outputs(gmx_nbnxn_cuda_t *nb,
627                        bool              computeVirial)
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 (computeVirial)
633     {
634         nbnxn_cuda_clear_e_fshift(nb);
635     }
636 }
637
638 void 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[InteractionLocality::Local];
648
649     natoms    = nbat->numAtoms();
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->params().lj_comb.data(),
704                           natoms*sizeof(*d_atdat->lj_comb), ls);
705     }
706     else
707     {
708         cu_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(),
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 {
720     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
721     {
722         destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
723     }
724 }
725
726 void gpu_free(gmx_nbnxn_cuda_t *nb)
727 {
728     cudaError_t      stat;
729     cu_atomdata_t   *atdat;
730     cu_nbparam_t    *nbparam;
731
732     if (nb == nullptr)
733     {
734         return;
735     }
736
737     atdat       = nb->atdat;
738     nbparam     = nb->nbparam;
739
740     nbnxn_cuda_free_nbparam_table(nbparam);
741
742     stat = cudaEventDestroy(nb->nonlocal_done);
743     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
744     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
745     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
746
747     delete nb->timers;
748     if (nb->bDoTime)
749     {
750         /* The non-local counters/stream (second in the array) are needed only with DD. */
751         for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
752         {
753             stat = cudaStreamDestroy(nb->stream[i]);
754             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
755         }
756     }
757
758     if (!useLjCombRule(nb->nbparam))
759     {
760         destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj);
761
762     }
763
764     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
765     {
766         destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
767     }
768
769     stat = cudaFree(atdat->shift_vec);
770     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
771     stat = cudaFree(atdat->fshift);
772     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
773
774     stat = cudaFree(atdat->e_lj);
775     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
776     stat = cudaFree(atdat->e_el);
777     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
778
779     freeDeviceBuffer(&atdat->f);
780     freeDeviceBuffer(&atdat->xq);
781     freeDeviceBuffer(&atdat->atom_types);
782     freeDeviceBuffer(&atdat->lj_comb);
783
784     /* Free plist */
785     auto *plist = nb->plist[InteractionLocality::Local];
786     freeDeviceBuffer(&plist->sci);
787     freeDeviceBuffer(&plist->cj4);
788     freeDeviceBuffer(&plist->imask);
789     freeDeviceBuffer(&plist->excl);
790     sfree(plist);
791     if (nb->bUseTwoStreams)
792     {
793         auto *plist_nl = nb->plist[InteractionLocality::NonLocal];
794         freeDeviceBuffer(&plist_nl->sci);
795         freeDeviceBuffer(&plist_nl->cj4);
796         freeDeviceBuffer(&plist_nl->imask);
797         freeDeviceBuffer(&plist_nl->excl);
798         sfree(plist_nl);
799     }
800
801     /* Free nbst */
802     pfree(nb->nbst.e_lj);
803     nb->nbst.e_lj = nullptr;
804
805     pfree(nb->nbst.e_el);
806     nb->nbst.e_el = nullptr;
807
808     pfree(nb->nbst.fshift);
809     nb->nbst.fshift = nullptr;
810
811     sfree(atdat);
812     sfree(nbparam);
813     sfree(nb->timings);
814     sfree(nb);
815
816     if (debug)
817     {
818         fprintf(debug, "Cleaned up CUDA data structures.\n");
819     }
820 }
821
822 //! This function is documented in the header file
823 gmx_wallclock_gpu_nbnxn_t *gpu_get_timings(gmx_nbnxn_cuda_t *nb)
824 {
825     return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
826 }
827
828 void gpu_reset_timings(nonbonded_verlet_t* nbv)
829 {
830     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
831     {
832         init_timings(nbv->gpu_nbv->timings);
833     }
834 }
835
836 int gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
837 {
838     return nb != nullptr ?
839            gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
840
841 }
842
843 gmx_bool gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
844 {
845     return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
846             (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
847 }
848
849 void *gpu_get_command_stream(gmx_nbnxn_gpu_t           *nb,
850                              const InteractionLocality  iloc)
851 {
852     assert(nb);
853
854     return static_cast<void *>(&nb->stream[iloc]);
855 }
856
857 void *gpu_get_xq(gmx_nbnxn_gpu_t *nb)
858 {
859     assert(nb);
860
861     return static_cast<void *>(nb->atdat->xq);
862 }
863
864 void *gpu_get_f(gmx_nbnxn_gpu_t *nb)
865 {
866     assert(nb);
867
868     return static_cast<void *>(nb->atdat->f);
869 }
870
871 rvec *gpu_get_fshift(gmx_nbnxn_gpu_t *nb)
872 {
873     assert(nb);
874
875     return reinterpret_cast<rvec *>(nb->atdat->fshift);
876 }
877
878 /* Initialization for X buffer operations on GPU. */
879 /* TODO  Remove explicit pinning from host arrays from here and manage in a more natural way*/
880 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet            &gridSet,
881                                 gmx_nbnxn_gpu_t                 *gpu_nbv)
882 {
883     cudaStream_t                     stream    = gpu_nbv->stream[InteractionLocality::Local];
884     bool                             bDoTime   = gpu_nbv->bDoTime;
885     const int maxNumColumns                    = gridSet.numColumnsMax();
886
887     reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns*gridSet.grids().size(),
888                            &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, nullptr);
889     reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns*gridSet.grids().size(),
890                            &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, nullptr);
891
892     for (unsigned int g = 0; g < gridSet.grids().size(); g++)
893     {
894
895         const Nbnxm::Grid  &grid       = gridSet.grids()[g];
896
897         const int           numColumns        = grid.numColumns();
898         const int          *atomIndices       = gridSet.atomIndices().data();
899         const int           atomIndicesSize   = gridSet.atomIndices().size();
900         const int          *cxy_na            = grid.cxy_na().data();
901         const int          *cxy_ind           = grid.cxy_ind().data();
902
903         reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize, &gpu_nbv->atomIndicesSize_alloc, nullptr);
904
905         if (atomIndicesSize > 0)
906         {
907
908             if (bDoTime)
909             {
910                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
911             }
912
913             copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, stream, GpuApiCallBehavior::Async, nullptr);
914
915             if (bDoTime)
916             {
917                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
918             }
919
920         }
921
922         if (numColumns > 0)
923         {
924             if (bDoTime)
925             {
926                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
927             }
928
929             int* destPtr = &gpu_nbv->cxy_na[maxNumColumns*g];
930             copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
931
932             if (bDoTime)
933             {
934                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
935             }
936
937             if (bDoTime)
938             {
939                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
940             }
941
942             destPtr = &gpu_nbv->cxy_ind[maxNumColumns*g];
943             copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
944
945             if (bDoTime)
946             {
947                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
948             }
949
950         }
951     }
952
953     // The above data is transferred on the local stream but is a
954     // dependency of the nonlocal stream (specifically the nonlocal X
955     // buf ops kernel).  We therefore set a dependency to ensure
956     // that the nonlocal stream waits on the local stream here.
957     // This call records an event in the local stream:
958     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
959     // ...and this call instructs the nonlocal stream to wait on that event:
960     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
961
962     return;
963 }
964
965 /* Initialization for F buffer operations on GPU. */
966 void nbnxn_gpu_init_add_nbat_f_to_f(const int                *cell,
967                                     gmx_nbnxn_gpu_t          *gpu_nbv,
968                                     int                       natoms_total)
969 {
970
971     cudaStream_t         stream  = gpu_nbv->stream[InteractionLocality::Local];
972
973     if (natoms_total > 0)
974     {
975         reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc, nullptr);
976         copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, stream, GpuApiCallBehavior::Async, nullptr);
977     }
978
979     return;
980 }
981
982 } // namespace Nbnxm