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