2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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 * \author Szilárd Páll <pszilard@kth.se>
42 * \ingroup module_mdlib
45 #ifndef NBNXN_OPENCL_TYPES_H
46 #define NBNXN_OPENCL_TYPES_H
48 #include "gromacs/gpu_utils/gmxopencl.h"
49 #include "gromacs/gpu_utils/gputraits_ocl.h"
50 #include "gromacs/gpu_utils/oclutils.h"
51 #include "gromacs/mdlib/nbnxn_gpu_types_common.h"
52 #include "gromacs/mdlib/nbnxn_pairlist.h"
53 #include "gromacs/mdtypes/interaction_const.h"
54 #include "gromacs/utility/real.h"
56 /* kernel does #include "gromacs/math/utilities.h" */
57 /* Move the actual useful stuff here: */
60 #define M_FLOAT_1_SQRTPI 0.564189583547756f
62 /*! \brief Macros defining platform-dependent defaults for the prune kernel's j4 processing concurrency.
64 * The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
67 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
68 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_AMD 4
69 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_NVIDIA 4
70 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_DEFAULT 4
72 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_AMD GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
73 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_NVIDIA GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
74 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_DEFAULT GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
77 /*! \brief Constants for platform-dependent defaults for the prune kernel's j4 processing concurrency.
79 * Initialized using macros that can be overridden at compile-time (using #GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY).
82 const int c_oclPruneKernelJ4ConcurrencyAMD = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_AMD;
83 const int c_oclPruneKernelJ4ConcurrencyNVIDIA = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_NVIDIA;
84 const int c_oclPruneKernelJ4ConcurrencyDefault = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_DEFAULT;
87 /*! \brief Returns the j4 processing concurrency parameter for the vendor \p vendorId
88 * \param vendorId takes values from #ocl_vendor_id_t.
90 static inline int getOclPruneKernelJ4Concurrency(int vendorId)
92 assert(vendorId < OCL_VENDOR_UNKNOWN);
95 case OCL_VENDOR_AMD: return c_oclPruneKernelJ4ConcurrencyAMD; break;
96 case OCL_VENDOR_NVIDIA: return c_oclPruneKernelJ4ConcurrencyNVIDIA; break;
97 default: return c_oclPruneKernelJ4ConcurrencyDefault; break;
106 /*! \brief Electrostatic OpenCL kernel flavors.
108 * Types of electrostatics implementations available in the OpenCL non-bonded
109 * force kernels. These represent both the electrostatics types implemented
110 * by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
111 * enums.h) as well as encode implementation details analytical/tabulated
112 * and single or twin cut-off (for Ewald kernels).
113 * Note that the cut-off and RF kernels have only analytical flavor and unlike
114 * in the CPU kernels, the tabulated kernels are ATM Ewald-only.
116 * The row-order of pointers to different electrostatic kernels defined in
117 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
118 * should match the order of enumerated types below.
121 eelOclCUT, eelOclRF, eelOclEWALD_TAB, eelOclEWALD_TAB_TWIN, eelOclEWALD_ANA, eelOclEWALD_ANA_TWIN, eelOclNR
124 /*! \brief VdW OpenCL kernel flavors.
126 * The enumerates values correspond to the LJ implementations in the OpenCL non-bonded
129 * The column-order of pointers to different electrostatic kernels defined in
130 * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
131 * should match the order of enumerated types below.
134 evdwOclCUT, evdwOclCUTCOMBGEOM, evdwOclCUTCOMBLB, evdwOclFSWITCH, evdwOclPSWITCH, evdwOclEWALDGEOM, evdwOclEWALDLB, evdwOclNR
137 /*! \brief Pruning kernel flavors.
139 * The values correspond to the first call of the pruning post-list generation
140 * and the rolling pruning, respectively.
143 epruneFirst, epruneRolling, ePruneNR
147 * \brief Staging area for temporary data downloaded from the GPU.
149 * The energies/shift forces get downloaded here first, before getting added
150 * to the CPU-side aggregate values.
152 typedef struct cl_nb_staging
154 float *e_lj; /**< LJ energy */
155 float *e_el; /**< electrostatic energy */
156 float (*fshift)[3]; /**< float3 buffer with shift forces */
160 * \brief Nonbonded atom data - both inputs and outputs.
162 typedef struct cl_atomdata
164 int natoms; /**< number of atoms */
165 int natoms_local; /**< number of local atoms */
166 int nalloc; /**< allocation size for the atom data (xq, f) */
168 cl_mem xq; /**< float4 buffer with atom coordinates + charges, size natoms */
170 cl_mem f; /**< float3 buffer with force output array, size natoms */
171 size_t f_elem_size; /**< Size in bytes for one element of f buffer */
173 cl_mem e_lj; /**< LJ energy output, size 1 */
174 cl_mem e_el; /**< Electrostatics energy input, size 1 */
176 cl_mem fshift; /**< float3 buffer with shift forces */
177 size_t fshift_elem_size; /**< Size in bytes for one element of fshift buffer */
179 int ntypes; /**< number of atom types */
180 cl_mem atom_types; /**< int buffer with atom type indices, size natoms */
181 cl_mem lj_comb; /**< float2 buffer with sqrt(c6),sqrt(c12), size natoms */
183 cl_mem shift_vec; /**< float3 buffer with shifts values */
184 size_t shift_vec_elem_size; /**< Size in bytes for one element of shift_vec buffer */
186 cl_bool bShiftVecUploaded; /**< true if the shift vector has been uploaded */
190 * \brief Parameters required for the OpenCL nonbonded calculations.
192 typedef struct cl_nbparam
195 int eeltype; /**< type of electrostatics, takes values from #eelOcl */
196 int vdwtype; /**< type of VdW impl., takes values from #evdwOcl */
198 float epsfac; /**< charge multiplication factor */
199 float c_rf; /**< Reaction-field/plain cutoff electrostatics const. */
200 float two_k_rf; /**< Reaction-field electrostatics constant */
201 float ewald_beta; /**< Ewald/PME parameter */
202 float sh_ewald; /**< Ewald/PME correction term substracted from the direct-space potential */
203 float sh_lj_ewald; /**< LJ-Ewald/PME correction term added to the correction potential */
204 float ewaldcoeff_lj; /**< LJ-Ewald/PME coefficient */
206 float rcoulomb_sq; /**< Coulomb cut-off squared */
208 float rvdw_sq; /**< VdW cut-off squared */
209 float rvdw_switch; /**< VdW switched cut-off */
210 float rlistOuter_sq; /**< Full, outer pair-list cut-off squared */
211 float rlistInner_sq; /**< Inner, dynamic pruned pair-list cut-off squared */
212 bool useDynamicPruning; /**< True if we use dynamic pair-list pruning */
214 shift_consts_t dispersion_shift; /**< VdW shift dispersion constants */
215 shift_consts_t repulsion_shift; /**< VdW shift repulsion constants */
216 switch_consts_t vdw_switch; /**< VdW switch constants */
218 /* LJ non-bonded parameters - accessed through texture memory */
219 cl_mem nbfp_climg2d; /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
220 cl_mem nbfp_comb_climg2d; /**< nonbonded parameter table per atom type, 2*ntype elements */
222 /* Ewald Coulomb force table data - accessed through texture memory */
223 float coulomb_tab_scale; /**< table scale/spacing */
224 cl_mem coulomb_tab_climg2d; /**< pointer to the table in the device memory */
228 * \brief Data structure shared between the OpenCL device code and OpenCL host code
230 * Must not contain OpenCL objects (buffers)
231 * TODO: review, improve */
232 typedef struct cl_nbparam_params
235 int eeltype; /**< type of electrostatics, takes values from #eelCu */
236 int vdwtype; /**< type of VdW impl., takes values from #evdwCu */
238 float epsfac; /**< charge multiplication factor */
239 float c_rf; /**< Reaction-field/plain cutoff electrostatics const. */
240 float two_k_rf; /**< Reaction-field electrostatics constant */
241 float ewald_beta; /**< Ewald/PME parameter */
242 float sh_ewald; /**< Ewald/PME correction term substracted from the direct-space potential */
243 float sh_lj_ewald; /**< LJ-Ewald/PME correction term added to the correction potential */
244 float ewaldcoeff_lj; /**< LJ-Ewald/PME coefficient */
246 float rcoulomb_sq; /**< Coulomb cut-off squared */
248 float rvdw_sq; /**< VdW cut-off squared */
249 float rvdw_switch; /**< VdW switched cut-off */
250 float rlistOuter_sq; /**< Full, outer pair-list cut-off squared */
251 float rlistInner_sq; /**< Inner, dynamic pruned pair-list cut-off squared */
253 shift_consts_t dispersion_shift; /**< VdW shift dispersion constants */
254 shift_consts_t repulsion_shift; /**< VdW shift repulsion constants */
255 switch_consts_t vdw_switch; /**< VdW switch constants */
257 /* Ewald Coulomb force table data - accessed through texture memory */
258 float coulomb_tab_scale; /**< table scale/spacing */
259 } cl_nbparam_params_t;
263 * \brief Pair list data.
265 using cl_plist_t = gpu_plist;
268 * \brief Typedef of actual timer type.
270 typedef struct nbnxn_gpu_timers_t cl_timers_t;
273 * \brief Main data structure for OpenCL nonbonded force calculations.
275 struct gmx_nbnxn_ocl_t
277 const gmx_device_info_t *dev_info; /**< OpenCL device information */
278 struct gmx_device_runtime_data_t *dev_rundata; /**< OpenCL runtime data (context, kernels) */
280 /**< Pointers to non-bonded kernel functions
281 * organized similar with nb_kfunc_xxx arrays in nbnxn_ocl.cpp */
283 cl_kernel kernel_noener_noprune_ptr[eelOclNR][evdwOclNR];
284 cl_kernel kernel_ener_noprune_ptr[eelOclNR][evdwOclNR];
285 cl_kernel kernel_noener_prune_ptr[eelOclNR][evdwOclNR];
286 cl_kernel kernel_ener_prune_ptr[eelOclNR][evdwOclNR];
288 cl_kernel kernel_pruneonly[ePruneNR]; /**< prune kernels, ePruneKind defined the kernel kinds */
290 bool bPrefetchLjParam; /**< true if prefetching fg i-atom LJ parameters should be used in the kernels */
292 /**< auxiliary kernels implementing memset-like functions */
294 cl_kernel kernel_memset_f;
295 cl_kernel kernel_memset_f2;
296 cl_kernel kernel_memset_f3;
297 cl_kernel kernel_zero_e_fshift;
300 cl_bool bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU */
301 cl_bool bNonLocalStreamActive; /**< true indicates that the nonlocal_done event was enqueued */
303 cl_atomdata_t *atdat; /**< atom data */
304 cl_nbparam_t *nbparam; /**< parameters required for the non-bonded calc. */
305 cl_plist_t *plist[2]; /**< pair-list data structures (local and non-local) */
306 cl_nb_staging_t nbst; /**< staging area where fshift/energies get downloaded */
308 cl_command_queue stream[2]; /**< local and non-local GPU queues */
310 /** events used for synchronization */
311 cl_event nonlocal_done; /**< event triggered when the non-local non-bonded kernel
312 is done (and the local transfer can proceed) */
313 cl_event misc_ops_and_local_H2D_done; /**< event triggered when the tasks issued in
314 the local stream that need to precede the
315 non-local force calculations are done
316 (e.g. f buffer 0-ing, local x/q H2D) */
318 cl_bool bDoTime; /**< True if event-based timing is enabled. */
319 cl_timers_t *timers; /**< OpenCL event-based timers. */
320 struct gmx_wallclock_gpu_nbnxn_t *timings; /**< Timing data. TODO: deprecate this and query timers for accumulated data instead */
327 #endif /* NBNXN_OPENCL_TYPES_H */