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