2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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.
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.
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.
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.
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.
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.
38 * Data types used internally in the nbnxn_ocl module.
40 * \author Anca Hamuraru <anca@streamcomputing.eu>
41 * \ingroup module_mdlib
44 #ifndef NBNXN_OPENCL_TYPES_H
45 #define NBNXN_OPENCL_TYPES_H
47 /*! \brief Declare to OpenCL SDKs that we intend to use OpenCL API
48 features that were deprecated in 2.0, so that they don't warn about
50 #define CL_USE_DEPRECATED_OPENCL_2_0_APIS
52 # include <OpenCL/opencl.h>
54 # include <CL/opencl.h>
57 #include "gromacs/mdlib/nbnxn_pairlist.h"
58 #include "gromacs/mdtypes/interaction_const.h"
59 #include "gromacs/utility/real.h"
61 /* kernel does #include "gromacs/math/utilities.h" */
62 /* Move the actual useful stuff here: */
65 #define M_FLOAT_1_SQRTPI 0.564189583547756f
72 /*! \brief Electrostatic OpenCL kernel flavors.
74 * Types of electrostatics implementations available in the OpenCL non-bonded
75 * force kernels. These represent both the electrostatics types implemented
76 * by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
77 * enums.h) as well as encode implementation details analytical/tabulated
78 * and single or twin cut-off (for Ewald kernels).
79 * Note that the cut-off and RF kernels have only analytical flavor and unlike
80 * in the CPU kernels, the tabulated kernels are ATM Ewald-only.
82 * The row-order of pointers to different electrostatic kernels defined in
83 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
84 * should match the order of enumerated types below.
87 eelOclCUT, eelOclRF, eelOclEWALD_TAB, eelOclEWALD_TAB_TWIN, eelOclEWALD_ANA, eelOclEWALD_ANA_TWIN, eelOclNR
90 /*! \brief VdW OpenCL kernel flavors.
92 * The enumerates values correspond to the LJ implementations in the OpenCL non-bonded
95 * The column-order of pointers to different electrostatic kernels defined in
96 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
97 * should match the order of enumerated types below.
100 evdwOclCUT, evdwOclFSWITCH, evdwOclPSWITCH, evdwOclEWALDGEOM, evdwOclEWALDLB, evdwOclNR
104 * \brief Staging area for temporary data downloaded from the GPU.
106 * The energies/shift forces get downloaded here first, before getting added
107 * to the CPU-side aggregate values.
109 typedef struct cl_nb_staging
111 float *e_lj; /**< LJ energy */
112 float *e_el; /**< electrostatic energy */
113 float (*fshift)[3]; /**< float3 buffer with shift forces */
117 * \brief Nonbonded atom data - both inputs and outputs.
119 typedef struct cl_atomdata
121 int natoms; /**< number of atoms */
122 int natoms_local; /**< number of local atoms */
123 int nalloc; /**< allocation size for the atom data (xq, f) */
125 cl_mem xq; /**< float4 buffer with atom coordinates + charges, size natoms */
127 cl_mem f; /**< float3 buffer with force output array, size natoms */
128 size_t f_elem_size; /**< Size in bytes for one element of f buffer */
130 cl_mem e_lj; /**< LJ energy output, size 1 */
131 cl_mem e_el; /**< Electrostatics energy input, size 1 */
133 cl_mem fshift; /**< float3 buffer with shift forces */
134 size_t fshift_elem_size; /**< Size in bytes for one element of fshift buffer */
136 int ntypes; /**< number of atom types */
137 cl_mem atom_types; /**< int buffer with atom type indices, size natoms */
139 cl_mem shift_vec; /**< float3 buffer with shifts values */
140 size_t shift_vec_elem_size; /**< Size in bytes for one element of shift_vec buffer */
142 cl_bool bShiftVecUploaded; /**< true if the shift vector has been uploaded */
146 * \brief Parameters required for the OpenCL nonbonded calculations.
148 typedef struct cl_nbparam
151 int eeltype; /**< type of electrostatics, takes values from #eelOcl */
152 int vdwtype; /**< type of VdW impl., takes values from #evdwOcl */
154 float epsfac; /**< charge multiplication factor */
155 float c_rf; /**< Reaction-field/plain cutoff electrostatics const. */
156 float two_k_rf; /**< Reaction-field electrostatics constant */
157 float ewald_beta; /**< Ewald/PME parameter */
158 float sh_ewald; /**< Ewald/PME correction term substracted from the direct-space potential */
159 float sh_lj_ewald; /**< LJ-Ewald/PME correction term added to the correction potential */
160 float ewaldcoeff_lj; /**< LJ-Ewald/PME coefficient */
162 float rcoulomb_sq; /**< Coulomb cut-off squared */
164 float rvdw_sq; /**< VdW cut-off squared */
165 float rvdw_switch; /**< VdW switched cut-off */
166 float rlist_sq; /**< pair-list cut-off squared */
168 shift_consts_t dispersion_shift; /**< VdW shift dispersion constants */
169 shift_consts_t repulsion_shift; /**< VdW shift repulsion constants */
170 switch_consts_t vdw_switch; /**< VdW switch constants */
172 /* LJ non-bonded parameters - accessed through texture memory */
173 cl_mem nbfp_climg2d; /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
174 cl_mem nbfp_comb_climg2d; /**< nonbonded parameter table per atom type, 2*ntype elements */
176 /* Ewald Coulomb force table data - accessed through texture memory */
177 int coulomb_tab_size; /**< table size (s.t. it fits in texture cache) */
178 float coulomb_tab_scale; /**< table scale/spacing */
179 cl_mem coulomb_tab_climg2d; /**< pointer to the table in the device memory */
183 * \brief Data structure shared between the OpenCL device code and OpenCL host code
185 * Must not contain OpenCL objects (buffers)
186 * TODO: review, improve */
187 typedef struct cl_nbparam_params
190 int eeltype; /**< type of electrostatics, takes values from #eelCu */
191 int vdwtype; /**< type of VdW impl., takes values from #evdwCu */
193 float epsfac; /**< charge multiplication factor */
194 float c_rf; /**< Reaction-field/plain cutoff electrostatics const. */
195 float two_k_rf; /**< Reaction-field electrostatics constant */
196 float ewald_beta; /**< Ewald/PME parameter */
197 float sh_ewald; /**< Ewald/PME correction term substracted from the direct-space potential */
198 float sh_lj_ewald; /**< LJ-Ewald/PME correction term added to the correction potential */
199 float ewaldcoeff_lj; /**< LJ-Ewald/PME coefficient */
201 float rcoulomb_sq; /**< Coulomb cut-off squared */
203 float rvdw_sq; /**< VdW cut-off squared */
204 float rvdw_switch; /**< VdW switched cut-off */
205 float rlist_sq; /**< pair-list cut-off squared */
207 shift_consts_t dispersion_shift; /**< VdW shift dispersion constants */
208 shift_consts_t repulsion_shift; /**< VdW shift repulsion constants */
209 switch_consts_t vdw_switch; /**< VdW switch constants */
211 /* Ewald Coulomb force table data - accessed through texture memory */
212 int coulomb_tab_size; /**< table size (s.t. it fits in texture cache) */
213 float coulomb_tab_scale; /**< table scale/spacing */
214 } cl_nbparam_params_t;
218 * \brief Pair list data.
220 typedef struct cl_plist
222 int na_c; /**< number of atoms per cluster */
224 int nsci; /**< size of sci, # of i clusters in the list */
225 int sci_nalloc; /**< allocation size of sci */
226 cl_mem sci; /**< list of i-cluster ("super-clusters").
227 It contains elements of type nbnxn_sci_t */
229 int ncj4; /**< total # of 4*j clusters */
230 int cj4_nalloc; /**< allocation size of cj4 */
231 cl_mem cj4; /**< 4*j cluster list, contains j cluster number and
232 index into the i cluster list.
233 It contains elements of type nbnxn_cj4_t */
234 cl_mem excl; /**< atom interaction bits
235 It contains elements of type nbnxn_excl_t */
236 int nexcl; /**< count for excl */
237 int excl_nalloc; /**< allocation size of excl */
239 cl_bool bDoPrune; /**< true if pair-list pruning needs to be
240 done during the current step */
245 * \brief OpenCL events used for timing GPU kernels and H2D/D2H transfers.
247 * The two-sized arrays hold the local and non-local values and should always
248 * be indexed with eintLocal/eintNonlocal.
250 typedef struct cl_timers
252 cl_event atdat; /**< event for atom data transfer (every PS step) */
254 cl_event nb_h2d[2]; /**< events for x/q H2D transfers (l/nl, every step) */
256 cl_event nb_d2h_f[2]; /**< events for f D2H transfer (l/nl, every step) */
257 cl_event nb_d2h_fshift[2]; /**< events for fshift D2H transfer (l/nl, every step) */
258 cl_event nb_d2h_e_el[2]; /**< events for e_el D2H transfer (l/nl, every step) */
259 cl_event nb_d2h_e_lj[2]; /**< events for e_lj D2H transfer (l/nl, every step) */
261 cl_event pl_h2d_sci[2]; /**< events for pair-list sci H2D transfers (l/nl, every PS step) */
262 cl_event pl_h2d_cj4[2]; /**< events for pair-list cj4 H2D transfers (l/nl, every PS step) */
263 cl_event pl_h2d_excl[2]; /**< events for pair-list excl H2D transfers (l/nl, every PS step)*/
265 cl_event nb_k[2]; /**< event for non-bonded kernels (l/nl, every step) */
269 * \brief Main data structure for OpenCL nonbonded force calculations.
271 struct gmx_nbnxn_ocl_t
273 struct gmx_device_info_t *dev_info; /**< OpenCL device information */
275 /**< Pointers to non-bonded kernel functions
276 * organized similar with nb_kfunc_xxx arrays in nbnxn_ocl.cpp */
278 cl_kernel kernel_noener_noprune_ptr[eelOclNR][evdwOclNR];
279 cl_kernel kernel_ener_noprune_ptr[eelOclNR][evdwOclNR];
280 cl_kernel kernel_noener_prune_ptr[eelOclNR][evdwOclNR];
281 cl_kernel kernel_ener_prune_ptr[eelOclNR][evdwOclNR];
284 /**< auxiliary kernels implementing memset-like functions */
286 cl_kernel kernel_memset_f;
287 cl_kernel kernel_memset_f2;
288 cl_kernel kernel_memset_f3;
289 cl_kernel kernel_zero_e_fshift;
292 cl_bool bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU */
293 cl_bool bNonLocalStreamActive; /**< true indicates that the nonlocal_done event was enqueued */
295 cl_atomdata_t *atdat; /**< atom data */
296 cl_nbparam_t *nbparam; /**< parameters required for the non-bonded calc. */
297 cl_plist_t *plist[2]; /**< pair-list data structures (local and non-local) */
298 cl_nb_staging_t nbst; /**< staging area where fshift/energies get downloaded */
300 cl_mem debug_buffer; /**< debug buffer */
302 cl_command_queue stream[2]; /**< local and non-local GPU queues */
304 /** events used for synchronization */
305 cl_event nonlocal_done; /**< event triggered when the non-local non-bonded kernel
306 is done (and the local transfer can proceed) */
307 cl_event misc_ops_and_local_H2D_done; /**< event triggered when the tasks issued in
308 the local stream that need to precede the
309 non-local force calculations are done
310 (e.g. f buffer 0-ing, local x/q H2D) */
312 cl_bool bDoTime; /**< True if event-based timing is enabled. */
313 cl_timers_t *timers; /**< OpenCL event-based timers. */
314 struct gmx_wallclock_gpu_t *timings; /**< Timing data. */
321 #endif /* NBNXN_OPENCL_TYPES_H */