Make cl_nbparam into a struct
[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/nbnxm_gpu_data_mgmt.h"
69 #include "gromacs/nbnxm/pairlistsets.h"
70 #include "gromacs/pbcutil/ishift.h"
71 #include "gromacs/timing/gpu_timing.h"
72 #include "gromacs/utility/basedefinitions.h"
73 #include "gromacs/utility/cstringutil.h"
74 #include "gromacs/utility/fatalerror.h"
75 #include "gromacs/utility/real.h"
76 #include "gromacs/utility/smalloc.h"
77
78 #include "nbnxm_cuda.h"
79
80 namespace Nbnxm
81 {
82
83 /* This is a heuristically determined parameter for the Kepler
84  * and Maxwell architectures for the minimum size of ci lists by multiplying
85  * this constant with the # of multiprocessors on the current device.
86  * Since the maximum number of blocks per multiprocessor is 16, the ideal
87  * count for small systems is 32 or 48 blocks per multiprocessor. Because
88  * there is a bit of fluctuations in the generated block counts, we use
89  * a target of 44 instead of the ideal value of 48.
90  */
91 static unsigned int gpu_min_ci_balanced_factor = 44;
92
93 /* Fw. decl. */
94 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb);
95
96 /*! Initializes the atomdata structure first time, it only gets filled at
97     pair-search. */
98 static void init_atomdata_first(cu_atomdata_t* ad, int ntypes, const DeviceContext& deviceContext)
99 {
100     ad->ntypes = ntypes;
101     allocateDeviceBuffer(&ad->shift_vec, SHIFTS, deviceContext);
102     ad->bShiftVecUploaded = false;
103
104     allocateDeviceBuffer(&ad->fshift, SHIFTS, deviceContext);
105     allocateDeviceBuffer(&ad->e_lj, 1, deviceContext);
106     allocateDeviceBuffer(&ad->e_el, 1, deviceContext);
107
108     /* initialize to nullptr poiters to data that is not allocated here and will
109        need reallocation in nbnxn_cuda_init_atomdata */
110     ad->xq = nullptr;
111     ad->f  = nullptr;
112
113     /* size -1 indicates that the respective array hasn't been initialized yet */
114     ad->natoms = -1;
115     ad->nalloc = -1;
116 }
117
118 /*! Initializes the nonbonded parameter data structure. */
119 static void init_nbparam(NBParamGpu*                     nbp,
120                          const interaction_const_t*      ic,
121                          const PairlistParams&           listParams,
122                          const nbnxn_atomdata_t::Params& nbatParams,
123                          const DeviceContext&            deviceContext)
124 {
125     int ntypes;
126
127     ntypes = nbatParams.numTypes;
128
129     set_cutoff_parameters(nbp, ic, listParams);
130
131     /* The kernel code supports LJ combination rules (geometric and LB) for
132      * all kernel types, but we only generate useful combination rule kernels.
133      * We currently only use LJ combination rule (geometric and LB) kernels
134      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
135      * with PME and 20% with RF, the other kernels speed up about half as much.
136      * For LJ force-switch the geometric rule would give 7% speed-up, but this
137      * combination is rarely used. LJ force-switch with LB rule is more common,
138      * but gives only 1% speed-up.
139      */
140     if (ic->vdwtype == evdwCUT)
141     {
142         switch (ic->vdw_modifier)
143         {
144             case eintmodNONE:
145             case eintmodPOTSHIFT:
146                 switch (nbatParams.comb_rule)
147                 {
148                     case ljcrNONE: nbp->vdwtype = evdwTypeCUT; break;
149                     case ljcrGEOM: nbp->vdwtype = evdwTypeCUTCOMBGEOM; break;
150                     case ljcrLB: nbp->vdwtype = evdwTypeCUTCOMBLB; break;
151                     default:
152                         gmx_incons(
153                                 "The requested LJ combination rule is not implemented in the CUDA "
154                                 "GPU accelerated kernels!");
155                 }
156                 break;
157             case eintmodFORCESWITCH: nbp->vdwtype = evdwTypeFSWITCH; break;
158             case eintmodPOTSWITCH: nbp->vdwtype = evdwTypePSWITCH; break;
159             default:
160                 gmx_incons(
161                         "The requested VdW interaction modifier is not implemented in the CUDA GPU "
162                         "accelerated kernels!");
163         }
164     }
165     else if (ic->vdwtype == evdwPME)
166     {
167         if (ic->ljpme_comb_rule == ljcrGEOM)
168         {
169             assert(nbatParams.comb_rule == ljcrGEOM);
170             nbp->vdwtype = evdwTypeEWALDGEOM;
171         }
172         else
173         {
174             assert(nbatParams.comb_rule == ljcrLB);
175             nbp->vdwtype = evdwTypeEWALDLB;
176         }
177     }
178     else
179     {
180         gmx_incons(
181                 "The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
182     }
183
184     if (ic->eeltype == eelCUT)
185     {
186         nbp->eeltype = eelTypeCUT;
187     }
188     else if (EEL_RF(ic->eeltype))
189     {
190         nbp->eeltype = eelTypeRF;
191     }
192     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
193     {
194         nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
195     }
196     else
197     {
198         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
199         gmx_incons(
200                 "The requested electrostatics type is not implemented in the CUDA GPU accelerated "
201                 "kernels!");
202     }
203
204     /* generate table for PME */
205     nbp->coulomb_tab = nullptr;
206     if (nbp->eeltype == eelTypeEWALD_TAB || nbp->eeltype == eelTypeEWALD_TAB_TWIN)
207     {
208         GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
209         init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, deviceContext);
210     }
211
212     /* set up LJ parameter lookup table */
213     if (!useLjCombRule(nbp->vdwtype))
214     {
215         initParamLookupTable(&nbp->nbfp, &nbp->nbfp_texobj, nbatParams.nbfp.data(),
216                              2 * ntypes * ntypes, deviceContext);
217     }
218
219     /* set up LJ-PME parameter lookup table */
220     if (ic->vdwtype == evdwPME)
221     {
222         initParamLookupTable(&nbp->nbfp_comb, &nbp->nbfp_comb_texobj, nbatParams.nbfp_comb.data(),
223                              2 * ntypes, deviceContext);
224     }
225 }
226
227 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
228  *  electrostatic type switching to twin cut-off (or back) if needed. */
229 void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interaction_const_t* ic)
230 {
231     if (!nbv || !nbv->useGpu())
232     {
233         return;
234     }
235     NbnxmGpu*   nb  = nbv->gpu_nbv;
236     NBParamGpu* nbp = nbv->gpu_nbv->nbparam;
237
238     set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
239
240     nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
241
242     GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
243     init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, *nb->deviceContext_);
244 }
245
246 /*! Initializes simulation constant data. */
247 static void cuda_init_const(NbnxmGpu*                       nb,
248                             const interaction_const_t*      ic,
249                             const PairlistParams&           listParams,
250                             const nbnxn_atomdata_t::Params& nbatParams)
251 {
252     init_atomdata_first(nb->atdat, nbatParams.numTypes, *nb->deviceContext_);
253     init_nbparam(nb->nbparam, ic, listParams, nbatParams, *nb->deviceContext_);
254
255     /* clear energy and shift force outputs */
256     nbnxn_cuda_clear_e_fshift(nb);
257 }
258
259 NbnxmGpu* gpu_init(const gmx::DeviceStreamManager& deviceStreamManager,
260                    const interaction_const_t*      ic,
261                    const PairlistParams&           listParams,
262                    const nbnxn_atomdata_t*         nbat,
263                    bool                            bLocalAndNonlocal)
264 {
265     cudaError_t stat;
266
267     auto nb            = new NbnxmGpu();
268     nb->deviceContext_ = &deviceStreamManager.context();
269     snew(nb->atdat, 1);
270     snew(nb->nbparam, 1);
271     snew(nb->plist[InteractionLocality::Local], 1);
272     if (bLocalAndNonlocal)
273     {
274         snew(nb->plist[InteractionLocality::NonLocal], 1);
275     }
276
277     nb->bUseTwoStreams = bLocalAndNonlocal;
278
279     nb->timers = new cu_timers_t();
280     snew(nb->timings, 1);
281
282     /* init nbst */
283     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
284     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
285     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
286
287     init_plist(nb->plist[InteractionLocality::Local]);
288
289     /* local/non-local GPU streams */
290     GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedLocal),
291                        "Local non-bonded stream should be initialized to use GPU for non-bonded.");
292     nb->deviceStreams[InteractionLocality::Local] =
293             &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedLocal);
294     if (nb->bUseTwoStreams)
295     {
296         init_plist(nb->plist[InteractionLocality::NonLocal]);
297
298         /* Note that the device we're running on does not have to support
299          * priorities, because we are querying the priority range which in this
300          * case will be a single value.
301          */
302         GMX_RELEASE_ASSERT(deviceStreamManager.streamIsValid(gmx::DeviceStreamType::NonBondedNonLocal),
303                            "Non-local non-bonded stream should be initialized to use GPU for "
304                            "non-bonded with domain decomposition.");
305         nb->deviceStreams[InteractionLocality::NonLocal] =
306                 &deviceStreamManager.stream(gmx::DeviceStreamType::NonBondedNonLocal);
307         ;
308     }
309
310     /* init events for sychronization (timing disabled for performance reasons!) */
311     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
312     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
313     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
314     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
315
316     nb->xNonLocalCopyD2HDone = new GpuEventSynchronizer();
317
318     /* WARNING: CUDA timings are incorrect with multiple streams.
319      *          This is the main reason why they are disabled by default.
320      */
321     // TODO: Consider turning on by default when we can detect nr of streams.
322     nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
323
324     if (nb->bDoTime)
325     {
326         init_timings(nb->timings);
327     }
328
329     /* set the kernel type for the current GPU */
330     /* pick L1 cache configuration */
331     cuda_set_cacheconfig();
332
333     cuda_init_const(nb, ic, listParams, nbat->params());
334
335     nb->atomIndicesSize       = 0;
336     nb->atomIndicesSize_alloc = 0;
337     nb->ncxy_na               = 0;
338     nb->ncxy_na_alloc         = 0;
339     nb->ncxy_ind              = 0;
340     nb->ncxy_ind_alloc        = 0;
341     nb->ncell                 = 0;
342     nb->ncell_alloc           = 0;
343
344     if (debug)
345     {
346         fprintf(debug, "Initialized CUDA data structures.\n");
347     }
348
349     return nb;
350 }
351
352 void gpu_init_pairlist(NbnxmGpu* nb, const NbnxnPairlistGpu* h_plist, const InteractionLocality iloc)
353 {
354     char                sbuf[STRLEN];
355     bool                bDoTime      = (nb->bDoTime && !h_plist->sci.empty());
356     const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
357     gpu_plist*          d_plist      = nb->plist[iloc];
358
359     if (d_plist->na_c < 0)
360     {
361         d_plist->na_c = h_plist->na_ci;
362     }
363     else
364     {
365         if (d_plist->na_c != h_plist->na_ci)
366         {
367             sprintf(sbuf, "In init_plist: the #atoms per cell has changed (from %d to %d)",
368                     d_plist->na_c, h_plist->na_ci);
369             gmx_incons(sbuf);
370         }
371     }
372
373     gpu_timers_t::Interaction& iTimers = nb->timers->interaction[iloc];
374
375     if (bDoTime)
376     {
377         iTimers.pl_h2d.openTimingRegion(deviceStream);
378         iTimers.didPairlistH2D = true;
379     }
380
381     const DeviceContext& deviceContext = *nb->deviceContext_;
382
383     reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(), &d_plist->nsci, &d_plist->sci_nalloc,
384                            deviceContext);
385     copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(), deviceStream,
386                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
387
388     reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(), &d_plist->ncj4, &d_plist->cj4_nalloc,
389                            deviceContext);
390     copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(), deviceStream,
391                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
392
393     reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size() * c_nbnxnGpuClusterpairSplit,
394                            &d_plist->nimask, &d_plist->imask_nalloc, deviceContext);
395
396     reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(), &d_plist->nexcl,
397                            &d_plist->excl_nalloc, deviceContext);
398     copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(), deviceStream,
399                        GpuApiCallBehavior::Async, bDoTime ? iTimers.pl_h2d.fetchNextEvent() : nullptr);
400
401     if (bDoTime)
402     {
403         iTimers.pl_h2d.closeTimingRegion(deviceStream);
404     }
405
406     /* the next use of thist list we be the first one, so we need to prune */
407     d_plist->haveFreshList = true;
408 }
409
410 void gpu_upload_shiftvec(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom)
411 {
412     cu_atomdata_t*      adat        = nb->atdat;
413     const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
414
415     /* only if we have a dynamic box */
416     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
417     {
418         static_assert(sizeof(adat->shift_vec[0]) == sizeof(nbatom->shift_vec[0]),
419                       "Sizes of host- and device-side shift vectors should be the same.");
420         copyToDeviceBuffer(&adat->shift_vec, reinterpret_cast<const float3*>(nbatom->shift_vec.data()),
421                            0, SHIFTS, localStream, GpuApiCallBehavior::Async, nullptr);
422         adat->bShiftVecUploaded = true;
423     }
424 }
425
426 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
427 static void nbnxn_cuda_clear_f(NbnxmGpu* nb, int natoms_clear)
428 {
429     cu_atomdata_t*      adat        = nb->atdat;
430     const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
431     clearDeviceBufferAsync(&adat->f, 0, natoms_clear, localStream);
432 }
433
434 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
435 static void nbnxn_cuda_clear_e_fshift(NbnxmGpu* nb)
436 {
437     cu_atomdata_t*      adat        = nb->atdat;
438     const DeviceStream& localStream = *nb->deviceStreams[InteractionLocality::Local];
439
440     clearDeviceBufferAsync(&adat->fshift, 0, SHIFTS, localStream);
441     clearDeviceBufferAsync(&adat->e_lj, 0, 1, localStream);
442     clearDeviceBufferAsync(&adat->e_el, 0, 1, localStream);
443 }
444
445 void gpu_clear_outputs(NbnxmGpu* nb, bool computeVirial)
446 {
447     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
448     /* clear shift force array and energies if the outputs were
449        used in the current step */
450     if (computeVirial)
451     {
452         nbnxn_cuda_clear_e_fshift(nb);
453     }
454 }
455
456 void gpu_init_atomdata(NbnxmGpu* nb, const nbnxn_atomdata_t* nbat)
457 {
458     int                  nalloc, natoms;
459     bool                 realloced;
460     bool                 bDoTime       = nb->bDoTime;
461     cu_timers_t*         timers        = nb->timers;
462     cu_atomdata_t*       d_atdat       = nb->atdat;
463     const DeviceContext& deviceContext = *nb->deviceContext_;
464     const DeviceStream&  localStream   = *nb->deviceStreams[InteractionLocality::Local];
465
466     natoms    = nbat->numAtoms();
467     realloced = false;
468
469     if (bDoTime)
470     {
471         /* time async copy */
472         timers->atdat.openTimingRegion(localStream);
473     }
474
475     /* need to reallocate if we have to copy more atoms than the amount of space
476        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
477     if (natoms > d_atdat->nalloc)
478     {
479         nalloc = over_alloc_small(natoms);
480
481         /* free up first if the arrays have already been initialized */
482         if (d_atdat->nalloc != -1)
483         {
484             freeDeviceBuffer(&d_atdat->f);
485             freeDeviceBuffer(&d_atdat->xq);
486             freeDeviceBuffer(&d_atdat->atom_types);
487             freeDeviceBuffer(&d_atdat->lj_comb);
488         }
489
490         allocateDeviceBuffer(&d_atdat->f, nalloc, deviceContext);
491         allocateDeviceBuffer(&d_atdat->xq, nalloc, deviceContext);
492         if (useLjCombRule(nb->nbparam->vdwtype))
493         {
494             allocateDeviceBuffer(&d_atdat->lj_comb, nalloc, deviceContext);
495         }
496         else
497         {
498             allocateDeviceBuffer(&d_atdat->atom_types, nalloc, deviceContext);
499         }
500
501         d_atdat->nalloc = nalloc;
502         realloced       = true;
503     }
504
505     d_atdat->natoms       = natoms;
506     d_atdat->natoms_local = nbat->natoms_local;
507
508     /* need to clear GPU f output if realloc happened */
509     if (realloced)
510     {
511         nbnxn_cuda_clear_f(nb, nalloc);
512     }
513
514     if (useLjCombRule(nb->nbparam->vdwtype))
515     {
516         static_assert(sizeof(d_atdat->lj_comb[0]) == sizeof(float2),
517                       "Size of the LJ parameters element should be equal to the size of float2.");
518         copyToDeviceBuffer(&d_atdat->lj_comb,
519                            reinterpret_cast<const float2*>(nbat->params().lj_comb.data()), 0,
520                            natoms, localStream, GpuApiCallBehavior::Async, nullptr);
521     }
522     else
523     {
524         static_assert(sizeof(d_atdat->atom_types[0]) == sizeof(nbat->params().type[0]),
525                       "Sizes of host- and device-side atom types should be the same.");
526         copyToDeviceBuffer(&d_atdat->atom_types, nbat->params().type.data(), 0, natoms, localStream,
527                            GpuApiCallBehavior::Async, nullptr);
528     }
529
530     if (bDoTime)
531     {
532         timers->atdat.closeTimingRegion(localStream);
533     }
534 }
535
536 void gpu_free(NbnxmGpu* nb)
537 {
538     cudaError_t    stat;
539     cu_atomdata_t* atdat;
540     NBParamGpu*    nbparam;
541
542     if (nb == nullptr)
543     {
544         return;
545     }
546
547     atdat   = nb->atdat;
548     nbparam = nb->nbparam;
549
550     if ((!nbparam->coulomb_tab)
551         && (nbparam->eeltype == eelTypeEWALD_TAB || nbparam->eeltype == eelTypeEWALD_TAB_TWIN))
552     {
553         destroyParamLookupTable(&nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
554     }
555
556     stat = cudaEventDestroy(nb->nonlocal_done);
557     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
558     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
559     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
560
561     delete nb->timers;
562
563     if (!useLjCombRule(nb->nbparam->vdwtype))
564     {
565         destroyParamLookupTable(&nbparam->nbfp, nbparam->nbfp_texobj);
566     }
567
568     if (nbparam->vdwtype == evdwTypeEWALDGEOM || nbparam->vdwtype == evdwTypeEWALDLB)
569     {
570         destroyParamLookupTable(&nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
571     }
572
573     freeDeviceBuffer(&atdat->shift_vec);
574     freeDeviceBuffer(&atdat->fshift);
575
576     freeDeviceBuffer(&atdat->e_lj);
577     freeDeviceBuffer(&atdat->e_el);
578
579     freeDeviceBuffer(&atdat->f);
580     freeDeviceBuffer(&atdat->xq);
581     freeDeviceBuffer(&atdat->atom_types);
582     freeDeviceBuffer(&atdat->lj_comb);
583
584     /* Free plist */
585     auto* plist = nb->plist[InteractionLocality::Local];
586     freeDeviceBuffer(&plist->sci);
587     freeDeviceBuffer(&plist->cj4);
588     freeDeviceBuffer(&plist->imask);
589     freeDeviceBuffer(&plist->excl);
590     sfree(plist);
591     if (nb->bUseTwoStreams)
592     {
593         auto* plist_nl = nb->plist[InteractionLocality::NonLocal];
594         freeDeviceBuffer(&plist_nl->sci);
595         freeDeviceBuffer(&plist_nl->cj4);
596         freeDeviceBuffer(&plist_nl->imask);
597         freeDeviceBuffer(&plist_nl->excl);
598         sfree(plist_nl);
599     }
600
601     /* Free nbst */
602     pfree(nb->nbst.e_lj);
603     nb->nbst.e_lj = nullptr;
604
605     pfree(nb->nbst.e_el);
606     nb->nbst.e_el = nullptr;
607
608     pfree(nb->nbst.fshift);
609     nb->nbst.fshift = nullptr;
610
611     sfree(atdat);
612     sfree(nbparam);
613     sfree(nb->timings);
614     delete nb;
615
616     if (debug)
617     {
618         fprintf(debug, "Cleaned up CUDA data structures.\n");
619     }
620 }
621
622 //! This function is documented in the header file
623 gmx_wallclock_gpu_nbnxn_t* gpu_get_timings(NbnxmGpu* nb)
624 {
625     return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
626 }
627
628 void gpu_reset_timings(nonbonded_verlet_t* nbv)
629 {
630     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
631     {
632         init_timings(nbv->gpu_nbv->timings);
633     }
634 }
635
636 int gpu_min_ci_balanced(NbnxmGpu* nb)
637 {
638     return nb != nullptr ? gpu_min_ci_balanced_factor * nb->deviceContext_->deviceInfo().prop.multiProcessorCount
639                          : 0;
640 }
641
642 gmx_bool gpu_is_kernel_ewald_analytical(const NbnxmGpu* nb)
643 {
644     return ((nb->nbparam->eeltype == eelTypeEWALD_ANA) || (nb->nbparam->eeltype == eelTypeEWALD_ANA_TWIN));
645 }
646
647 void* gpu_get_xq(NbnxmGpu* nb)
648 {
649     assert(nb);
650
651     return static_cast<void*>(nb->atdat->xq);
652 }
653
654 DeviceBuffer<gmx::RVec> gpu_get_f(NbnxmGpu* nb)
655 {
656     assert(nb);
657
658     return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->f);
659 }
660
661 DeviceBuffer<gmx::RVec> gpu_get_fshift(NbnxmGpu* nb)
662 {
663     assert(nb);
664
665     return reinterpret_cast<DeviceBuffer<gmx::RVec>>(nb->atdat->fshift);
666 }
667
668 /* Initialization for X buffer operations on GPU. */
669 /* TODO  Remove explicit pinning from host arrays from here and manage in a more natural way*/
670 void nbnxn_gpu_init_x_to_nbat_x(const Nbnxm::GridSet& gridSet, NbnxmGpu* gpu_nbv)
671 {
672     const DeviceStream& deviceStream  = *gpu_nbv->deviceStreams[InteractionLocality::Local];
673     bool                bDoTime       = gpu_nbv->bDoTime;
674     const int           maxNumColumns = gridSet.numColumnsMax();
675
676     reallocateDeviceBuffer(&gpu_nbv->cxy_na, maxNumColumns * gridSet.grids().size(),
677                            &gpu_nbv->ncxy_na, &gpu_nbv->ncxy_na_alloc, *gpu_nbv->deviceContext_);
678     reallocateDeviceBuffer(&gpu_nbv->cxy_ind, maxNumColumns * gridSet.grids().size(),
679                            &gpu_nbv->ncxy_ind, &gpu_nbv->ncxy_ind_alloc, *gpu_nbv->deviceContext_);
680
681     for (unsigned int g = 0; g < gridSet.grids().size(); g++)
682     {
683
684         const Nbnxm::Grid& grid = gridSet.grids()[g];
685
686         const int  numColumns      = grid.numColumns();
687         const int* atomIndices     = gridSet.atomIndices().data();
688         const int  atomIndicesSize = gridSet.atomIndices().size();
689         const int* cxy_na          = grid.cxy_na().data();
690         const int* cxy_ind         = grid.cxy_ind().data();
691
692         reallocateDeviceBuffer(&gpu_nbv->atomIndices, atomIndicesSize, &gpu_nbv->atomIndicesSize,
693                                &gpu_nbv->atomIndicesSize_alloc, *gpu_nbv->deviceContext_);
694
695         if (atomIndicesSize > 0)
696         {
697
698             if (bDoTime)
699             {
700                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
701             }
702
703             copyToDeviceBuffer(&gpu_nbv->atomIndices, atomIndices, 0, atomIndicesSize, deviceStream,
704                                GpuApiCallBehavior::Async, nullptr);
705
706             if (bDoTime)
707             {
708                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
709             }
710         }
711
712         if (numColumns > 0)
713         {
714             if (bDoTime)
715             {
716                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
717             }
718
719             int* destPtr = &gpu_nbv->cxy_na[maxNumColumns * g];
720             copyToDeviceBuffer(&destPtr, cxy_na, 0, numColumns, deviceStream,
721                                GpuApiCallBehavior::Async, nullptr);
722
723             if (bDoTime)
724             {
725                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
726             }
727
728             if (bDoTime)
729             {
730                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.openTimingRegion(deviceStream);
731             }
732
733             destPtr = &gpu_nbv->cxy_ind[maxNumColumns * g];
734             copyToDeviceBuffer(&destPtr, cxy_ind, 0, numColumns, deviceStream,
735                                GpuApiCallBehavior::Async, nullptr);
736
737             if (bDoTime)
738             {
739                 gpu_nbv->timers->xf[AtomLocality::Local].nb_h2d.closeTimingRegion(deviceStream);
740             }
741         }
742     }
743
744     // The above data is transferred on the local stream but is a
745     // dependency of the nonlocal stream (specifically the nonlocal X
746     // buf ops kernel).  We therefore set a dependency to ensure
747     // that the nonlocal stream waits on the local stream here.
748     // This call records an event in the local stream:
749     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::Local);
750     // ...and this call instructs the nonlocal stream to wait on that event:
751     nbnxnInsertNonlocalGpuDependency(gpu_nbv, Nbnxm::InteractionLocality::NonLocal);
752
753     return;
754 }
755
756 /* Initialization for F buffer operations on GPU. */
757 void nbnxn_gpu_init_add_nbat_f_to_f(const int*                  cell,
758                                     NbnxmGpu*                   gpu_nbv,
759                                     int                         natoms_total,
760                                     GpuEventSynchronizer* const localReductionDone)
761 {
762
763     const DeviceStream& deviceStream = *gpu_nbv->deviceStreams[InteractionLocality::Local];
764
765     GMX_ASSERT(localReductionDone, "localReductionDone should be a valid pointer");
766     gpu_nbv->localFReductionDone = localReductionDone;
767
768     if (natoms_total > 0)
769     {
770         reallocateDeviceBuffer(&gpu_nbv->cell, natoms_total, &gpu_nbv->ncell, &gpu_nbv->ncell_alloc,
771                                *gpu_nbv->deviceContext_);
772         copyToDeviceBuffer(&gpu_nbv->cell, cell, 0, natoms_total, deviceStream,
773                            GpuApiCallBehavior::Async, nullptr);
774     }
775
776     return;
777 }
778
779 } // namespace Nbnxm