b2fc758417a814784fcc3ee1b9abf8c7881aa5de
[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 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \file
37  *  \brief Define CUDA implementation of nbnxn_gpu_data_mgmt.h
38  *
39  *  \author Szilard Pall <pall.szilard@gmail.com>
40  */
41 #include "gmxpre.h"
42
43 #include <assert.h>
44 #include <stdarg.h>
45 #include <stdio.h>
46 #include <stdlib.h>
47
48 // TODO We would like to move this down, but the way NbnxmGpu
49 //      is currently declared means this has to be before gpu_types.h
50 #include "nbnxm_cuda_types.h"
51
52 // TODO Remove this comment when the above order issue is resolved
53 #include "gromacs/gpu_utils/cudautils.cuh"
54 #include "gromacs/gpu_utils/gpu_utils.h"
55 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
56 #include "gromacs/gpu_utils/pmalloc_cuda.h"
57 #include "gromacs/hardware/gpu_hw_info.h"
58 #include "gromacs/math/vectypes.h"
59 #include "gromacs/mdlib/force_flags.h"
60 #include "gromacs/mdtypes/interaction_const.h"
61 #include "gromacs/mdtypes/md_enums.h"
62 #include "gromacs/nbnxm/atomdata.h"
63 #include "gromacs/nbnxm/gpu_data_mgmt.h"
64 #include "gromacs/nbnxm/gridset.h"
65 #include "gromacs/nbnxm/nbnxm.h"
66 #include "gromacs/nbnxm/nbnxm_gpu.h"
67 #include "gromacs/nbnxm/pairlistsets.h"
68 #include "gromacs/pbcutil/ishift.h"
69 #include "gromacs/timing/gpu_timing.h"
70 #include "gromacs/utility/basedefinitions.h"
71 #include "gromacs/utility/cstringutil.h"
72 #include "gromacs/utility/fatalerror.h"
73 #include "gromacs/utility/real.h"
74 #include "gromacs/utility/smalloc.h"
75
76 #include "nbnxm_cuda.h"
77
78 namespace Nbnxm
79 {
80
81 /* This is a heuristically determined parameter for the Kepler
82  * and Maxwell architectures for the minimum size of ci lists by multiplying
83  * this constant with the # of multiprocessors on the current device.
84  * Since the maximum number of blocks per multiprocessor is 16, the ideal
85  * count for small systems is 32 or 48 blocks per multiprocessor. Because
86  * there is a bit of fluctuations in the generated block counts, we use
87  * a target of 44 instead of the ideal value of 48.
88  */
89 static unsigned int gpu_min_ci_balanced_factor = 44;
90
91 /* Fw. decl. */
92 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
93
94 /* Fw. decl, */
95 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t* nbparam);
96
97 /*! \brief Return whether combination rules are used.
98  *
99  * \param[in]   pointer to nonbonded paramter struct
100  * \return      true if combination rules are used in this run, false otherwise
101  */
102 static inline bool useLjCombRule(const cu_nbparam_t* nbparam)
103 {
104     return (nbparam->vdwtype == evdwCuCUTCOMBGEOM || 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, cu_nbparam_t* nbp)
114 {
115     if (nbp->coulomb_tab != nullptr)
116     {
117         nbnxn_cuda_free_nbparam_table(nbp);
118     }
119
120     nbp->coulomb_tab_scale = tables.scale;
121     initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj, tables.tableF.data(),
122                          tables.tableF.size());
123 }
124
125
126 /*! Initializes the atomdata structure first time, it only gets filled at
127     pair-search. */
128 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes)
129 {
130     cudaError_t stat;
131
132     ad->ntypes = ntypes;
133     stat       = cudaMalloc((void**)&ad->shift_vec, SHIFTS * sizeof(*ad->shift_vec));
134     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
135     ad->bShiftVecUploaded = false;
136
137     stat = cudaMalloc((void**)&ad->fshift, SHIFTS * sizeof(*ad->fshift));
138     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
139
140     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
141     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
142     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
143     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
144
145     /* initialize to nullptr poiters to data that is not allocated here and will
146        need reallocation in nbnxn_cuda_init_atomdata */
147     ad->xq = nullptr;
148     ad->f  = nullptr;
149
150     /* size -1 indicates that the respective array hasn't been initialized yet */
151     ad->natoms = -1;
152     ad->nalloc = -1;
153 }
154
155 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
156     earlier GPUs, single or twin cut-off. */
157 static int pick_ewald_kernel_type(const interaction_const_t& ic)
158 {
159     bool bTwinCut = (ic.rcoulomb != ic.rvdw);
160     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
161     int  kernel_type;
162
163     /* Benchmarking/development environment variables to force the use of
164        analytical or tabulated Ewald kernel. */
165     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
166     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
167
168     if (bForceAnalyticalEwald && bForceTabulatedEwald)
169     {
170         gmx_incons(
171                 "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, const interaction_const_t* ic, const PairlistParams& listParams)
210 {
211     nbp->ewald_beta        = ic->ewaldcoeff_q;
212     nbp->sh_ewald          = ic->sh_ewald;
213     nbp->epsfac            = ic->epsfac;
214     nbp->two_k_rf          = 2.0 * ic->k_rf;
215     nbp->c_rf              = ic->c_rf;
216     nbp->rvdw_sq           = ic->rvdw * ic->rvdw;
217     nbp->rcoulomb_sq       = ic->rcoulomb * ic->rcoulomb;
218     nbp->rlistOuter_sq     = listParams.rlistOuter * listParams.rlistOuter;
219     nbp->rlistInner_sq     = listParams.rlistInner * listParams.rlistInner;
220     nbp->useDynamicPruning = listParams.useDynamicPruning;
221
222     nbp->sh_lj_ewald   = ic->sh_lj_ewald;
223     nbp->ewaldcoeff_lj = ic->ewaldcoeff_lj;
224
225     nbp->rvdw_switch      = ic->rvdw_switch;
226     nbp->dispersion_shift = ic->dispersion_shift;
227     nbp->repulsion_shift  = ic->repulsion_shift;
228     nbp->vdw_switch       = ic->vdw_switch;
229 }
230
231 /*! Initializes the nonbonded parameter data structure. */
232 static void init_nbparam(cu_nbparam_t*                   nbp,
233                          const interaction_const_t*      ic,
234                          const PairlistParams&           listParams,
235                          const nbnxn_atomdata_t::Params& nbatParams)
236 {
237     int ntypes;
238
239     ntypes = nbatParams.numTypes;
240
241     set_cutoff_parameters(nbp, ic, listParams);
242
243     /* The kernel code supports LJ combination rules (geometric and LB) for
244      * all kernel types, but we only generate useful combination rule kernels.
245      * We currently only use LJ combination rule (geometric and LB) kernels
246      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
247      * with PME and 20% with RF, the other kernels speed up about half as much.
248      * For LJ force-switch the geometric rule would give 7% speed-up, but this
249      * combination is rarely used. LJ force-switch with LB rule is more common,
250      * but gives only 1% speed-up.
251      */
252     if (ic->vdwtype == evdwCUT)
253     {
254         switch (ic->vdw_modifier)
255         {
256             case eintmodNONE:
257             case eintmodPOTSHIFT:
258                 switch (nbatParams.comb_rule)
259                 {
260                     case ljcrNONE: nbp->vdwtype = evdwCuCUT; break;
261                     case ljcrGEOM: nbp->vdwtype = evdwCuCUTCOMBGEOM; break;
262                     case ljcrLB: nbp->vdwtype = evdwCuCUTCOMBLB; break;
263                     default:
264                         gmx_incons(
265                                 "The requested LJ combination rule is not implemented in the CUDA "
266                                 "GPU accelerated kernels!");
267                 }
268                 break;
269             case eintmodFORCESWITCH: nbp->vdwtype = evdwCuFSWITCH; break;
270             case eintmodPOTSWITCH: nbp->vdwtype = evdwCuPSWITCH; break;
271             default:
272                 gmx_incons(
273                         "The requested VdW interaction modifier is not implemented in the CUDA GPU "
274                         "accelerated kernels!");
275         }
276     }
277     else if (ic->vdwtype == evdwPME)
278     {
279         if (ic->ljpme_comb_rule == ljcrGEOM)
280         {
281             assert(nbatParams.comb_rule == ljcrGEOM);
282             nbp->vdwtype = evdwCuEWALDGEOM;
283         }
284         else
285         {
286             assert(nbatParams.comb_rule == ljcrLB);
287             nbp->vdwtype = evdwCuEWALDLB;
288         }
289     }
290     else
291     {
292         gmx_incons(
293                 "The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
294     }
295
296     if (ic->eeltype == eelCUT)
297     {
298         nbp->eeltype = eelCuCUT;
299     }
300     else if (EEL_RF(ic->eeltype))
301     {
302         nbp->eeltype = eelCuRF;
303     }
304     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
305     {
306         nbp->eeltype = pick_ewald_kernel_type(*ic);
307     }
308     else
309     {
310         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
311         gmx_incons(
312                 "The requested electrostatics type is not implemented in the CUDA GPU accelerated "
313                 "kernels!");
314     }
315
316     /* generate table for PME */
317     nbp->coulomb_tab = nullptr;
318     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
319     {
320         GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
321         init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
322     }
323
324     /* set up LJ parameter lookup table */
325     if (!useLjCombRule(nbp))
326     {
327         initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj, nbatParams.nbfp.data(), 2 * ntypes * ntypes);
328     }
329
330     /* set up LJ-PME parameter lookup table */
331     if (ic->vdwtype == evdwPME)
332     {
333         initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(), 2 * ntypes);
334     }
335 }
336
337 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
338  *  electrostatic type switching to twin cut-off (or back) if needed. */
339 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
340 {
341     if (!nbv || !nbv->useGpu())
342     {
343         return;
344     }
345     cu_nbparam_t* nbp = nbv->gpu_nbv->nbparam;
346
347     set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
348
349     nbp->eeltype = pick_ewald_kernel_type(*ic);
350
351     GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
352     init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp);
353 }
354
355 /*! Initializes the pair list data structure. */
356 static void init_plist(cu_plist_t* pl)
357 {
358     /* initialize to nullptr pointers to data that is not allocated here and will
359        need reallocation in nbnxn_gpu_init_pairlist */
360     pl->sci   = nullptr;
361     pl->cj4   = nullptr;
362     pl->imask = nullptr;
363     pl->excl  = nullptr;
364
365     /* size -1 indicates that the respective array hasn't been initialized yet */
366     pl->na_c          = -1;
367     pl->nsci          = -1;
368     pl->sci_nalloc    = -1;
369     pl->ncj4          = -1;
370     pl->cj4_nalloc    = -1;
371     pl->nimask        = -1;
372     pl->imask_nalloc  = -1;
373     pl->nexcl         = -1;
374     pl->excl_nalloc   = -1;
375     pl->haveFreshList = false;
376 }
377
378 /*! Initializes the timings data structure. */
379 static void init_timings(gmx_wallclock_gpu_nbnxn_t* t)
380 {
381     int i, j;
382
383     t->nb_h2d_t = 0.0;
384     t->nb_d2h_t = 0.0;
385     t->nb_c     = 0;
386     t->pl_h2d_t = 0.0;
387     t->pl_h2d_c = 0;
388     for (i = 0; i < 2; i++)
389     {
390         for (j = 0; j < 2; j++)
391         {
392             t->ktime[i][j].t = 0.0;
393             t->ktime[i][j].c = 0;
394         }
395     }
396     t->pruneTime.c        = 0;
397     t->pruneTime.t        = 0.0;
398     t->dynamicPruneTime.c = 0;
399     t->dynamicPruneTime.t = 0.0;
400 }
401
402 /*! Initializes simulation constant data. */
403 static void cuda_init_const(NbnxmGpu*                       nb,
404                             const interaction_const_t*      ic,
405                             const PairlistParams&           listParams,
406                             const nbnxn_atomdata_t::Params& nbatParams)
407 {
408     init_atomdata_first(nb->atdat, nbatParams.numTypes);
409     init_nbparam(nb->nbparam, ic, listParams, nbatParams);
410
411     /* clear energy and shift force outputs */
412     nbnxn_cuda_clear_e_fshift(nb);
413 }
414
415 NbnxmGpu* gpu_init(const gmx_device_info_t*   deviceInfo,
416                    const interaction_const_t* ic,
417                    const PairlistParams&      listParams,
418                    const nbnxn_atomdata_t*    nbat,
419                    int /*rank*/,
420                    bool bLocalAndNonlocal)
421 {
422     cudaError_t stat;
423
424     auto nb = new NbnxmGpu;
425     snew(nb->atdat, 1);
426     snew(nb->nbparam, 1);
427     snew(nb->plist[InteractionLocality::Local], 1);
428     if (bLocalAndNonlocal)
429     {
430         snew(nb->plist[InteractionLocality::NonLocal], 1);
431     }
432
433     nb->bUseTwoStreams = bLocalAndNonlocal;
434
435     nb->timers = new cu_timers_t();
436     snew(nb->timings, 1);
437
438     /* init nbst */
439     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
440     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
441     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
442
443     init_plist(nb->plist[InteractionLocality::Local]);
444
445     /* set device info, just point it to the right GPU among the detected ones */
446     nb->dev_info = deviceInfo;
447
448     /* local/non-local GPU streams */
449     stat = cudaStreamCreate(&nb->stream[InteractionLocality::Local]);
450     CU_RET_ERR(stat, "cudaStreamCreate on stream[InterationLocality::Local] failed");
451     if (nb->bUseTwoStreams)
452     {
453         init_plist(nb->plist[InteractionLocality::NonLocal]);
454
455         /* Note that the device we're running on does not have to support
456          * priorities, because we are querying the priority range which in this
457          * case will be a single value.
458          */
459         int highest_priority;
460         stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
461         CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
462
463         stat = cudaStreamCreateWithPriority(&nb->stream[InteractionLocality::NonLocal],
464                                             cudaStreamDefault, highest_priority);
465         CU_RET_ERR(stat,
466                    "cudaStreamCreateWithPriority on stream[InteractionLocality::NonLocal] failed");
467     }
468
469     /* init events for sychronization (timing disabled for performance reasons!) */
470     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
471     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
472     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
473     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
474
475     nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
476
477     /* WARNING: CUDA timings are incorrect with multiple streams.
478      *          This is the main reason why they are disabled by default.
479      */
480     // TODO: Consider turning on by default when we can detect nr of streams.
481     nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
482
483     if (nb->bDoTime)
484     {
485         init_timings(nb->timings);
486     }
487
488     /* set the kernel type for the current GPU */
489     /* pick L1 cache configuration */
490     cuda_set_cacheconfig();
491
492     cuda_init_const(nb, ic, listParams, nbat->params());
493
494     nb->atomIndicesSize       = 0;
495     nb->atomIndicesSize_alloc = 0;
496     nb->ncxy_na               = 0;
497     nb->ncxy_na_alloc         = 0;
498     nb->ncxy_ind              = 0;
499     nb->ncxy_ind_alloc        = 0;
500     nb->ncell                 = 0;
501     nb->ncell_alloc           = 0;
502
503     if (debug)
504     {
505         fprintf(debug, "Initialized CUDA data structures.\n");
506     }
507
508     return nb;
509 }
510
511 void gpu_init_pairlist(NbnxmGpu* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
512 {
513     char         sbuf[STRLEN];
514     bool         bDoTime = (nb->bDoTime && !h_plist->sci.empty());
515     cudaStream_t stream  = nb->stream[iloc];
516     cu_plist_t*  d_plist = nb->plist[iloc];
517
518     if (d_plist->na_c < 0)
519     {
520         d_plist->na_c = h_plist->na_ci;
521     }
522     else
523     {
524         if (d_plist->na_c != h_plist->na_ci)
525         {
526             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
527                     d_plist->na_c, h_plist->na_ci);
528             gmx_incons(sbuf);
529         }
530     }
531
532     gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
533
534     if (bDoTime)
535     {
536         iTimers.pl_h2d.openTimingRegion(stream);
537         iTimers.didPairlistH2D = true;
538     }
539
540     DeviceContext context = nullptr;
541
542     reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc, context);
543     copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), stream,
544                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
545
546     reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc, context);
547     copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), stream,
548                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
549
550     reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
551                            &d_plist->nimask, &d_plist->imask_nalloc, context);
552
553     reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
554                            &d_plist->excl_nalloc, context);
555     copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), stream,
556                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
557
558     if (bDoTime)
559     {
560         iTimers.pl_h2d.closeTimingRegion(stream);
561     }
562
563     /* the next use of thist list we be the first one, so we need to prune */
564     d_plist->haveFreshList = true;
565 }
566
567 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
568 {
569     cu_atomdata_t* adat = nb->atdat;
570     cudaStream_t   ls   = nb->stream[InteractionLocality::Local];
571
572     /* only if we have a dynamic box */
573     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
574     {
575         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(), SHIFTS * sizeof(*adat->shift_vec), ls);
576         adat->bShiftVecUploaded = true;
577     }
578 }
579
580 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
581 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
582 {
583     cudaError_t    stat;
584     cu_atomdata_t* adat = nb->atdat;
585     cudaStream_t   ls   = nb->stream[InteractionLocality::Local];
586
587     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
588     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
589 }
590
591 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
592 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
593 {
594     cudaError_t    stat;
595     cu_atomdata_t* adat = nb->atdat;
596     cudaStream_t   ls   = nb->stream[InteractionLocality::Local];
597
598     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
599     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
600     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
601     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
602     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
603     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
604 }
605
606 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
607 {
608     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
609     /* clear shift force array and energies if the outputs were
610        used in the current step */
611     if (computeVirial)
612     {
613         nbnxn_cuda_clear_e_fshift(nb);
614     }
615 }
616
617 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
618 {
619     cudaError_t    stat;
620     int            nalloc, natoms;
621     bool           realloced;
622     bool           bDoTime = nb->bDoTime;
623     cu_timers_t*   timers  = nb->timers;
624     cu_atomdata_t* d_atdat = nb->atdat;
625     cudaStream_t   ls      = nb->stream[InteractionLocality::Local];
626
627     natoms    = nbat->numAtoms();
628     realloced = false;
629
630     if (bDoTime)
631     {
632         /* time async copy */
633         timers->atdat.openTimingRegion(ls);
634     }
635
636     /* need to reallocate if we have to copy more atoms than the amount of space
637        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
638     if (natoms > d_atdat->nalloc)
639     {
640         nalloc = over_alloc_small(natoms);
641
642         /* free up first if the arrays have already been initialized */
643         if (d_atdat->nalloc != -1)
644         {
645             freeDeviceBuffer(&d_atdat->f);
646             freeDeviceBuffer(&d_atdat->xq);
647             freeDeviceBuffer(&d_atdat->atom_types);
648             freeDeviceBuffer(&d_atdat->lj_comb);
649         }
650
651         stat = cudaMalloc((void**)&d_atdat->f, nalloc * sizeof(*d_atdat->f));
652         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
653         stat = cudaMalloc((void**)&d_atdat->xq, nalloc * sizeof(*d_atdat->xq));
654         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
655         if (useLjCombRule(nb->nbparam))
656         {
657             stat = cudaMalloc((void**)&d_atdat->lj_comb, nalloc * sizeof(*d_atdat->lj_comb));
658             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
659         }
660         else
661         {
662             stat = cudaMalloc((void**)&d_atdat->atom_types, nalloc * sizeof(*d_atdat->atom_types));
663             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
664         }
665
666         d_atdat->nalloc = nalloc;
667         realloced       = true;
668     }
669
670     d_atdat->natoms       = natoms;
671     d_atdat->natoms_local = nbat->natoms_local;
672
673     /* need to clear GPU f output if realloc happened */
674     if (realloced)
675     {
676         nbnxn_cuda_clear_f(nb, nalloc);
677     }
678
679     if (useLjCombRule(nb->nbparam))
680     {
681         cu_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(),
682                           natoms * sizeof(*d_atdat->lj_comb), ls);
683     }
684     else
685     {
686         cu_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(),
687                           natoms * sizeof(*d_atdat->atom_types), ls);
688     }
689
690     if (bDoTime)
691     {
692         timers->atdat.closeTimingRegion(ls);
693     }
694 }
695
696 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t* nbparam)
697 {
698     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
699     {
700         destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
701     }
702 }
703
704 void gpu_free(NbnxmGpu* nb)
705 {
706     cudaError_t    stat;
707     cu_atomdata_t* atdat;
708     cu_nbparam_t*  nbparam;
709
710     if (nb == nullptr)
711     {
712         return;
713     }
714
715     atdat   = nb->atdat;
716     nbparam = nb->nbparam;
717
718     nbnxn_cuda_free_nbparam_table(nbparam);
719
720     stat = cudaEventDestroy(nb->nonlocal_done);
721     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
722     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
723     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
724
725     delete nb->timers;
726     if (nb->bDoTime)
727     {
728         /* The non-local counters/stream (second in the array) are needed only with DD. */
729         for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
730         {
731             stat = cudaStreamDestroy(nb->stream[i]);
732             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
733         }
734     }
735
736     if (!useLjCombRule(nb->nbparam))
737     {
738         destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj);
739     }
740
741     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
742     {
743         destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
744     }
745
746     stat = cudaFree(atdat->shift_vec);
747     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
748     stat = cudaFree(atdat->fshift);
749     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
750
751     stat = cudaFree(atdat->e_lj);
752     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
753     stat = cudaFree(atdat->e_el);
754     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
755
756     freeDeviceBuffer(&atdat->f);
757     freeDeviceBuffer(&atdat->xq);
758     freeDeviceBuffer(&atdat->atom_types);
759     freeDeviceBuffer(&atdat->lj_comb);
760
761     /* Free plist */
762     auto* plist = nb->plist[InteractionLocality::Local];
763     freeDeviceBuffer(&plist->sci);
764     freeDeviceBuffer(&plist->cj4);
765     freeDeviceBuffer(&plist->imask);
766     freeDeviceBuffer(&plist->excl);
767     sfree(plist);
768     if (nb->bUseTwoStreams)
769     {
770         auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
771         freeDeviceBuffer(&plist_nl->sci);
772         freeDeviceBuffer(&plist_nl->cj4);
773         freeDeviceBuffer(&plist_nl->imask);
774         freeDeviceBuffer(&plist_nl->excl);
775         sfree(plist_nl);
776     }
777
778     /* Free nbst */
779     pfree(nb->nbst.e_lj);
780     nb->nbst.e_lj = nullptr;
781
782     pfree(nb->nbst.e_el);
783     nb->nbst.e_el = nullptr;
784
785     pfree(nb->nbst.fshift);
786     nb->nbst.fshift = nullptr;
787
788     sfree(atdat);
789     sfree(nbparam);
790     sfree(nb->timings);
791     delete nb;
792
793     if (debug)
794     {
795         fprintf(debug, "Cleaned up CUDA data structures.\n");
796     }
797 }
798
799 //! This function is documented in the header file
800 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(NbnxmGpu* nb)
801 {
802     return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
803 }
804
805 void gpu_reset_timings(nonbonded_verlet_t* nbv)
806 {
807     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
808     {
809         init_timings(nbv->gpu_nbv->timings);
810     }
811 }
812
813 int gpu_min_ci_balanced(NbnxmGpu* nb)
814 {
815     return nb != nullptr ? gpu_min_ci_balanced_factor * nb->dev_info->prop.multiProcessorCount : 0;
816 }
817
818 gmx_bool gpu_is_kernel_ewald_analytical(const NbnxmGpu* nb)
819 {
820     return ((nb->nbparam->eeltype == eelCuEWALD_ANA) || (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
821 }
822
823 void* gpu_get_command_stream(NbnxmGpu* nb, const InteractionLocality iloc)
824 {
825     assert(nb);
826
827     return static_cast<void*>(&nb->stream[iloc]);
828 }
829
830 void* gpu_get_xq(NbnxmGpu* nb)
831 {
832     assert(nb);
833
834     return static_cast<void*>(nb->atdat->xq);
835 }
836
837 void* gpu_get_f(NbnxmGpu* nb)
838 {
839     assert(nb);
840
841     return static_cast<void*>(nb->atdat->f);
842 }
843
844 rvec* gpu_get_fshift(NbnxmGpu* nb)
845 {
846     assert(nb);
847
848     return reinterpret_cast<rvec*>(nb->atdat->fshift);
849 }
850
851 /* Initialization for X buffer operations on GPU. */
852 /* TODO  Remove explicit pinning from host arrays from here and manage in a more natural way*/
853 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
854 {
855     cudaStream_t stream        = gpu_nbv->stream[InteractionLocality::Local];
856     bool         bDoTime       = gpu_nbv->bDoTime;
857     const int    maxNumColumns = gridSet.numColumnsMax();
858
859     reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
860                            &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, nullptr);
861     reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
862                            &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, nullptr);
863
864     for (unsigned int g = 0; g < gridSet.grids().size(); g++)
865     {
866
867         const Nbnxm::Grid& grid = gridSet.grids()[g];
868
869         const int  numColumns      = grid.numColumns();
870         const int* atomIndices     = gridSet.atomIndices().data();
871         const int  atomIndicesSize = gridSet.atomIndices().size();
872         const int* cxy_na          = grid.cxy_na().data();
873         const int* cxy_ind         = grid.cxy_ind().data();
874
875         reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
876                                &gpu_nbv->atomIndicesSize_alloc, nullptr);
877
878         if (atomIndicesSize > 0)
879         {
880
881             if (bDoTime)
882             {
883                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
884             }
885
886             copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, stream,
887                                GpuApiCallBehavior::Async, nullptr);
888
889             if (bDoTime)
890             {
891                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
892             }
893         }
894
895         if (numColumns > 0)
896         {
897             if (bDoTime)
898             {
899                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
900             }
901
902             int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
903             copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
904
905             if (bDoTime)
906             {
907                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
908             }
909
910             if (bDoTime)
911             {
912                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(stream);
913             }
914
915             destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
916             copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, stream, GpuApiCallBehavior::Async, nullptr);
917
918             if (bDoTime)
919             {
920                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(stream);
921             }
922         }
923     }
924
925     // The above data is transferred on the local stream but is a
926     // dependency of the nonlocal stream (specifically the nonlocal X
927     // buf ops kernel).  We therefore set a dependency to ensure
928     // that the nonlocal stream waits on the local stream here.
929     // This call records an event in the local stream:
930     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
931     // ...and this call instructs the nonlocal stream to wait on that event:
932     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
933
934     return;
935 }
936
937 /* Initialization for F buffer operations on GPU. */
938 void nbnxn_gpu_init_add_nbat_f_to_f(const int*                  cell,
939                                     NbnxmGpu*                   gpu_nbv,
940                                     int                         natoms_total,
941                                     GpuEventSynchronizer* const localReductionDone)
942 {
943
944     cudaStream_t stream = gpu_nbv->stream[InteractionLocality::Local];
945
946     GMX_ASSERT(localReductionDone, "localReductionDone should be a valid pointer");
947     gpu_nbv->localFReductionDone = localReductionDone;
948
949     if (natoms_total > 0)
950     {
951         reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc, nullptr);
952         copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, stream, GpuApiCallBehavior::Async, nullptr);
953     }
954
955     return;
956 }
957
958 } // namespace Nbnxm