2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2012, The GROMACS development team.
6 * Copyright (c) 2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
40 * Data types used internally in the nbnxn_cuda module.
42 * \author Szilárd Páll <pall.szilard@gmail.com>
43 * \ingroup module_mdlib
46 #ifndef NBNXN_CUDA_TYPES_H
47 #define NBNXN_CUDA_TYPES_H
49 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/gpu_utils/gpuregiontimer.cuh"
52 #include "gromacs/mdlib/nbnxn_consts.h"
53 #include "gromacs/mdlib/nbnxn_pairlist.h"
54 #include "gromacs/mdtypes/interaction_const.h"
55 #include "gromacs/timing/gpu_timing.h"
58 /*! \brief Macro definining default for the prune kernel's j4 processing concurrency.
60 * The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
62 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
63 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY 4
65 /*! \brief Default for the prune kernel's j4 processing concurrency.
67 * Initialized using the #GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro which allows compile-time override.
69 const int c_cudaPruneKernelJ4Concurrency = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY;
71 /* TODO: consider moving this to kernel_utils */
72 /* Convenience defines */
73 /*! \brief number of clusters per supercluster. */
74 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
75 /*! \brief cluster size = number of atoms per cluster. */
76 static const int c_clSize = c_nbnxnGpuClusterSize;
82 /*! \brief Electrostatic CUDA kernel flavors.
84 * Types of electrostatics implementations available in the CUDA non-bonded
85 * force kernels. These represent both the electrostatics types implemented
86 * by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
87 * enums.h) as well as encode implementation details analytical/tabulated
88 * and single or twin cut-off (for Ewald kernels).
89 * Note that the cut-off and RF kernels have only analytical flavor and unlike
90 * in the CPU kernels, the tabulated kernels are ATM Ewald-only.
92 * The row-order of pointers to different electrostatic kernels defined in
93 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
94 * should match the order of enumerated types below.
97 eelCuCUT, eelCuRF, eelCuEWALD_TAB, eelCuEWALD_TAB_TWIN, eelCuEWALD_ANA, eelCuEWALD_ANA_TWIN, eelCuNR
100 /*! \brief VdW CUDA kernel flavors.
102 * The enumerates values correspond to the LJ implementations in the CUDA non-bonded
105 * The column-order of pointers to different electrostatic kernels defined in
106 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
107 * should match the order of enumerated types below.
110 evdwCuCUT, evdwCuCUTCOMBGEOM, evdwCuCUTCOMBLB, evdwCuFSWITCH, evdwCuPSWITCH, evdwCuEWALDGEOM, evdwCuEWALDLB, evdwCuNR
113 /* All structs prefixed with "cu_" hold data used in GPU calculations and
114 * are passed to the kernels, except cu_timers_t. */
116 typedef struct cu_plist cu_plist_t;
117 typedef struct cu_atomdata cu_atomdata_t;
118 typedef struct cu_nbparam cu_nbparam_t;
119 typedef struct cu_timers cu_timers_t;
120 typedef struct nb_staging nb_staging_t;
125 * \brief Staging area for temporary data downloaded from the GPU.
127 * The energies/shift forces get downloaded here first, before getting added
128 * to the CPU-side aggregate values.
132 float *e_lj; /**< LJ energy */
133 float *e_el; /**< electrostatic energy */
134 float3 *fshift; /**< shift forces */
138 * \brief Nonbonded atom data - both inputs and outputs.
142 int natoms; /**< number of atoms */
143 int natoms_local; /**< number of local atoms */
144 int nalloc; /**< allocation size for the atom data (xq, f) */
146 float4 *xq; /**< atom coordinates + charges, size natoms */
147 float3 *f; /**< force output array, size natoms */
149 float *e_lj; /**< LJ energy output, size 1 */
150 float *e_el; /**< Electrostatics energy input, size 1 */
152 float3 *fshift; /**< shift forces */
154 int ntypes; /**< number of atom types */
155 int *atom_types; /**< atom type indices, size natoms */
156 float2 *lj_comb; /**< sqrt(c6),sqrt(c12) size natoms */
158 float3 *shift_vec; /**< shifts */
159 bool bShiftVecUploaded; /**< true if the shift vector has been uploaded */
163 * \brief Parameters required for the CUDA nonbonded calculations.
168 int eeltype; /**< type of electrostatics, takes values from #eelCu */
169 int vdwtype; /**< type of VdW impl., takes values from #evdwCu */
171 float epsfac; /**< charge multiplication factor */
172 float c_rf; /**< Reaction-field/plain cutoff electrostatics const. */
173 float two_k_rf; /**< Reaction-field electrostatics constant */
174 float ewald_beta; /**< Ewald/PME parameter */
175 float sh_ewald; /**< Ewald/PME correction term substracted from the direct-space potential */
176 float sh_lj_ewald; /**< LJ-Ewald/PME correction term added to the correction potential */
177 float ewaldcoeff_lj; /**< LJ-Ewald/PME coefficient */
179 float rcoulomb_sq; /**< Coulomb cut-off squared */
181 float rvdw_sq; /**< VdW cut-off squared */
182 float rvdw_switch; /**< VdW switched cut-off */
183 float rlistOuter_sq; /**< Full, outer pair-list cut-off squared */
184 float rlistInner_sq; /**< Inner, dynamic pruned pair-list cut-off squared */
185 bool useDynamicPruning; /**< True if we use dynamic pair-list pruning */
187 shift_consts_t dispersion_shift; /**< VdW shift dispersion constants */
188 shift_consts_t repulsion_shift; /**< VdW shift repulsion constants */
189 switch_consts_t vdw_switch; /**< VdW switch constants */
191 /* LJ non-bonded parameters - accessed through texture memory */
192 float *nbfp; /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
193 cudaTextureObject_t nbfp_texobj; /**< texture object bound to nbfp */
194 float *nbfp_comb; /**< nonbonded parameter table per atom type, 2*ntype elements */
195 cudaTextureObject_t nbfp_comb_texobj; /**< texture object bound to nbfp_texobj */
197 /* Ewald Coulomb force table data - accessed through texture memory */
198 float coulomb_tab_scale; /**< table scale/spacing */
199 float *coulomb_tab; /**< pointer to the table in the device memory */
200 cudaTextureObject_t coulomb_tab_texobj; /**< texture object bound to coulomb_tab */
204 * \brief Pair list data.
208 int na_c; /**< number of atoms per cluster */
210 int nsci; /**< size of sci, # of i clusters in the list */
211 int sci_nalloc; /**< allocation size of sci */
212 nbnxn_sci_t *sci; /**< list of i-cluster ("super-clusters") */
214 int ncj4; /**< total # of 4*j clusters */
215 int cj4_nalloc; /**< allocation size of cj4 */
216 nbnxn_cj4_t *cj4; /**< 4*j cluster list, contains j cluster number
217 and index into the i cluster list */
218 int nimask; /**< # of 4*j clusters * # of warps */
219 int imask_nalloc; /**< allocation size of imask */
220 unsigned int *imask; /**< imask for 2 warps for each 4*j cluster group */
221 nbnxn_excl_t *excl; /**< atom interaction bits */
222 int nexcl; /**< count for excl */
223 int excl_nalloc; /**< allocation size of excl */
225 /* parameter+variables for normal and rolling pruning */
226 bool haveFreshList; /**< true after search, indictes that initial pruning with outer prunning is needed */
227 int rollingPruningNumParts; /**< the number of parts/steps over which one cyle of roling pruning takes places */
228 int rollingPruningPart; /**< the next part to which the roling pruning needs to be applied */
232 * \brief CUDA timers used for timing GPU kernels and H2D/D2H transfers.
234 * The two-sized arrays hold the local and non-local values and should always
235 * be indexed with eintLocal/eintNonlocal.
239 GpuRegionTimer atdat; /**< timer for atom data transfer (every PS step) */
240 GpuRegionTimer nb_h2d[2]; /**< timer for x/q H2D transfers (l/nl, every step) */
241 GpuRegionTimer nb_d2h[2]; /**< timer for f D2H transfer (l/nl, every step) */
242 GpuRegionTimer pl_h2d[2]; /**< timer for pair-list H2D transfers (l/nl, every PS step) */
243 bool didPairlistH2D[2]; /**< true when a pair-list transfer has been done at this step */
244 GpuRegionTimer nb_k[2]; /**< timer for non-bonded kernels (l/nl, every step) */
245 GpuRegionTimer prune_k[2]; /**< timer for the 1st pass list pruning kernel (l/nl, every PS step) */
246 bool didPrune[2]; /**< true when we timed pruning and the timings need to be accounted for */
247 GpuRegionTimer rollingPrune_k[2]; /**< timer for rolling pruning kernels (l/nl, frequency depends on chunk size) */
248 bool didRollingPrune[2]; /**< true when we timed rolling pruning (at the previous step) and the timings need to be accounted for */
252 * \brief Main data structure for CUDA nonbonded force calculations.
254 struct gmx_nbnxn_cuda_t
256 const gmx_device_info_t *dev_info; /**< CUDA device information */
257 bool bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU */
258 cu_atomdata_t *atdat; /**< atom data */
259 cu_nbparam_t *nbparam; /**< parameters required for the non-bonded calc. */
260 cu_plist_t *plist[2]; /**< pair-list data structures (local and non-local) */
261 nb_staging_t nbst; /**< staging area where fshift/energies get downloaded */
263 cudaStream_t stream[2]; /**< local and non-local GPU streams */
265 /** events used for synchronization */
266 cudaEvent_t nonlocal_done; /**< event triggered when the non-local non-bonded kernel
267 is done (and the local transfer can proceed) */
268 cudaEvent_t misc_ops_and_local_H2D_done; /**< event triggered when the tasks issued in
269 the local stream that need to precede the
270 non-local force calculations are done
271 (e.g. f buffer 0-ing, local x/q H2D) */
273 /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
274 * concurrent streams, so we won't time if both l/nl work is done on GPUs.
275 * Timer init/uninit is still done even with timing off so only the condition
276 * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
277 bool bDoTime; /**< True if event-based timing is enabled. */
278 cu_timers_t *timers; /**< CUDA event-based timers. */
279 gmx_wallclock_gpu_t *timings; /**< Timing data. */
286 #endif /* NBNXN_CUDA_TYPES_H */