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/mdlib/nbnxn_consts.h"
52 #include "gromacs/mdlib/nbnxn_pairlist.h"
53 #include "gromacs/mdtypes/interaction_const.h"
54 #include "gromacs/timing/gpu_timing.h"
57 /*! \brief Macro definining default for the prune kernel's j4 processing concurrency.
59 * The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
61 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
62 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY 4
64 /*! \brief Default for the prune kernel's j4 processing concurrency.
66 * Initialized using the #GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro which allows compile-time override.
68 const int c_cudaPruneKernelJ4Concurrency = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY;
70 /* TODO: consider moving this to kernel_utils */
71 /* Convenience defines */
72 /*! \brief number of clusters per supercluster. */
73 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
74 /*! \brief cluster size = number of atoms per cluster. */
75 static const int c_clSize = c_nbnxnGpuClusterSize;
81 /*! \brief Electrostatic CUDA kernel flavors.
83 * Types of electrostatics implementations available in the CUDA non-bonded
84 * force kernels. These represent both the electrostatics types implemented
85 * by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
86 * enums.h) as well as encode implementation details analytical/tabulated
87 * and single or twin cut-off (for Ewald kernels).
88 * Note that the cut-off and RF kernels have only analytical flavor and unlike
89 * in the CPU kernels, the tabulated kernels are ATM Ewald-only.
91 * The row-order of pointers to different electrostatic kernels defined in
92 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
93 * should match the order of enumerated types below.
96 eelCuCUT, eelCuRF, eelCuEWALD_TAB, eelCuEWALD_TAB_TWIN, eelCuEWALD_ANA, eelCuEWALD_ANA_TWIN, eelCuNR
99 /*! \brief VdW CUDA kernel flavors.
101 * The enumerates values correspond to the LJ implementations in the CUDA non-bonded
104 * The column-order of pointers to different electrostatic kernels defined in
105 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
106 * should match the order of enumerated types below.
109 evdwCuCUT, evdwCuCUTCOMBGEOM, evdwCuCUTCOMBLB, evdwCuFSWITCH, evdwCuPSWITCH, evdwCuEWALDGEOM, evdwCuEWALDLB, evdwCuNR
112 /* All structs prefixed with "cu_" hold data used in GPU calculations and
113 * are passed to the kernels, except cu_timers_t. */
115 typedef struct cu_plist cu_plist_t;
116 typedef struct cu_atomdata cu_atomdata_t;
117 typedef struct cu_nbparam cu_nbparam_t;
118 typedef struct cu_timers cu_timers_t;
119 typedef struct nb_staging nb_staging_t;
124 * \brief Staging area for temporary data downloaded from the GPU.
126 * The energies/shift forces get downloaded here first, before getting added
127 * to the CPU-side aggregate values.
131 float *e_lj; /**< LJ energy */
132 float *e_el; /**< electrostatic energy */
133 float3 *fshift; /**< shift forces */
137 * \brief Nonbonded atom data - both inputs and outputs.
141 int natoms; /**< number of atoms */
142 int natoms_local; /**< number of local atoms */
143 int nalloc; /**< allocation size for the atom data (xq, f) */
145 float4 *xq; /**< atom coordinates + charges, size natoms */
146 float3 *f; /**< force output array, size natoms */
148 float *e_lj; /**< LJ energy output, size 1 */
149 float *e_el; /**< Electrostatics energy input, size 1 */
151 float3 *fshift; /**< shift forces */
153 int ntypes; /**< number of atom types */
154 int *atom_types; /**< atom type indices, size natoms */
155 float2 *lj_comb; /**< sqrt(c6),sqrt(c12) size natoms */
157 float3 *shift_vec; /**< shifts */
158 bool bShiftVecUploaded; /**< true if the shift vector has been uploaded */
162 * \brief Parameters required for the CUDA nonbonded calculations.
167 int eeltype; /**< type of electrostatics, takes values from #eelCu */
168 int vdwtype; /**< type of VdW impl., takes values from #evdwCu */
170 float epsfac; /**< charge multiplication factor */
171 float c_rf; /**< Reaction-field/plain cutoff electrostatics const. */
172 float two_k_rf; /**< Reaction-field electrostatics constant */
173 float ewald_beta; /**< Ewald/PME parameter */
174 float sh_ewald; /**< Ewald/PME correction term substracted from the direct-space potential */
175 float sh_lj_ewald; /**< LJ-Ewald/PME correction term added to the correction potential */
176 float ewaldcoeff_lj; /**< LJ-Ewald/PME coefficient */
178 float rcoulomb_sq; /**< Coulomb cut-off squared */
180 float rvdw_sq; /**< VdW cut-off squared */
181 float rvdw_switch; /**< VdW switched cut-off */
182 float rlistOuter_sq; /**< Full, outer pair-list cut-off squared */
183 float rlistInner_sq; /**< Inner, dynamic pruned pair-list cut-off squared */
184 bool useDynamicPruning; /**< True if we use dynamic pair-list pruning */
186 shift_consts_t dispersion_shift; /**< VdW shift dispersion constants */
187 shift_consts_t repulsion_shift; /**< VdW shift repulsion constants */
188 switch_consts_t vdw_switch; /**< VdW switch constants */
190 /* LJ non-bonded parameters - accessed through texture memory */
191 float *nbfp; /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
192 cudaTextureObject_t nbfp_texobj; /**< texture object bound to nbfp */
193 float *nbfp_comb; /**< nonbonded parameter table per atom type, 2*ntype elements */
194 cudaTextureObject_t nbfp_comb_texobj; /**< texture object bound to nbfp_texobj */
196 /* Ewald Coulomb force table data - accessed through texture memory */
197 float coulomb_tab_scale; /**< table scale/spacing */
198 float *coulomb_tab; /**< pointer to the table in the device memory */
199 cudaTextureObject_t coulomb_tab_texobj; /**< texture object bound to coulomb_tab */
203 * \brief Pair list data.
207 int na_c; /**< number of atoms per cluster */
209 int nsci; /**< size of sci, # of i clusters in the list */
210 int sci_nalloc; /**< allocation size of sci */
211 nbnxn_sci_t *sci; /**< list of i-cluster ("super-clusters") */
213 int ncj4; /**< total # of 4*j clusters */
214 int cj4_nalloc; /**< allocation size of cj4 */
215 nbnxn_cj4_t *cj4; /**< 4*j cluster list, contains j cluster number
216 and index into the i cluster list */
217 int nimask; /**< # of 4*j clusters * # of warps */
218 int imask_nalloc; /**< allocation size of imask */
219 unsigned int *imask; /**< imask for 2 warps for each 4*j cluster group */
220 nbnxn_excl_t *excl; /**< atom interaction bits */
221 int nexcl; /**< count for excl */
222 int excl_nalloc; /**< allocation size of excl */
224 /* parameter+variables for normal and rolling pruning */
225 bool haveFreshList; /**< true after search, indictes that initial pruning with outer prunning is needed */
226 int rollingPruningNumParts; /**< the number of parts/steps over which one cyle of roling pruning takes places */
227 int rollingPruningPart; /**< the next part to which the roling pruning needs to be applied */
231 * \brief CUDA events used for timing GPU kernels and H2D/D2H transfers.
233 * The two-sized arrays hold the local and non-local values and should always
234 * be indexed with eintLocal/eintNonlocal.
238 cudaEvent_t start_atdat; /**< start event for atom data transfer (every PS step) */
239 cudaEvent_t stop_atdat; /**< stop event for atom data transfer (every PS step) */
240 cudaEvent_t start_nb_h2d[2]; /**< start events for x/q H2D transfers (l/nl, every step) */
241 cudaEvent_t stop_nb_h2d[2]; /**< stop events for x/q H2D transfers (l/nl, every step) */
242 cudaEvent_t start_nb_d2h[2]; /**< start events for f D2H transfer (l/nl, every step) */
243 cudaEvent_t stop_nb_d2h[2]; /**< stop events for f D2H transfer (l/nl, every step) */
244 cudaEvent_t start_pl_h2d[2]; /**< start events for pair-list H2D transfers (l/nl, every PS step) */
245 cudaEvent_t stop_pl_h2d[2]; /**< start events for pair-list H2D transfers (l/nl, every PS step) */
246 bool didPairlistH2D[2]; /**< true when a pair-list transfer has been done at this step */
247 cudaEvent_t start_nb_k[2]; /**< start event for non-bonded kernels (l/nl, every step) */
248 cudaEvent_t stop_nb_k[2]; /**< stop event non-bonded kernels (l/nl, every step) */
249 cudaEvent_t start_prune_k[2]; /**< start event for the 1st pass list pruning kernel (l/nl, every PS step) */
250 cudaEvent_t stop_prune_k[2]; /**< stop event for the 1st pass list pruning kernel (l/nl, every PS step) */
251 bool didPrune[2]; /**< true when we timed pruning and the timings need to be accounted for */
252 cudaEvent_t start_rollingPrune_k[2]; /**< start event for rolling pruning kernels (l/nl, frequency depends on chunk size) */
253 cudaEvent_t stop_rollingPrune_k[2]; /**< stop event for rolling pruning kernels (l/nl, frequency depends on chunk size) */
254 bool didRollingPrune[2]; /**< true when we timed rolling pruning (at the previous step) and the timings need to be accounted for */
258 * \brief Main data structure for CUDA nonbonded force calculations.
260 struct gmx_nbnxn_cuda_t
262 const gmx_device_info_t *dev_info; /**< CUDA device information */
263 bool bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU */
264 cu_atomdata_t *atdat; /**< atom data */
265 cu_nbparam_t *nbparam; /**< parameters required for the non-bonded calc. */
266 cu_plist_t *plist[2]; /**< pair-list data structures (local and non-local) */
267 nb_staging_t nbst; /**< staging area where fshift/energies get downloaded */
269 cudaStream_t stream[2]; /**< local and non-local GPU streams */
271 /** events used for synchronization */
272 cudaEvent_t nonlocal_done; /**< event triggered when the non-local non-bonded kernel
273 is done (and the local transfer can proceed) */
274 cudaEvent_t misc_ops_and_local_H2D_done; /**< event triggered when the tasks issued in
275 the local stream that need to precede the
276 non-local force calculations are done
277 (e.g. f buffer 0-ing, local x/q H2D) */
279 /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
280 * concurrent streams, so we won't time if both l/nl work is done on GPUs.
281 * Timer init/uninit is still done even with timing off so only the condition
282 * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
283 bool bDoTime; /**< True if event-based timing is enabled. */
284 cu_timers_t *timers; /**< CUDA event-based timers. */
285 gmx_wallclock_gpu_t *timings; /**< Timing data. */
292 #endif /* NBNXN_CUDA_TYPES_H */