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