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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
39 #ifndef NBNXN_CUDA_TYPES_H
40 #define NBNXN_CUDA_TYPES_H
42 #include "types/nbnxn_pairlist.h"
43 #include "types/nbnxn_cuda_types_ext.h"
44 #include "../../gmxlib/cuda_tools/cudautils.cuh"
46 /* CUDA versions from 5.0 above support texture objects. */
47 #if CUDA_VERSION >= 5000
48 #define TEXOBJ_SUPPORTED
49 #else /* CUDA_VERSION */
50 /* This typedef allows us to define only one version of struct cu_nbparam */
51 typedef int cudaTextureObject_t;
52 #endif /* CUDA_VERSION */
58 /** Types of electrostatics implementations available in the CUDA non-bonded
59 * force kernels. These represent both the electrostatics types implemented
60 * by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
61 * enums.h) as well as encode implementation details analytical/tabulated
62 * and single or twin cut-off (for Ewald kernels).
63 * Note that the cut-off and RF kernels have only analytical flavor and unlike
64 * in the CPU kernels, the tabulated kernels are ATM Ewald-only.
66 * The order of pointers to different electrostatic kernels defined in
67 * nbnxn_cuda.cu by the nb_default_kfunc_ptr and nb_legacy_kfunc_ptr arrays
68 * should match the order of enumerated types below. */
70 eelCuCUT, eelCuRF, eelCuEWALD_TAB, eelCuEWALD_TAB_TWIN, eelCuEWALD_ANA, eelCuEWALD_ANA_TWIN, eelCuNR
73 /** Kernel flavors with different set of optimizations: default for CUDA <=v4.1
74 * compilers and legacy for earlier, 3.2 and 4.0 CUDA compilers. */
76 eNbnxnCuKDefault, eNbnxnCuKLegacy, eNbnxnCuKNR
79 #define NBNXN_KVER_OLD(k) (k == eNbnxnCuKOld)
80 #define NBNXN_KVER_LEGACY(k) (k == eNbnxnCuKLegacy)
81 #define NBNXN_KVER_DEFAULT(k) (k == eNbnxnCuKDefault)
83 /* Non-bonded kernel versions. */
85 /* All structs prefixed with "cu_" hold data used in GPU calculations and
86 * are passed to the kernels, except cu_timers_t. */
87 typedef struct cu_plist cu_plist_t;
88 typedef struct cu_atomdata cu_atomdata_t;
89 typedef struct cu_nbparam cu_nbparam_t;
90 typedef struct cu_timers cu_timers_t;
91 typedef struct nb_staging nb_staging_t;
94 /** Staging area for temporary data. The energies get downloaded here first,
95 * before getting added to the CPU-side aggregate values.
99 float *e_lj; /**< LJ energy */
100 float *e_el; /**< electrostatic energy */
101 float3 *fshift; /**< shift forces */
104 /** Nonbonded atom data -- both inputs and outputs. */
107 int natoms; /**< number of atoms */
108 int natoms_local; /**< number of local atoms */
109 int nalloc; /**< allocation size for the atom data (xq, f) */
111 float4 *xq; /**< atom coordinates + charges, size natoms */
112 float3 *f; /**< force output array, size natoms */
113 /* TODO: try float2 for the energies */
114 float *e_lj, /**< LJ energy output, size 1 */
115 *e_el; /**< Electrostatics energy input, size 1 */
117 float3 *fshift; /**< shift forces */
119 int ntypes; /**< number of atom types */
120 int *atom_types; /**< atom type indices, size natoms */
122 float3 *shift_vec; /**< shifts */
123 bool bShiftVecUploaded; /**< true if the shift vector has been uploaded */
126 /** Parameters required for the CUDA nonbonded calculations. */
129 int eeltype; /**< type of electrostatics */
131 float epsfac; /**< charge multiplication factor */
132 float c_rf; /**< Reaction-field/plain cutoff electrostatics const. */
133 float two_k_rf; /**< Reaction-field electrostatics constant */
134 float ewald_beta; /**< Ewald/PME parameter */
135 float sh_ewald; /**< Ewald/PME correction term */
136 float rvdw_sq; /**< VdW cut-off */
137 float rcoulomb_sq; /**< Coulomb cut-off */
138 float rlist_sq; /**< pair-list cut-off */
139 float sh_invrc6; /**< LJ potential correction term */
141 /* Non-bonded parameters - accessed through texture memory */
142 float *nbfp; /**< nonbonded parameter table with C6/C12 pairs */
143 cudaTextureObject_t nbfp_texobj; /**< texture object bound to nbfp */
145 /* Ewald Coulomb force table data - accessed through texture memory */
146 int coulomb_tab_size; /**< table size (s.t. it fits in texture cache) */
147 float coulomb_tab_scale; /**< table scale/spacing */
148 float *coulomb_tab; /**< pointer to the table in the device memory */
149 cudaTextureObject_t coulomb_tab_texobj; /**< texture object bound to coulomb_tab */
152 /** Pair list data */
155 int na_c; /**< number of atoms per cluster */
157 int nsci; /**< size of sci, # of i clusters in the list */
158 int sci_nalloc; /**< allocation size of sci */
159 nbnxn_sci_t *sci; /**< list of i-cluster ("super-clusters") */
161 int ncj4; /**< total # of 4*j clusters */
162 int cj4_nalloc; /**< allocation size of cj4 */
163 nbnxn_cj4_t *cj4; /**< 4*j cluster list, contains j cluster number
164 and index into the i cluster list */
165 nbnxn_excl_t *excl; /**< atom interaction bits */
166 int nexcl; /**< count for excl */
167 int excl_nalloc; /**< allocation size of excl */
169 bool bDoPrune; /**< true if pair-list pruning needs to be
170 done during the current step */
173 /** CUDA events used for timing GPU kernels and H2D/D2H transfers.
174 * The two-sized arrays hold the local and non-local values and should always
175 * be indexed with eintLocal/eintNonlocal.
179 cudaEvent_t start_atdat; /**< start event for atom data transfer (every PS step) */
180 cudaEvent_t stop_atdat; /**< stop event for atom data transfer (every PS step) */
181 cudaEvent_t start_nb_h2d[2]; /**< start events for x/q H2D transfers (l/nl, every step) */
182 cudaEvent_t stop_nb_h2d[2]; /**< stop events for x/q H2D transfers (l/nl, every step) */
183 cudaEvent_t start_nb_d2h[2]; /**< start events for f D2H transfer (l/nl, every step) */
184 cudaEvent_t stop_nb_d2h[2]; /**< stop events for f D2H transfer (l/nl, every step) */
185 cudaEvent_t start_pl_h2d[2]; /**< start events for pair-list H2D transfers (l/nl, every PS step) */
186 cudaEvent_t stop_pl_h2d[2]; /**< start events for pair-list H2D transfers (l/nl, every PS step) */
187 cudaEvent_t start_nb_k[2]; /**< start event for non-bonded kernels (l/nl, every step) */
188 cudaEvent_t stop_nb_k[2]; /**< stop event non-bonded kernels (l/nl, every step) */
191 /** Main data structure for CUDA nonbonded force calculations. */
194 cuda_dev_info_t *dev_info; /**< CUDA device information */
195 int kernel_ver; /**< The version of the kernel to be executed on the
196 device in use, possible values: eNbnxnCuK* */
197 bool bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU */
198 bool bUseStreamSync; /**< true if the standard cudaStreamSynchronize is used
199 and not memory polling-based waiting */
200 cu_atomdata_t *atdat; /**< atom data */
201 cu_nbparam_t *nbparam; /**< parameters required for the non-bonded calc. */
202 cu_plist_t *plist[2]; /**< pair-list data structures (local and non-local) */
203 nb_staging_t nbst; /**< staging area where fshift/energies get downloaded */
205 cudaStream_t stream[2]; /**< local and non-local GPU streams */
207 /** events used for synchronization */
208 cudaEvent_t nonlocal_done; /**< event triggered when the non-local non-bonded kernel
209 is done (and the local transfer can proceed) */
210 cudaEvent_t misc_ops_done; /**< event triggered when the operations that precede the
211 main force calculations are done (e.g. buffer 0-ing) */
213 /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
214 * concurrent streams, so we won't time if both l/nl work is done on GPUs.
215 * Timer init/uninit is still done even with timing off so only the condition
216 * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
217 bool bDoTime; /**< True if event-based timing is enabled. */
218 cu_timers_t *timers; /**< CUDA event-based timers. */
219 wallclock_gpu_t *timings; /**< Timing data. */
226 #endif /* NBNXN_CUDA_TYPES_H */