Tests of restrained listed potentials.
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_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 #include "gromacs/gpu_utils/cudautils.cuh"
48 #include "gromacs/gpu_utils/gpu_utils.h"
49 #include "gromacs/gpu_utils/pmalloc_cuda.h"
50 #include "gromacs/hardware/gpu_hw_info.h"
51 #include "gromacs/math/vectypes.h"
52 #include "gromacs/mdlib/force_flags.h"
53 #include "gromacs/mdlib/nb_verlet.h"
54 #include "gromacs/mdlib/nbnxn_consts.h"
55 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/md_enums.h"
58 #include "gromacs/pbcutil/ishift.h"
59 #include "gromacs/timing/gpu_timing.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/real.h"
64 #include "gromacs/utility/smalloc.h"
65
66 #include "nbnxn_cuda.h"
67 #include "nbnxn_cuda_types.h"
68
69 /* This is a heuristically determined parameter for the Kepler
70  * and Maxwell architectures for the minimum size of ci lists by multiplying
71  * this constant with the # of multiprocessors on the current device.
72  * Since the maximum number of blocks per multiprocessor is 16, the ideal
73  * count for small systems is 32 or 48 blocks per multiprocessor. Because
74  * there is a bit of fluctuations in the generated block counts, we use
75  * a target of 44 instead of the ideal value of 48.
76  */
77 static unsigned int gpu_min_ci_balanced_factor = 44;
78
79 /* Fw. decl. */
80 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
81
82 /* Fw. decl, */
83 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam);
84
85 /*! \brief Return whether combination rules are used.
86  *
87  * \param[in]   pointer to nonbonded paramter struct
88  * \return      true if combination rules are used in this run, false otherwise
89  */
90 static inline bool useLjCombRule(const cu_nbparam_t  *nbparam)
91 {
92     return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
93             nbparam->vdwtype == evdwCuCUTCOMBLB);
94 }
95
96 /*! \brief Initialized the Ewald Coulomb correction GPU table.
97
98     Tabulates the Ewald Coulomb force and initializes the size/scale
99     and the table GPU array. If called with an already allocated table,
100     it just re-uploads the table.
101  */
102 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
103                                            cu_nbparam_t              *nbp)
104 {
105     if (nbp->coulomb_tab != nullptr)
106     {
107         nbnxn_cuda_free_nbparam_table(nbp);
108     }
109
110     nbp->coulomb_tab_scale = ic->tabq_scale;
111     initParamLookupTable(nbp->coulomb_tab, nbp->coulomb_tab_texobj,
112                          ic->tabq_coul_F, ic->tabq_size);
113 }
114
115
116 /*! Initializes the atomdata structure first time, it only gets filled at
117     pair-search. */
118 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
119 {
120     cudaError_t stat;
121
122     ad->ntypes  = ntypes;
123     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
124     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
125     ad->bShiftVecUploaded = false;
126
127     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
128     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
129
130     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
131     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
132     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
133     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
134
135     /* initialize to nullptr poiters to data that is not allocated here and will
136        need reallocation in nbnxn_cuda_init_atomdata */
137     ad->xq = nullptr;
138     ad->f  = nullptr;
139
140     /* size -1 indicates that the respective array hasn't been initialized yet */
141     ad->natoms = -1;
142     ad->nalloc = -1;
143 }
144
145 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
146     earlier GPUs, single or twin cut-off. */
147 static int pick_ewald_kernel_type(bool                     bTwinCut)
148 {
149     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
150     int  kernel_type;
151
152     /* Benchmarking/development environment variables to force the use of
153        analytical or tabulated Ewald kernel. */
154     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != nullptr);
155     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != nullptr);
156
157     if (bForceAnalyticalEwald && bForceTabulatedEwald)
158     {
159         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
160                    "requested through environment variables.");
161     }
162
163     /* By default use analytical Ewald. */
164     bUseAnalyticalEwald = true;
165     if (bForceAnalyticalEwald)
166     {
167         if (debug)
168         {
169             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
170         }
171     }
172     else if (bForceTabulatedEwald)
173     {
174         bUseAnalyticalEwald = false;
175
176         if (debug)
177         {
178             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
179         }
180     }
181
182     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
183        forces it (use it for debugging/benchmarking only). */
184     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == nullptr))
185     {
186         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
187     }
188     else
189     {
190         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
191     }
192
193     return kernel_type;
194 }
195
196 /*! Copies all parameters related to the cut-off from ic to nbp */
197 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
198                                   const interaction_const_t *ic,
199                                   const NbnxnListParameters *listParams)
200 {
201     nbp->ewald_beta        = ic->ewaldcoeff_q;
202     nbp->sh_ewald          = ic->sh_ewald;
203     nbp->epsfac            = ic->epsfac;
204     nbp->two_k_rf          = 2.0 * ic->k_rf;
205     nbp->c_rf              = ic->c_rf;
206     nbp->rvdw_sq           = ic->rvdw * ic->rvdw;
207     nbp->rcoulomb_sq       = ic->rcoulomb * ic->rcoulomb;
208     nbp->rlistOuter_sq     = listParams->rlistOuter * listParams->rlistOuter;
209     nbp->rlistInner_sq     = listParams->rlistInner * listParams->rlistInner;
210     nbp->useDynamicPruning = listParams->useDynamicPruning;
211
212     nbp->sh_lj_ewald       = ic->sh_lj_ewald;
213     nbp->ewaldcoeff_lj     = ic->ewaldcoeff_lj;
214
215     nbp->rvdw_switch       = ic->rvdw_switch;
216     nbp->dispersion_shift  = ic->dispersion_shift;
217     nbp->repulsion_shift   = ic->repulsion_shift;
218     nbp->vdw_switch        = ic->vdw_switch;
219 }
220
221 /*! Initializes the nonbonded parameter data structure. */
222 static void init_nbparam(cu_nbparam_t                   *nbp,
223                          const interaction_const_t      *ic,
224                          const NbnxnListParameters      *listParams,
225                          const nbnxn_atomdata_t::Params &nbatParams)
226 {
227     int         ntypes;
228
229     ntypes  = nbatParams.numTypes;
230
231     set_cutoff_parameters(nbp, ic, listParams);
232
233     /* The kernel code supports LJ combination rules (geometric and LB) for
234      * all kernel types, but we only generate useful combination rule kernels.
235      * We currently only use LJ combination rule (geometric and LB) kernels
236      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
237      * with PME and 20% with RF, the other kernels speed up about half as much.
238      * For LJ force-switch the geometric rule would give 7% speed-up, but this
239      * combination is rarely used. LJ force-switch with LB rule is more common,
240      * but gives only 1% speed-up.
241      */
242     if (ic->vdwtype == evdwCUT)
243     {
244         switch (ic->vdw_modifier)
245         {
246             case eintmodNONE:
247             case eintmodPOTSHIFT:
248                 switch (nbatParams.comb_rule)
249                 {
250                     case ljcrNONE:
251                         nbp->vdwtype = evdwCuCUT;
252                         break;
253                     case ljcrGEOM:
254                         nbp->vdwtype = evdwCuCUTCOMBGEOM;
255                         break;
256                     case ljcrLB:
257                         nbp->vdwtype = evdwCuCUTCOMBLB;
258                         break;
259                     default:
260                         gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
261                 }
262                 break;
263             case eintmodFORCESWITCH:
264                 nbp->vdwtype = evdwCuFSWITCH;
265                 break;
266             case eintmodPOTSWITCH:
267                 nbp->vdwtype = evdwCuPSWITCH;
268                 break;
269             default:
270                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
271         }
272     }
273     else if (ic->vdwtype == evdwPME)
274     {
275         if (ic->ljpme_comb_rule == ljcrGEOM)
276         {
277             assert(nbatParams.comb_rule == ljcrGEOM);
278             nbp->vdwtype = evdwCuEWALDGEOM;
279         }
280         else
281         {
282             assert(nbatParams.comb_rule == ljcrLB);
283             nbp->vdwtype = evdwCuEWALDLB;
284         }
285     }
286     else
287     {
288         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
289     }
290
291     if (ic->eeltype == eelCUT)
292     {
293         nbp->eeltype = eelCuCUT;
294     }
295     else if (EEL_RF(ic->eeltype))
296     {
297         nbp->eeltype = eelCuRF;
298     }
299     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
300     {
301         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
302         nbp->eeltype = pick_ewald_kernel_type(false);
303     }
304     else
305     {
306         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
307         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
308     }
309
310     /* generate table for PME */
311     nbp->coulomb_tab = nullptr;
312     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
313     {
314         init_ewald_coulomb_force_table(ic, nbp);
315     }
316
317     /* set up LJ parameter lookup table */
318     if (!useLjCombRule(nbp))
319     {
320         initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
321                              nbatParams.nbfp.data(), 2*ntypes*ntypes);
322     }
323
324     /* set up LJ-PME parameter lookup table */
325     if (ic->vdwtype == evdwPME)
326     {
327         initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
328                              nbatParams.nbfp_comb.data(), 2*ntypes);
329     }
330 }
331
332 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
333  *  electrostatic type switching to twin cut-off (or back) if needed. */
334 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t    *nbv,
335                                         const interaction_const_t   *ic,
336                                         const NbnxnListParameters   *listParams)
337 {
338     if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
339     {
340         return;
341     }
342     gmx_nbnxn_cuda_t *nb    = nbv->gpu_nbv;
343     cu_nbparam_t     *nbp   = nb->nbparam;
344
345     set_cutoff_parameters(nbp, ic, listParams);
346
347     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw);
348
349     init_ewald_coulomb_force_table(ic, nb->nbparam);
350 }
351
352 /*! Initializes the pair list data structure. */
353 static void init_plist(cu_plist_t *pl)
354 {
355     /* initialize to nullptr pointers to data that is not allocated here and will
356        need reallocation in nbnxn_gpu_init_pairlist */
357     pl->sci      = nullptr;
358     pl->cj4      = nullptr;
359     pl->imask    = nullptr;
360     pl->excl     = nullptr;
361
362     /* size -1 indicates that the respective array hasn't been initialized yet */
363     pl->na_c           = -1;
364     pl->nsci           = -1;
365     pl->sci_nalloc     = -1;
366     pl->ncj4           = -1;
367     pl->cj4_nalloc     = -1;
368     pl->nimask         = -1;
369     pl->imask_nalloc   = -1;
370     pl->nexcl          = -1;
371     pl->excl_nalloc    = -1;
372     pl->haveFreshList  = false;
373 }
374
375 /*! Initializes the timer data structure. */
376 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
377 {
378     /* The non-local counters/stream (second in the array) are needed only with DD. */
379     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
380     {
381         t->didPairlistH2D[i]  = false;
382         t->didPrune[i]        = false;
383         t->didRollingPrune[i] = false;
384     }
385 }
386
387 /*! Initializes the timings data structure. */
388 static void init_timings(gmx_wallclock_gpu_nbnxn_t *t)
389 {
390     int i, j;
391
392     t->nb_h2d_t = 0.0;
393     t->nb_d2h_t = 0.0;
394     t->nb_c     = 0;
395     t->pl_h2d_t = 0.0;
396     t->pl_h2d_c = 0;
397     for (i = 0; i < 2; i++)
398     {
399         for (j = 0; j < 2; j++)
400         {
401             t->ktime[i][j].t = 0.0;
402             t->ktime[i][j].c = 0;
403         }
404     }
405     t->pruneTime.c        = 0;
406     t->pruneTime.t        = 0.0;
407     t->dynamicPruneTime.c = 0;
408     t->dynamicPruneTime.t = 0.0;
409 }
410
411 /*! Initializes simulation constant data. */
412 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t               *nb,
413                                   const interaction_const_t      *ic,
414                                   const NbnxnListParameters      *listParams,
415                                   const nbnxn_atomdata_t::Params &nbatParams)
416 {
417     init_atomdata_first(nb->atdat, nbatParams.numTypes);
418     init_nbparam(nb->nbparam, ic, listParams, nbatParams);
419
420     /* clear energy and shift force outputs */
421     nbnxn_cuda_clear_e_fshift(nb);
422 }
423
424 void nbnxn_gpu_init(gmx_nbnxn_cuda_t         **p_nb,
425                     const gmx_device_info_t   *deviceInfo,
426                     const interaction_const_t *ic,
427                     const NbnxnListParameters *listParams,
428                     const nbnxn_atomdata_t    *nbat,
429                     int                        /*rank*/,
430                     gmx_bool                   bLocalAndNonlocal)
431 {
432     cudaError_t       stat;
433     gmx_nbnxn_cuda_t *nb;
434
435     if (p_nb == nullptr)
436     {
437         return;
438     }
439
440     snew(nb, 1);
441     snew(nb->atdat, 1);
442     snew(nb->nbparam, 1);
443     snew(nb->plist[eintLocal], 1);
444     if (bLocalAndNonlocal)
445     {
446         snew(nb->plist[eintNonlocal], 1);
447     }
448
449     nb->bUseTwoStreams = bLocalAndNonlocal;
450
451     nb->timers = new cu_timers_t();
452     snew(nb->timings, 1);
453
454     /* init nbst */
455     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
456     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
457     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
458
459     init_plist(nb->plist[eintLocal]);
460
461     /* set device info, just point it to the right GPU among the detected ones */
462     nb->dev_info = deviceInfo;
463
464     /* local/non-local GPU streams */
465     stat = cudaStreamCreate(&nb->stream[eintLocal]);
466     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
467     if (nb->bUseTwoStreams)
468     {
469         init_plist(nb->plist[eintNonlocal]);
470
471         /* Note that the device we're running on does not have to support
472          * priorities, because we are querying the priority range which in this
473          * case will be a single value.
474          */
475         int highest_priority;
476         stat = cudaDeviceGetStreamPriorityRange(nullptr, &highest_priority);
477         CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
478
479         stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
480                                             cudaStreamDefault,
481                                             highest_priority);
482         CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
483     }
484
485     /* init events for sychronization (timing disabled for performance reasons!) */
486     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
487     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
488     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
489     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
490
491     /* WARNING: CUDA timings are incorrect with multiple streams.
492      *          This is the main reason why they are disabled by default.
493      */
494     // TODO: Consider turning on by default when we can detect nr of streams.
495     nb->bDoTime = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
496
497     if (nb->bDoTime)
498     {
499         init_timers(nb->timers, nb->bUseTwoStreams);
500         init_timings(nb->timings);
501     }
502
503     /* set the kernel type for the current GPU */
504     /* pick L1 cache configuration */
505     nbnxn_cuda_set_cacheconfig();
506
507     nbnxn_cuda_init_const(nb, ic, listParams, nbat->params());
508
509     *p_nb = nb;
510
511     if (debug)
512     {
513         fprintf(debug, "Initialized CUDA data structures.\n");
514     }
515 }
516
517 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t       *nb,
518                              const NbnxnPairlistGpu *h_plist,
519                              int                     iloc)
520 {
521     char          sbuf[STRLEN];
522     bool          bDoTime    =  (nb->bDoTime && !h_plist->sci.empty());
523     cudaStream_t  stream     = nb->stream[iloc];
524     cu_plist_t   *d_plist    = nb->plist[iloc];
525
526     if (d_plist->na_c < 0)
527     {
528         d_plist->na_c = h_plist->na_ci;
529     }
530     else
531     {
532         if (d_plist->na_c != h_plist->na_ci)
533         {
534             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
535                     d_plist->na_c, h_plist->na_ci);
536             gmx_incons(sbuf);
537         }
538     }
539
540     if (bDoTime)
541     {
542         nb->timers->pl_h2d[iloc].openTimingRegion(stream);
543         nb->timers->didPairlistH2D[iloc] = true;
544     }
545
546     Context context = nullptr;
547
548     reallocateDeviceBuffer(&d_plist->sci, h_plist->sci.size(),
549                            &d_plist->nsci, &d_plist->sci_nalloc, context);
550     copyToDeviceBuffer(&d_plist->sci, h_plist->sci.data(), 0, h_plist->sci.size(),
551                        stream, GpuApiCallBehavior::Async,
552                        bDoTime ? nb->timers->pl_h2d[iloc].fetchNextEvent() : nullptr);
553
554     reallocateDeviceBuffer(&d_plist->cj4, h_plist->cj4.size(),
555                            &d_plist->ncj4, &d_plist->cj4_nalloc, context);
556     copyToDeviceBuffer(&d_plist->cj4, h_plist->cj4.data(), 0, h_plist->cj4.size(),
557                        stream, GpuApiCallBehavior::Async,
558                        bDoTime ? nb->timers->pl_h2d[iloc].fetchNextEvent() : nullptr);
559
560     reallocateDeviceBuffer(&d_plist->imask, h_plist->cj4.size()*c_nbnxnGpuClusterpairSplit,
561                            &d_plist->nimask, &d_plist->imask_nalloc, context);
562
563     reallocateDeviceBuffer(&d_plist->excl, h_plist->excl.size(),
564                            &d_plist->nexcl, &d_plist->excl_nalloc, context);
565     copyToDeviceBuffer(&d_plist->excl, h_plist->excl.data(), 0, h_plist->excl.size(),
566                        stream, GpuApiCallBehavior::Async,
567                        bDoTime ? nb->timers->pl_h2d[iloc].fetchNextEvent() : nullptr);
568
569     if (bDoTime)
570     {
571         nb->timers->pl_h2d[iloc].closeTimingRegion(stream);
572     }
573
574     /* the next use of thist list we be the first one, so we need to prune */
575     d_plist->haveFreshList = true;
576 }
577
578 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t       *nb,
579                                const nbnxn_atomdata_t *nbatom)
580 {
581     cu_atomdata_t *adat  = nb->atdat;
582     cudaStream_t   ls    = nb->stream[eintLocal];
583
584     /* only if we have a dynamic box */
585     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
586     {
587         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec.data(),
588                           SHIFTS * sizeof(*adat->shift_vec), ls);
589         adat->bShiftVecUploaded = true;
590     }
591 }
592
593 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
594 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
595 {
596     cudaError_t    stat;
597     cu_atomdata_t *adat  = nb->atdat;
598     cudaStream_t   ls    = nb->stream[eintLocal];
599
600     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
601     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
602 }
603
604 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
605 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
606 {
607     cudaError_t    stat;
608     cu_atomdata_t *adat  = nb->atdat;
609     cudaStream_t   ls    = nb->stream[eintLocal];
610
611     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
612     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
613     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
614     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
615     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
616     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
617 }
618
619 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
620 {
621     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
622     /* clear shift force array and energies if the outputs were
623        used in the current step */
624     if (flags & GMX_FORCE_VIRIAL)
625     {
626         nbnxn_cuda_clear_e_fshift(nb);
627     }
628 }
629
630 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t              *nb,
631                              const nbnxn_atomdata_t        *nbat)
632 {
633     cudaError_t    stat;
634     int            nalloc, natoms;
635     bool           realloced;
636     bool           bDoTime   = nb->bDoTime;
637     cu_timers_t   *timers    = nb->timers;
638     cu_atomdata_t *d_atdat   = nb->atdat;
639     cudaStream_t   ls        = nb->stream[eintLocal];
640
641     natoms    = nbat->numAtoms();
642     realloced = false;
643
644     if (bDoTime)
645     {
646         /* time async copy */
647         timers->atdat.openTimingRegion(ls);
648     }
649
650     /* need to reallocate if we have to copy more atoms than the amount of space
651        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
652     if (natoms > d_atdat->nalloc)
653     {
654         nalloc = over_alloc_small(natoms);
655
656         /* free up first if the arrays have already been initialized */
657         if (d_atdat->nalloc != -1)
658         {
659             freeDeviceBuffer(&d_atdat->f);
660             freeDeviceBuffer(&d_atdat->xq);
661             freeDeviceBuffer(&d_atdat->atom_types);
662             freeDeviceBuffer(&d_atdat->lj_comb);
663         }
664
665         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
666         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
667         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
668         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
669         if (useLjCombRule(nb->nbparam))
670         {
671             stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
672             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
673         }
674         else
675         {
676             stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
677             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
678         }
679
680         d_atdat->nalloc = nalloc;
681         realloced       = true;
682     }
683
684     d_atdat->natoms       = natoms;
685     d_atdat->natoms_local = nbat->natoms_local;
686
687     /* need to clear GPU f output if realloc happened */
688     if (realloced)
689     {
690         nbnxn_cuda_clear_f(nb, nalloc);
691     }
692
693     if (useLjCombRule(nb->nbparam))
694     {
695         cu_copy_H2D_async(d_atdat->lj_comb, nbat->params().lj_comb.data(),
696                           natoms*sizeof(*d_atdat->lj_comb), ls);
697     }
698     else
699     {
700         cu_copy_H2D_async(d_atdat->atom_types, nbat->params().type.data(),
701                           natoms*sizeof(*d_atdat->atom_types), ls);
702     }
703
704     if (bDoTime)
705     {
706         timers->atdat.closeTimingRegion(ls);
707     }
708 }
709
710 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam)
711 {
712     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
713     {
714         destroyParamLookupTable(nbparam->coulomb_tab, nbparam->coulomb_tab_texobj);
715     }
716 }
717
718 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
719 {
720     cudaError_t      stat;
721     cu_atomdata_t   *atdat;
722     cu_nbparam_t    *nbparam;
723
724     if (nb == nullptr)
725     {
726         return;
727     }
728
729     atdat       = nb->atdat;
730     nbparam     = nb->nbparam;
731
732     nbnxn_cuda_free_nbparam_table(nbparam);
733
734     stat = cudaEventDestroy(nb->nonlocal_done);
735     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
736     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
737     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
738
739     delete nb->timers;
740     if (nb->bDoTime)
741     {
742         /* The non-local counters/stream (second in the array) are needed only with DD. */
743         for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
744         {
745             stat = cudaStreamDestroy(nb->stream[i]);
746             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
747         }
748     }
749
750     if (!useLjCombRule(nb->nbparam))
751     {
752         destroyParamLookupTable(nbparam->nbfp, nbparam->nbfp_texobj);
753
754     }
755
756     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
757     {
758         destroyParamLookupTable(nbparam->nbfp_comb, nbparam->nbfp_comb_texobj);
759     }
760
761     stat = cudaFree(atdat->shift_vec);
762     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
763     stat = cudaFree(atdat->fshift);
764     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
765
766     stat = cudaFree(atdat->e_lj);
767     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
768     stat = cudaFree(atdat->e_el);
769     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
770
771     freeDeviceBuffer(&atdat->f);
772     freeDeviceBuffer(&atdat->xq);
773     freeDeviceBuffer(&atdat->atom_types);
774     freeDeviceBuffer(&atdat->lj_comb);
775
776     /* Free plist */
777     auto *plist = nb->plist[eintLocal];
778     freeDeviceBuffer(&plist->sci);
779     freeDeviceBuffer(&plist->cj4);
780     freeDeviceBuffer(&plist->imask);
781     freeDeviceBuffer(&plist->excl);
782     sfree(plist);
783     if (nb->bUseTwoStreams)
784     {
785         auto *plist_nl = nb->plist[eintNonlocal];
786         freeDeviceBuffer(&plist_nl->sci);
787         freeDeviceBuffer(&plist_nl->cj4);
788         freeDeviceBuffer(&plist_nl->imask);
789         freeDeviceBuffer(&plist_nl->excl);
790         sfree(plist_nl);
791     }
792
793     /* Free nbst */
794     pfree(nb->nbst.e_lj);
795     nb->nbst.e_lj = nullptr;
796
797     pfree(nb->nbst.e_el);
798     nb->nbst.e_el = nullptr;
799
800     pfree(nb->nbst.fshift);
801     nb->nbst.fshift = nullptr;
802
803     sfree(atdat);
804     sfree(nbparam);
805     sfree(nb->timings);
806     sfree(nb);
807
808     if (debug)
809     {
810         fprintf(debug, "Cleaned up CUDA data structures.\n");
811     }
812 }
813
814 //! This function is documented in the header file
815 gmx_wallclock_gpu_nbnxn_t *nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
816 {
817     return (nb != nullptr && nb->bDoTime) ? nb->timings : nullptr;
818 }
819
820 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
821 {
822     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
823     {
824         init_timings(nbv->gpu_nbv->timings);
825     }
826 }
827
828 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
829 {
830     return nb != nullptr ?
831            gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
832
833 }
834
835 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
836 {
837     return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
838             (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));
839 }
840
841 void *nbnxn_gpu_get_command_stream(gmx_nbnxn_gpu_t *nb,
842                                    int              iloc)
843 {
844     assert(nb);
845
846     return static_cast<void *>(&nb->stream[iloc]);
847 }
848
849 void *nbnxn_gpu_get_xq(gmx_nbnxn_gpu_t *nb)
850 {
851     assert(nb);
852
853     return static_cast<void *>(nb->atdat->xq);
854 }
855
856 void *nbnxn_gpu_get_f(gmx_nbnxn_gpu_t *nb)
857 {
858     assert(nb);
859
860     return static_cast<void *>(nb->atdat->f);
861 }
862
863 rvec *nbnxn_gpu_get_fshift(gmx_nbnxn_gpu_t *nb)
864 {
865     assert(nb);
866
867     return reinterpret_cast<rvec *>(nb->atdat->fshift);
868 }