Remove unused function cu_synchstream_atdat()
[gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_data_mgmt.cu
blob4cef612eb63ef416bbcf6442529c277ce74bcb0f
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017, 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"
42 #include "config.h"
44 #include <assert.h>
45 #include <stdarg.h>
46 #include <stdio.h>
47 #include <stdlib.h>
49 #include "gromacs/gpu_utils/cudautils.cuh"
50 #include "gromacs/gpu_utils/gpu_utils.h"
51 #include "gromacs/gpu_utils/pmalloc_cuda.h"
52 #include "gromacs/hardware/detecthardware.h"
53 #include "gromacs/hardware/gpu_hw_info.h"
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/mdlib/force_flags.h"
56 #include "gromacs/mdlib/nb_verlet.h"
57 #include "gromacs/mdlib/nbnxn_consts.h"
58 #include "gromacs/mdlib/nbnxn_gpu_data_mgmt.h"
59 #include "gromacs/mdtypes/interaction_const.h"
60 #include "gromacs/mdtypes/md_enums.h"
61 #include "gromacs/pbcutil/ishift.h"
62 #include "gromacs/timing/gpu_timing.h"
63 #include "gromacs/utility/basedefinitions.h"
64 #include "gromacs/utility/cstringutil.h"
65 #include "gromacs/utility/fatalerror.h"
66 #include "gromacs/utility/real.h"
67 #include "gromacs/utility/smalloc.h"
69 #include "nbnxn_cuda_types.h"
71 static bool bUseCudaEventBlockingSync = false; /* makes the CPU thread block */
73 /* This is a heuristically determined parameter for the Fermi, Kepler
74  * and Maxwell architectures for the minimum size of ci lists by multiplying
75  * this constant with the # of multiprocessors on the current device.
76  * Since the maximum number of blocks per multiprocessor is 16, the ideal
77  * count for small systems is 32 or 48 blocks per multiprocessor. Because
78  * there is a bit of fluctuations in the generated block counts, we use
79  * a target of 44 instead of the ideal value of 48.
80  */
81 static unsigned int gpu_min_ci_balanced_factor = 44;
83 /* Functions from nbnxn_cuda.cu */
84 extern void nbnxn_cuda_set_cacheconfig(gmx_device_info_t *devinfo);
85 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_texref();
86 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_nbfp_comb_texref();
87 extern const struct texture<float, 1, cudaReadModeElementType> &nbnxn_cuda_get_coulomb_tab_texref();
90 /* Fw. decl. */
91 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb);
93 /* Fw. decl, */
94 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
95                                           const gmx_device_info_t *dev_info);
98 /*! \brief Return whether texture objects are used on this device.
99  *
100  * \param[in]   pointer to the GPU device info structure to inspect for texture objects support
101  * \return      true if texture objects are used on this device
102  */
103 static bool use_texobj(const gmx_device_info_t *dev_info)
105     /* Only device CC >= 3.0 (Kepler and later) support texture objects */
106     return (dev_info->prop.major >= 3);
109 /*! \brief Return whether combination rules are used.
111  * \param[in]   pointer to nonbonded paramter struct
112  * \return      true if combination rules are used in this run, false otherwise
113  */
114 static inline bool useLjCombRule(const cu_nbparam_t  *nbparam)
116     return (nbparam->vdwtype == evdwCuCUTCOMBGEOM ||
117             nbparam->vdwtype == evdwCuCUTCOMBLB);
120 /*! \brief Set up float texture object.
122  * Set up texture object for float data and bind it to the device memory
123  * \p devPtr points to.
125  * \param[out] texObj   texture object to initialize
126  * \param[in]  devPtr   pointer to device global memory to bind \p texObj to
127  * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
128  */
129 static void setup1DFloatTexture(cudaTextureObject_t &texObj,
130                                 void                *devPtr,
131                                 size_t               sizeInBytes)
133     cudaError_t      stat;
134     cudaResourceDesc rd;
135     cudaTextureDesc  td;
137     memset(&rd, 0, sizeof(rd));
138     rd.resType                  = cudaResourceTypeLinear;
139     rd.res.linear.devPtr        = devPtr;
140     rd.res.linear.desc.f        = cudaChannelFormatKindFloat;
141     rd.res.linear.desc.x        = 32;
142     rd.res.linear.sizeInBytes   = sizeInBytes;
144     memset(&td, 0, sizeof(td));
145     td.readMode                 = cudaReadModeElementType;
146     stat = cudaCreateTextureObject(&texObj, &rd, &td, NULL);
147     CU_RET_ERR(stat, "cudaCreateTextureObject failed");
150 /*! \brief Set up float texture reference.
152  * Set up texture object for float data and bind it to the device memory
153  * \p devPtr points to.
155  * \param[out] texObj   texture reference to initialize
156  * \param[in]  devPtr   pointer to device global memory to bind \p texObj to
157  * \param[in]  sizeInBytes  size of memory area to bind \p texObj to
158  */
159 static void setup1DFloatTexture(const struct texture<float, 1, cudaReadModeElementType> *texRef,
160                                 const void                                              *devPtr,
161                                 size_t                                                   sizeInBytes)
163     cudaError_t           stat;
164     cudaChannelFormatDesc cd;
166     cd   = cudaCreateChannelDesc<float>();
167     stat = cudaBindTexture(NULL, texRef, devPtr, &cd, sizeInBytes);
168     CU_RET_ERR(stat, "cudaBindTexture failed");
172 /*! \brief Initialized the Ewald Coulomb correction GPU table.
174     Tabulates the Ewald Coulomb force and initializes the size/scale
175     and the table GPU array. If called with an already allocated table,
176     it just re-uploads the table.
177  */
178 static void init_ewald_coulomb_force_table(const interaction_const_t *ic,
179                                            cu_nbparam_t              *nbp,
180                                            const gmx_device_info_t   *dev_info)
182     float       *coul_tab;
183     cudaError_t  stat;
185     if (nbp->coulomb_tab != NULL)
186     {
187         nbnxn_cuda_free_nbparam_table(nbp, dev_info);
188     }
190     /* initialize table data in nbp and crete/copy into in global mem */
191     stat = cudaMalloc((void **)&coul_tab, ic->tabq_size*sizeof(*coul_tab));
192     CU_RET_ERR(stat, "cudaMalloc failed on coulumb_tab");
193     cu_copy_H2D(coul_tab, ic->tabq_coul_F, ic->tabq_size*sizeof(*coul_tab));
195     nbp->coulomb_tab       = coul_tab;
196     nbp->coulomb_tab_size  = ic->tabq_size;
197     nbp->coulomb_tab_scale = ic->tabq_scale;
199     if (use_texobj(dev_info))
200     {
201         setup1DFloatTexture(nbp->coulomb_tab_texobj, nbp->coulomb_tab,
202                             nbp->coulomb_tab_size*sizeof(*nbp->coulomb_tab));
203     }
204     else
205     {
206         setup1DFloatTexture(&nbnxn_cuda_get_coulomb_tab_texref(), nbp->coulomb_tab,
207                             nbp->coulomb_tab_size*sizeof(*nbp->coulomb_tab));
208     }
212 /*! Initializes the atomdata structure first time, it only gets filled at
213     pair-search. */
214 static void init_atomdata_first(cu_atomdata_t *ad, int ntypes)
216     cudaError_t stat;
218     ad->ntypes  = ntypes;
219     stat        = cudaMalloc((void**)&ad->shift_vec, SHIFTS*sizeof(*ad->shift_vec));
220     CU_RET_ERR(stat, "cudaMalloc failed on ad->shift_vec");
221     ad->bShiftVecUploaded = false;
223     stat = cudaMalloc((void**)&ad->fshift, SHIFTS*sizeof(*ad->fshift));
224     CU_RET_ERR(stat, "cudaMalloc failed on ad->fshift");
226     stat = cudaMalloc((void**)&ad->e_lj, sizeof(*ad->e_lj));
227     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_lj");
228     stat = cudaMalloc((void**)&ad->e_el, sizeof(*ad->e_el));
229     CU_RET_ERR(stat, "cudaMalloc failed on ad->e_el");
231     /* initialize to NULL poiters to data that is not allocated here and will
232        need reallocation in nbnxn_cuda_init_atomdata */
233     ad->xq = NULL;
234     ad->f  = NULL;
236     /* size -1 indicates that the respective array hasn't been initialized yet */
237     ad->natoms = -1;
238     ad->nalloc = -1;
241 /*! Selects the Ewald kernel type, analytical on SM 3.0 and later, tabulated on
242     earlier GPUs, single or twin cut-off. */
243 static int pick_ewald_kernel_type(bool                     bTwinCut,
244                                   const gmx_device_info_t *dev_info)
246     bool bUseAnalyticalEwald, bForceAnalyticalEwald, bForceTabulatedEwald;
247     int  kernel_type;
249     /* Benchmarking/development environment variables to force the use of
250        analytical or tabulated Ewald kernel. */
251     bForceAnalyticalEwald = (getenv("GMX_CUDA_NB_ANA_EWALD") != NULL);
252     bForceTabulatedEwald  = (getenv("GMX_CUDA_NB_TAB_EWALD") != NULL);
254     if (bForceAnalyticalEwald && bForceTabulatedEwald)
255     {
256         gmx_incons("Both analytical and tabulated Ewald CUDA non-bonded kernels "
257                    "requested through environment variables.");
258     }
260     /* By default, on SM 3.0 and later use analytical Ewald, on earlier tabulated. */
261     if ((dev_info->prop.major >= 3 || bForceAnalyticalEwald) && !bForceTabulatedEwald)
262     {
263         bUseAnalyticalEwald = true;
265         if (debug)
266         {
267             fprintf(debug, "Using analytical Ewald CUDA kernels\n");
268         }
269     }
270     else
271     {
272         bUseAnalyticalEwald = false;
274         if (debug)
275         {
276             fprintf(debug, "Using tabulated Ewald CUDA kernels\n");
277         }
278     }
280     /* Use twin cut-off kernels if requested by bTwinCut or the env. var.
281        forces it (use it for debugging/benchmarking only). */
282     if (!bTwinCut && (getenv("GMX_CUDA_NB_EWALD_TWINCUT") == NULL))
283     {
284         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA : eelCuEWALD_TAB;
285     }
286     else
287     {
288         kernel_type = bUseAnalyticalEwald ? eelCuEWALD_ANA_TWIN : eelCuEWALD_TAB_TWIN;
289     }
291     return kernel_type;
294 /*! Copies all parameters related to the cut-off from ic to nbp */
295 static void set_cutoff_parameters(cu_nbparam_t              *nbp,
296                                   const interaction_const_t *ic)
298     nbp->ewald_beta       = ic->ewaldcoeff_q;
299     nbp->sh_ewald         = ic->sh_ewald;
300     nbp->epsfac           = ic->epsfac;
301     nbp->two_k_rf         = 2.0 * ic->k_rf;
302     nbp->c_rf             = ic->c_rf;
303     nbp->rvdw_sq          = ic->rvdw * ic->rvdw;
304     nbp->rcoulomb_sq      = ic->rcoulomb * ic->rcoulomb;
305     nbp->rlist_sq         = ic->rlist * ic->rlist;
307     nbp->sh_lj_ewald      = ic->sh_lj_ewald;
308     nbp->ewaldcoeff_lj    = ic->ewaldcoeff_lj;
310     nbp->rvdw_switch      = ic->rvdw_switch;
311     nbp->dispersion_shift = ic->dispersion_shift;
312     nbp->repulsion_shift  = ic->repulsion_shift;
313     nbp->vdw_switch       = ic->vdw_switch;
316 /*! \brief Initialize LJ parameter lookup table.
318  * Initializes device memory and copies data from host an binds
319  * a texture to allocated device memory to be used for LJ parameter
320  * lookup.
322  * \param[out] devPtr    device pointer to the memory to be allocated
323  * \param[out] texObj    texture object to be initialized
324  * \param[out] texRef    texture reference to be initialized
325  * \param[in]  hostPtr   pointer to the host memory to be uploaded to the device
326  * \param[in]  numElem   number of elements in the hostPtr
327  * \param[in]  devInfo   pointer to the info struct of the device in use
328  */
329 static void initParamLookupTable(float                    * &devPtr,
330                                  cudaTextureObject_t       &texObj,
331                                  const struct texture<float, 1, cudaReadModeElementType> *texRef,
332                                  const float               *hostPtr,
333                                  int                        numElem,
334                                  const gmx_device_info_t   *devInfo)
336     cudaError_t stat;
338     size_t      sizeInBytes = numElem*sizeof(*devPtr);
340     stat  = cudaMalloc((void **)&devPtr, sizeInBytes);
341     CU_RET_ERR(stat, "cudaMalloc failed in initParamLookupTable");
342     cu_copy_H2D(devPtr, (void *)hostPtr, sizeInBytes);
344     if (use_texobj(devInfo))
345     {
346         setup1DFloatTexture(texObj, devPtr, sizeInBytes);
347     }
348     else
349     {
350         setup1DFloatTexture(texRef, devPtr, sizeInBytes);
351     }
354 /*! Initializes the nonbonded parameter data structure. */
355 static void init_nbparam(cu_nbparam_t              *nbp,
356                          const interaction_const_t *ic,
357                          const nbnxn_atomdata_t    *nbat,
358                          const gmx_device_info_t   *dev_info)
360     int         ntypes;
362     ntypes  = nbat->ntype;
364     set_cutoff_parameters(nbp, ic);
366     /* The kernel code supports LJ combination rules (geometric and LB) for
367      * all kernel types, but we only generate useful combination rule kernels.
368      * We currently only use LJ combination rule (geometric and LB) kernels
369      * for plain cut-off LJ. On Maxwell the force only kernels speed up 15%
370      * with PME and 20% with RF, the other kernels speed up about half as much.
371      * For LJ force-switch the geometric rule would give 7% speed-up, but this
372      * combination is rarely used. LJ force-switch with LB rule is more common,
373      * but gives only 1% speed-up.
374      */
375     if (ic->vdwtype == evdwCUT)
376     {
377         switch (ic->vdw_modifier)
378         {
379             case eintmodNONE:
380             case eintmodPOTSHIFT:
381                 switch (nbat->comb_rule)
382                 {
383                     case ljcrNONE:
384                         nbp->vdwtype = evdwCuCUT;
385                         break;
386                     case ljcrGEOM:
387                         nbp->vdwtype = evdwCuCUTCOMBGEOM;
388                         break;
389                     case ljcrLB:
390                         nbp->vdwtype = evdwCuCUTCOMBLB;
391                         break;
392                     default:
393                         gmx_incons("The requested LJ combination rule is not implemented in the CUDA GPU accelerated kernels!");
394                         break;
395                 }
396                 break;
397             case eintmodFORCESWITCH:
398                 nbp->vdwtype = evdwCuFSWITCH;
399                 break;
400             case eintmodPOTSWITCH:
401                 nbp->vdwtype = evdwCuPSWITCH;
402                 break;
403             default:
404                 gmx_incons("The requested VdW interaction modifier is not implemented in the CUDA GPU accelerated kernels!");
405                 break;
406         }
407     }
408     else if (ic->vdwtype == evdwPME)
409     {
410         if (ic->ljpme_comb_rule == ljcrGEOM)
411         {
412             assert(nbat->comb_rule == ljcrGEOM);
413             nbp->vdwtype = evdwCuEWALDGEOM;
414         }
415         else
416         {
417             assert(nbat->comb_rule == ljcrLB);
418             nbp->vdwtype = evdwCuEWALDLB;
419         }
420     }
421     else
422     {
423         gmx_incons("The requested VdW type is not implemented in the CUDA GPU accelerated kernels!");
424     }
426     if (ic->eeltype == eelCUT)
427     {
428         nbp->eeltype = eelCuCUT;
429     }
430     else if (EEL_RF(ic->eeltype))
431     {
432         nbp->eeltype = eelCuRF;
433     }
434     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
435     {
436         /* Initially rcoulomb == rvdw, so it's surely not twin cut-off. */
437         nbp->eeltype = pick_ewald_kernel_type(false, dev_info);
438     }
439     else
440     {
441         /* Shouldn't happen, as this is checked when choosing Verlet-scheme */
442         gmx_incons("The requested electrostatics type is not implemented in the CUDA GPU accelerated kernels!");
443     }
445     /* generate table for PME */
446     nbp->coulomb_tab = NULL;
447     if (nbp->eeltype == eelCuEWALD_TAB || nbp->eeltype == eelCuEWALD_TAB_TWIN)
448     {
449         init_ewald_coulomb_force_table(ic, nbp, dev_info);
450     }
452     /* set up LJ parameter lookup table */
453     if (!useLjCombRule(nbp))
454     {
455         initParamLookupTable(nbp->nbfp, nbp->nbfp_texobj,
456                              &nbnxn_cuda_get_nbfp_texref(),
457                              nbat->nbfp, 2*ntypes*ntypes, dev_info);
458     }
460     /* set up LJ-PME parameter lookup table */
461     if (ic->vdwtype == evdwPME)
462     {
463         initParamLookupTable(nbp->nbfp_comb, nbp->nbfp_comb_texobj,
464                              &nbnxn_cuda_get_nbfp_comb_texref(),
465                              nbat->nbfp_comb, 2*ntypes, dev_info);
466     }
469 /*! Re-generate the GPU Ewald force table, resets rlist, and update the
470  *  electrostatic type switching to twin cut-off (or back) if needed. */
471 void nbnxn_gpu_pme_loadbal_update_param(const nonbonded_verlet_t    *nbv,
472                                         const interaction_const_t   *ic)
474     if (!nbv || nbv->grp[0].kernel_type != nbnxnk8x8x8_GPU)
475     {
476         return;
477     }
478     gmx_nbnxn_cuda_t *nb    = nbv->gpu_nbv;
479     cu_nbparam_t     *nbp   = nb->nbparam;
481     set_cutoff_parameters(nbp, ic);
483     nbp->eeltype        = pick_ewald_kernel_type(ic->rcoulomb != ic->rvdw,
484                                                  nb->dev_info);
486     init_ewald_coulomb_force_table(ic, nb->nbparam, nb->dev_info);
489 /*! Initializes the pair list data structure. */
490 static void init_plist(cu_plist_t *pl)
492     /* initialize to NULL pointers to data that is not allocated here and will
493        need reallocation in nbnxn_gpu_init_pairlist */
494     pl->sci     = NULL;
495     pl->cj4     = NULL;
496     pl->excl    = NULL;
498     /* size -1 indicates that the respective array hasn't been initialized yet */
499     pl->na_c        = -1;
500     pl->nsci        = -1;
501     pl->sci_nalloc  = -1;
502     pl->ncj4        = -1;
503     pl->cj4_nalloc  = -1;
504     pl->nexcl       = -1;
505     pl->excl_nalloc = -1;
506     pl->bDoPrune    = false;
509 /*! Initializes the timer data structure. */
510 static void init_timers(cu_timers_t *t, bool bUseTwoStreams)
512     cudaError_t stat;
513     int         eventflags = ( bUseCudaEventBlockingSync ? cudaEventBlockingSync : cudaEventDefault );
515     stat = cudaEventCreateWithFlags(&(t->start_atdat), eventflags);
516     CU_RET_ERR(stat, "cudaEventCreate on start_atdat failed");
517     stat = cudaEventCreateWithFlags(&(t->stop_atdat), eventflags);
518     CU_RET_ERR(stat, "cudaEventCreate on stop_atdat failed");
520     /* The non-local counters/stream (second in the array) are needed only with DD. */
521     for (int i = 0; i <= (bUseTwoStreams ? 1 : 0); i++)
522     {
523         stat = cudaEventCreateWithFlags(&(t->start_nb_k[i]), eventflags);
524         CU_RET_ERR(stat, "cudaEventCreate on start_nb_k failed");
525         stat = cudaEventCreateWithFlags(&(t->stop_nb_k[i]), eventflags);
526         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_k failed");
529         stat = cudaEventCreateWithFlags(&(t->start_pl_h2d[i]), eventflags);
530         CU_RET_ERR(stat, "cudaEventCreate on start_pl_h2d failed");
531         stat = cudaEventCreateWithFlags(&(t->stop_pl_h2d[i]), eventflags);
532         CU_RET_ERR(stat, "cudaEventCreate on stop_pl_h2d failed");
534         stat = cudaEventCreateWithFlags(&(t->start_nb_h2d[i]), eventflags);
535         CU_RET_ERR(stat, "cudaEventCreate on start_nb_h2d failed");
536         stat = cudaEventCreateWithFlags(&(t->stop_nb_h2d[i]), eventflags);
537         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_h2d failed");
539         stat = cudaEventCreateWithFlags(&(t->start_nb_d2h[i]), eventflags);
540         CU_RET_ERR(stat, "cudaEventCreate on start_nb_d2h failed");
541         stat = cudaEventCreateWithFlags(&(t->stop_nb_d2h[i]), eventflags);
542         CU_RET_ERR(stat, "cudaEventCreate on stop_nb_d2h failed");
543     }
546 /*! Initializes the timings data structure. */
547 static void init_timings(gmx_wallclock_gpu_t *t)
549     int i, j;
551     t->nb_h2d_t = 0.0;
552     t->nb_d2h_t = 0.0;
553     t->nb_c     = 0;
554     t->pl_h2d_t = 0.0;
555     t->pl_h2d_c = 0;
556     for (i = 0; i < 2; i++)
557     {
558         for (j = 0; j < 2; j++)
559         {
560             t->ktime[i][j].t = 0.0;
561             t->ktime[i][j].c = 0;
562         }
563     }
566 /*! Initializes simulation constant data. */
567 static void nbnxn_cuda_init_const(gmx_nbnxn_cuda_t               *nb,
568                                   const interaction_const_t      *ic,
569                                   const nonbonded_verlet_group_t *nbv_group)
571     init_atomdata_first(nb->atdat, nbv_group[0].nbat->ntype);
572     init_nbparam(nb->nbparam, ic, nbv_group[0].nbat, nb->dev_info);
574     /* clear energy and shift force outputs */
575     nbnxn_cuda_clear_e_fshift(nb);
578 void nbnxn_gpu_init(gmx_nbnxn_cuda_t         **p_nb,
579                     const gmx_gpu_info_t      *gpu_info,
580                     const gmx_gpu_opt_t       *gpu_opt,
581                     const interaction_const_t *ic,
582                     nonbonded_verlet_group_t  *nbv_grp,
583                     int                        my_gpu_index,
584                     int                        /*rank*/,
585                     gmx_bool                   bLocalAndNonlocal)
587     cudaError_t       stat;
588     gmx_nbnxn_cuda_t *nb;
590     assert(gpu_info);
592     if (p_nb == NULL)
593     {
594         return;
595     }
597     snew(nb, 1);
598     snew(nb->atdat, 1);
599     snew(nb->nbparam, 1);
600     snew(nb->plist[eintLocal], 1);
601     if (bLocalAndNonlocal)
602     {
603         snew(nb->plist[eintNonlocal], 1);
604     }
606     nb->bUseTwoStreams = bLocalAndNonlocal;
608     snew(nb->timers, 1);
609     snew(nb->timings, 1);
611     /* init nbst */
612     pmalloc((void**)&nb->nbst.e_lj, sizeof(*nb->nbst.e_lj));
613     pmalloc((void**)&nb->nbst.e_el, sizeof(*nb->nbst.e_el));
614     pmalloc((void**)&nb->nbst.fshift, SHIFTS * sizeof(*nb->nbst.fshift));
616     init_plist(nb->plist[eintLocal]);
618     /* set device info, just point it to the right GPU among the detected ones */
619     nb->dev_info = &gpu_info->gpu_dev[get_gpu_device_id(gpu_info, gpu_opt, my_gpu_index)];
621     /* local/non-local GPU streams */
622     stat = cudaStreamCreate(&nb->stream[eintLocal]);
623     CU_RET_ERR(stat, "cudaStreamCreate on stream[eintLocal] failed");
624     if (nb->bUseTwoStreams)
625     {
626         init_plist(nb->plist[eintNonlocal]);
628         /* CUDA stream priority available in the CUDA RT 5.5 API.
629          * Note that the device we're running on does not have to support
630          * priorities, because we are querying the priority range which in this
631          * case will be a single value.
632          */
633 #if GMX_CUDA_VERSION >= 5050
634         {
635             int highest_priority;
636             stat = cudaDeviceGetStreamPriorityRange(NULL, &highest_priority);
637             CU_RET_ERR(stat, "cudaDeviceGetStreamPriorityRange failed");
639             stat = cudaStreamCreateWithPriority(&nb->stream[eintNonlocal],
640                                                 cudaStreamDefault,
641                                                 highest_priority);
642             CU_RET_ERR(stat, "cudaStreamCreateWithPriority on stream[eintNonlocal] failed");
643         }
644 #else
645         stat = cudaStreamCreate(&nb->stream[eintNonlocal]);
646         CU_RET_ERR(stat, "cudaStreamCreate on stream[eintNonlocal] failed");
647 #endif
648     }
650     /* init events for sychronization (timing disabled for performance reasons!) */
651     stat = cudaEventCreateWithFlags(&nb->nonlocal_done, cudaEventDisableTiming);
652     CU_RET_ERR(stat, "cudaEventCreate on nonlocal_done failed");
653     stat = cudaEventCreateWithFlags(&nb->misc_ops_and_local_H2D_done, cudaEventDisableTiming);
654     CU_RET_ERR(stat, "cudaEventCreate on misc_ops_and_local_H2D_done failed");
656     /* CUDA timing disabled as event timers don't work:
657        - with multiple streams = domain-decomposition;
658        - when turned off by GMX_DISABLE_CUDA_TIMING/GMX_DISABLE_GPU_TIMING.
659      */
660     nb->bDoTime = (!nb->bUseTwoStreams &&
661                    (getenv("GMX_DISABLE_CUDA_TIMING") == NULL) &&
662                    (getenv("GMX_DISABLE_GPU_TIMING") == NULL));
664     if (nb->bDoTime)
665     {
666         init_timers(nb->timers, nb->bUseTwoStreams);
667         init_timings(nb->timings);
668     }
670     /* set the kernel type for the current GPU */
671     /* pick L1 cache configuration */
672     nbnxn_cuda_set_cacheconfig(nb->dev_info);
674     nbnxn_cuda_init_const(nb, ic, nbv_grp);
676     *p_nb = nb;
678     if (debug)
679     {
680         fprintf(debug, "Initialized CUDA data structures.\n");
681     }
684 void nbnxn_gpu_init_pairlist(gmx_nbnxn_cuda_t       *nb,
685                              const nbnxn_pairlist_t *h_plist,
686                              int                     iloc)
688     char          sbuf[STRLEN];
689     cudaError_t   stat;
690     bool          bDoTime    = nb->bDoTime;
691     cudaStream_t  stream     = nb->stream[iloc];
692     cu_plist_t   *d_plist    = nb->plist[iloc];
694     if (d_plist->na_c < 0)
695     {
696         d_plist->na_c = h_plist->na_ci;
697     }
698     else
699     {
700         if (d_plist->na_c != h_plist->na_ci)
701         {
702             sprintf(sbuf, "In cu_init_plist: the #atoms per cell has changed (from %d to %d)",
703                     d_plist->na_c, h_plist->na_ci);
704             gmx_incons(sbuf);
705         }
706     }
708     if (bDoTime)
709     {
710         stat = cudaEventRecord(nb->timers->start_pl_h2d[iloc], stream);
711         CU_RET_ERR(stat, "cudaEventRecord failed");
712     }
714     cu_realloc_buffered((void **)&d_plist->sci, h_plist->sci, sizeof(*d_plist->sci),
715                         &d_plist->nsci, &d_plist->sci_nalloc,
716                         h_plist->nsci,
717                         stream, true);
719     cu_realloc_buffered((void **)&d_plist->cj4, h_plist->cj4, sizeof(*d_plist->cj4),
720                         &d_plist->ncj4, &d_plist->cj4_nalloc,
721                         h_plist->ncj4,
722                         stream, true);
724     cu_realloc_buffered((void **)&d_plist->excl, h_plist->excl, sizeof(*d_plist->excl),
725                         &d_plist->nexcl, &d_plist->excl_nalloc,
726                         h_plist->nexcl,
727                         stream, true);
729     if (bDoTime)
730     {
731         stat = cudaEventRecord(nb->timers->stop_pl_h2d[iloc], stream);
732         CU_RET_ERR(stat, "cudaEventRecord failed");
733     }
735     /* need to prune the pair list during the next step */
736     d_plist->bDoPrune = true;
739 void nbnxn_gpu_upload_shiftvec(gmx_nbnxn_cuda_t       *nb,
740                                const nbnxn_atomdata_t *nbatom)
742     cu_atomdata_t *adat  = nb->atdat;
743     cudaStream_t   ls    = nb->stream[eintLocal];
745     /* only if we have a dynamic box */
746     if (nbatom->bDynamicBox || !adat->bShiftVecUploaded)
747     {
748         cu_copy_H2D_async(adat->shift_vec, nbatom->shift_vec,
749                           SHIFTS * sizeof(*adat->shift_vec), ls);
750         adat->bShiftVecUploaded = true;
751     }
754 /*! Clears the first natoms_clear elements of the GPU nonbonded force output array. */
755 static void nbnxn_cuda_clear_f(gmx_nbnxn_cuda_t *nb, int natoms_clear)
757     cudaError_t    stat;
758     cu_atomdata_t *adat  = nb->atdat;
759     cudaStream_t   ls    = nb->stream[eintLocal];
761     stat = cudaMemsetAsync(adat->f, 0, natoms_clear * sizeof(*adat->f), ls);
762     CU_RET_ERR(stat, "cudaMemsetAsync on f falied");
765 /*! Clears nonbonded shift force output array and energy outputs on the GPU. */
766 static void nbnxn_cuda_clear_e_fshift(gmx_nbnxn_cuda_t *nb)
768     cudaError_t    stat;
769     cu_atomdata_t *adat  = nb->atdat;
770     cudaStream_t   ls    = nb->stream[eintLocal];
772     stat = cudaMemsetAsync(adat->fshift, 0, SHIFTS * sizeof(*adat->fshift), ls);
773     CU_RET_ERR(stat, "cudaMemsetAsync on fshift falied");
774     stat = cudaMemsetAsync(adat->e_lj, 0, sizeof(*adat->e_lj), ls);
775     CU_RET_ERR(stat, "cudaMemsetAsync on e_lj falied");
776     stat = cudaMemsetAsync(adat->e_el, 0, sizeof(*adat->e_el), ls);
777     CU_RET_ERR(stat, "cudaMemsetAsync on e_el falied");
780 void nbnxn_gpu_clear_outputs(gmx_nbnxn_cuda_t *nb, int flags)
782     nbnxn_cuda_clear_f(nb, nb->atdat->natoms);
783     /* clear shift force array and energies if the outputs were
784        used in the current step */
785     if (flags & GMX_FORCE_VIRIAL)
786     {
787         nbnxn_cuda_clear_e_fshift(nb);
788     }
791 void nbnxn_gpu_init_atomdata(gmx_nbnxn_cuda_t              *nb,
792                              const struct nbnxn_atomdata_t *nbat)
794     cudaError_t    stat;
795     int            nalloc, natoms;
796     bool           realloced;
797     bool           bDoTime   = nb->bDoTime;
798     cu_timers_t   *timers    = nb->timers;
799     cu_atomdata_t *d_atdat   = nb->atdat;
800     cudaStream_t   ls        = nb->stream[eintLocal];
802     natoms    = nbat->natoms;
803     realloced = false;
805     if (bDoTime)
806     {
807         /* time async copy */
808         stat = cudaEventRecord(timers->start_atdat, ls);
809         CU_RET_ERR(stat, "cudaEventRecord failed");
810     }
812     /* need to reallocate if we have to copy more atoms than the amount of space
813        available and only allocate if we haven't initialized yet, i.e d_atdat->natoms == -1 */
814     if (natoms > d_atdat->nalloc)
815     {
816         nalloc = over_alloc_small(natoms);
818         /* free up first if the arrays have already been initialized */
819         if (d_atdat->nalloc != -1)
820         {
821             cu_free_buffered(d_atdat->f, &d_atdat->natoms, &d_atdat->nalloc);
822             cu_free_buffered(d_atdat->xq);
823             cu_free_buffered(d_atdat->atom_types);
824             cu_free_buffered(d_atdat->lj_comb);
825         }
827         stat = cudaMalloc((void **)&d_atdat->f, nalloc*sizeof(*d_atdat->f));
828         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->f");
829         stat = cudaMalloc((void **)&d_atdat->xq, nalloc*sizeof(*d_atdat->xq));
830         CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->xq");
831         if (useLjCombRule(nb->nbparam))
832         {
833             stat = cudaMalloc((void **)&d_atdat->lj_comb, nalloc*sizeof(*d_atdat->lj_comb));
834             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->lj_comb");
835         }
836         else
837         {
838             stat = cudaMalloc((void **)&d_atdat->atom_types, nalloc*sizeof(*d_atdat->atom_types));
839             CU_RET_ERR(stat, "cudaMalloc failed on d_atdat->atom_types");
840         }
842         d_atdat->nalloc = nalloc;
843         realloced       = true;
844     }
846     d_atdat->natoms       = natoms;
847     d_atdat->natoms_local = nbat->natoms_local;
849     /* need to clear GPU f output if realloc happened */
850     if (realloced)
851     {
852         nbnxn_cuda_clear_f(nb, nalloc);
853     }
855     if (useLjCombRule(nb->nbparam))
856     {
857         cu_copy_H2D_async(d_atdat->lj_comb, nbat->lj_comb,
858                           natoms*sizeof(*d_atdat->lj_comb), ls);
859     }
860     else
861     {
862         cu_copy_H2D_async(d_atdat->atom_types, nbat->type,
863                           natoms*sizeof(*d_atdat->atom_types), ls);
864     }
866     if (bDoTime)
867     {
868         stat = cudaEventRecord(timers->stop_atdat, ls);
869         CU_RET_ERR(stat, "cudaEventRecord failed");
870     }
873 static void nbnxn_cuda_free_nbparam_table(cu_nbparam_t            *nbparam,
874                                           const gmx_device_info_t *dev_info)
876     cudaError_t stat;
878     if (nbparam->eeltype == eelCuEWALD_TAB || nbparam->eeltype == eelCuEWALD_TAB_TWIN)
879     {
880         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
881         if (use_texobj(dev_info))
882         {
883             stat = cudaDestroyTextureObject(nbparam->coulomb_tab_texobj);
884             CU_RET_ERR(stat, "cudaDestroyTextureObject on coulomb_tab_texobj failed");
885         }
886         else
887         {
888             GMX_UNUSED_VALUE(dev_info);
889             stat = cudaUnbindTexture(nbnxn_cuda_get_coulomb_tab_texref());
890             CU_RET_ERR(stat, "cudaUnbindTexture on coulomb_tab_texref failed");
891         }
892         cu_free_buffered(nbparam->coulomb_tab, &nbparam->coulomb_tab_size);
893     }
896 void nbnxn_gpu_free(gmx_nbnxn_cuda_t *nb)
898     cudaError_t      stat;
899     cu_atomdata_t   *atdat;
900     cu_nbparam_t    *nbparam;
901     cu_plist_t      *plist, *plist_nl;
902     cu_timers_t     *timers;
904     if (nb == NULL)
905     {
906         return;
907     }
909     atdat       = nb->atdat;
910     nbparam     = nb->nbparam;
911     plist       = nb->plist[eintLocal];
912     plist_nl    = nb->plist[eintNonlocal];
913     timers      = nb->timers;
915     nbnxn_cuda_free_nbparam_table(nbparam, nb->dev_info);
917     stat = cudaEventDestroy(nb->nonlocal_done);
918     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->nonlocal_done");
919     stat = cudaEventDestroy(nb->misc_ops_and_local_H2D_done);
920     CU_RET_ERR(stat, "cudaEventDestroy failed on timers->misc_ops_and_local_H2D_done");
922     if (nb->bDoTime)
923     {
924         stat = cudaEventDestroy(timers->start_atdat);
925         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_atdat");
926         stat = cudaEventDestroy(timers->stop_atdat);
927         CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_atdat");
929         /* The non-local counters/stream (second in the array) are needed only with DD. */
930         for (int i = 0; i <= (nb->bUseTwoStreams ? 1 : 0); i++)
931         {
932             stat = cudaEventDestroy(timers->start_nb_k[i]);
933             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_k");
934             stat = cudaEventDestroy(timers->stop_nb_k[i]);
935             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_k");
937             stat = cudaEventDestroy(timers->start_pl_h2d[i]);
938             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_pl_h2d");
939             stat = cudaEventDestroy(timers->stop_pl_h2d[i]);
940             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_pl_h2d");
942             stat = cudaStreamDestroy(nb->stream[i]);
943             CU_RET_ERR(stat, "cudaStreamDestroy failed on stream");
945             stat = cudaEventDestroy(timers->start_nb_h2d[i]);
946             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_h2d");
947             stat = cudaEventDestroy(timers->stop_nb_h2d[i]);
948             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_h2d");
950             stat = cudaEventDestroy(timers->start_nb_d2h[i]);
951             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->start_nb_d2h");
952             stat = cudaEventDestroy(timers->stop_nb_d2h[i]);
953             CU_RET_ERR(stat, "cudaEventDestroy failed on timers->stop_nb_d2h");
954         }
955     }
957     if (!useLjCombRule(nb->nbparam))
958     {
959         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
960         if (use_texobj(nb->dev_info))
961         {
962             stat = cudaDestroyTextureObject(nbparam->nbfp_texobj);
963             CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_texobj failed");
964         }
965         else
966         {
967             stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_texref());
968             CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_texref failed");
969         }
970         cu_free_buffered(nbparam->nbfp);
971     }
973     if (nbparam->vdwtype == evdwCuEWALDGEOM || nbparam->vdwtype == evdwCuEWALDLB)
974     {
975         /* Only device CC >= 3.0 (Kepler and later) support texture objects */
976         if (use_texobj(nb->dev_info))
977         {
978             stat = cudaDestroyTextureObject(nbparam->nbfp_comb_texobj);
979             CU_RET_ERR(stat, "cudaDestroyTextureObject on nbfp_comb_texobj failed");
980         }
981         else
982         {
983             stat = cudaUnbindTexture(nbnxn_cuda_get_nbfp_comb_texref());
984             CU_RET_ERR(stat, "cudaUnbindTexture on nbfp_comb_texref failed");
985         }
986         cu_free_buffered(nbparam->nbfp_comb);
987     }
989     stat = cudaFree(atdat->shift_vec);
990     CU_RET_ERR(stat, "cudaFree failed on atdat->shift_vec");
991     stat = cudaFree(atdat->fshift);
992     CU_RET_ERR(stat, "cudaFree failed on atdat->fshift");
994     stat = cudaFree(atdat->e_lj);
995     CU_RET_ERR(stat, "cudaFree failed on atdat->e_lj");
996     stat = cudaFree(atdat->e_el);
997     CU_RET_ERR(stat, "cudaFree failed on atdat->e_el");
999     cu_free_buffered(atdat->f, &atdat->natoms, &atdat->nalloc);
1000     cu_free_buffered(atdat->xq);
1001     cu_free_buffered(atdat->atom_types, &atdat->ntypes);
1002     cu_free_buffered(atdat->lj_comb);
1004     cu_free_buffered(plist->sci, &plist->nsci, &plist->sci_nalloc);
1005     cu_free_buffered(plist->cj4, &plist->ncj4, &plist->cj4_nalloc);
1006     cu_free_buffered(plist->excl, &plist->nexcl, &plist->excl_nalloc);
1007     if (nb->bUseTwoStreams)
1008     {
1009         cu_free_buffered(plist_nl->sci, &plist_nl->nsci, &plist_nl->sci_nalloc);
1010         cu_free_buffered(plist_nl->cj4, &plist_nl->ncj4, &plist_nl->cj4_nalloc);
1011         cu_free_buffered(plist_nl->excl, &plist_nl->nexcl, &plist->excl_nalloc);
1012     }
1014     sfree(atdat);
1015     sfree(nbparam);
1016     sfree(plist);
1017     if (nb->bUseTwoStreams)
1018     {
1019         sfree(plist_nl);
1020     }
1021     sfree(timers);
1022     sfree(nb->timings);
1023     sfree(nb);
1025     if (debug)
1026     {
1027         fprintf(debug, "Cleaned up CUDA data structures.\n");
1028     }
1031 gmx_wallclock_gpu_t * nbnxn_gpu_get_timings(gmx_nbnxn_cuda_t *nb)
1033     return (nb != NULL && nb->bDoTime) ? nb->timings : NULL;
1036 void nbnxn_gpu_reset_timings(nonbonded_verlet_t* nbv)
1038     if (nbv->gpu_nbv && nbv->gpu_nbv->bDoTime)
1039     {
1040         init_timings(nbv->gpu_nbv->timings);
1041     }
1044 int nbnxn_gpu_min_ci_balanced(gmx_nbnxn_cuda_t *nb)
1046     return nb != NULL ?
1047            gpu_min_ci_balanced_factor*nb->dev_info->prop.multiProcessorCount : 0;
1051 gmx_bool nbnxn_gpu_is_kernel_ewald_analytical(const gmx_nbnxn_cuda_t *nb)
1053     return ((nb->nbparam->eeltype == eelCuEWALD_ANA) ||
1054             (nb->nbparam->eeltype == eelCuEWALD_ANA_TWIN));