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