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