4d428afa84df431edaea965de04a53ce3cc8354d
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_ocl / nbnxn_ocl_types.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016,2017, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *  \brief
38  *  Data types used internally in the nbnxn_ocl module.
39  *
40  *  \author Anca Hamuraru <anca@streamcomputing.eu>
41  *  \author Szilárd Páll <pszilard@kth.se>
42  *  \ingroup module_mdlib
43  */
44
45 #ifndef NBNXN_OPENCL_TYPES_H
46 #define NBNXN_OPENCL_TYPES_H
47
48 #include "gromacs/gpu_utils/gmxopencl.h"
49 #include "gromacs/gpu_utils/oclutils.h"
50 #include "gromacs/mdlib/nbnxn_gpu_types_common.h"
51 #include "gromacs/mdlib/nbnxn_pairlist.h"
52 #include "gromacs/mdtypes/interaction_const.h"
53 #include "gromacs/utility/real.h"
54
55 /* kernel does #include "gromacs/math/utilities.h" */
56 /* Move the actual useful stuff here: */
57
58 //! Define 1/sqrt(pi)
59 #define M_FLOAT_1_SQRTPI 0.564189583547756f
60
61 /*! \brief Macros defining platform-dependent defaults for the prune kernel's j4 processing concurrency.
62  *
63  *  The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
64  */
65 /*! @{ */
66 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
67 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_AMD       4
68 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_NVIDIA    4
69 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_DEFAULT   4
70 #else
71 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_AMD       GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
72 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_NVIDIA    GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
73 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_DEFAULT   GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
74 #endif
75 /*! @} */
76 /*! \brief Constants for platform-dependent defaults for the prune kernel's j4 processing concurrency.
77  *
78  *  Initialized using macros that can be overridden at compile-time (using #GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY).
79  */
80 /*! @{ */
81 const int c_oclPruneKernelJ4ConcurrencyAMD     = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_AMD;
82 const int c_oclPruneKernelJ4ConcurrencyNVIDIA  = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_NVIDIA;
83 const int c_oclPruneKernelJ4ConcurrencyDefault = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY_DEFAULT;
84 /*! @} */
85
86 /*! \brief Returns the j4 processing concurrency parameter for the vendor \p vendorId
87  *  \param vendorId takes values from #ocl_vendor_id_t.
88  */
89 static inline int getOclPruneKernelJ4Concurrency(int vendorId)
90 {
91     assert(vendorId < OCL_VENDOR_UNKNOWN);
92     switch (vendorId)
93     {
94         case OCL_VENDOR_AMD:    return c_oclPruneKernelJ4ConcurrencyAMD;     break;
95         case OCL_VENDOR_NVIDIA: return c_oclPruneKernelJ4ConcurrencyNVIDIA;  break;
96         default:                return c_oclPruneKernelJ4ConcurrencyDefault; break;
97     }
98 }
99
100
101 #ifdef __cplusplus
102 extern "C" {
103 #endif
104
105 /*! \brief Electrostatic OpenCL kernel flavors.
106  *
107  *  Types of electrostatics implementations available in the OpenCL non-bonded
108  *  force kernels. These represent both the electrostatics types implemented
109  *  by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
110  *  enums.h) as well as encode implementation details analytical/tabulated
111  *  and single or twin cut-off (for Ewald kernels).
112  *  Note that the cut-off and RF kernels have only analytical flavor and unlike
113  *  in the CPU kernels, the tabulated kernels are ATM Ewald-only.
114  *
115  *  The row-order of pointers to different electrostatic kernels defined in
116  *  nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
117  *  should match the order of enumerated types below.
118  */
119 enum eelOcl {
120     eelOclCUT, eelOclRF, eelOclEWALD_TAB, eelOclEWALD_TAB_TWIN, eelOclEWALD_ANA, eelOclEWALD_ANA_TWIN, eelOclNR
121 };
122
123 /*! \brief VdW OpenCL kernel flavors.
124  *
125  * The enumerates values correspond to the LJ implementations in the OpenCL non-bonded
126  * kernels.
127  *
128  * The column-order of pointers to different electrostatic kernels defined in
129  * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
130  * should match the order of enumerated types below.
131  */
132 enum evdwOcl {
133     evdwOclCUT, evdwOclCUTCOMBGEOM, evdwOclCUTCOMBLB, evdwOclFSWITCH, evdwOclPSWITCH, evdwOclEWALDGEOM, evdwOclEWALDLB, evdwOclNR
134 };
135
136 /*! \brief Pruning kernel flavors.
137  *
138  * The values correspond to the first call of the pruning post-list generation
139  * and the rolling pruning, respectively.
140  */
141 enum ePruneKind {
142     epruneFirst, epruneRolling, ePruneNR
143 };
144
145 /*! \internal
146  * \brief Staging area for temporary data downloaded from the GPU.
147  *
148  *  The energies/shift forces get downloaded here first, before getting added
149  *  to the CPU-side aggregate values.
150  */
151 typedef struct cl_nb_staging
152 {
153     float    *e_lj;           /**< LJ energy                       */
154     float    *e_el;           /**< electrostatic energy            */
155     float   (*fshift)[3];     /**< float3 buffer with shift forces */
156 } cl_nb_staging_t;
157
158 /*! \internal
159  * \brief Nonbonded atom data - both inputs and outputs.
160  */
161 typedef struct cl_atomdata
162 {
163     int         natoms;              /**< number of atoms                              */
164     int         natoms_local;        /**< number of local atoms                        */
165     int         nalloc;              /**< allocation size for the atom data (xq, f)    */
166
167     cl_mem      xq;                  /**< float4 buffer with atom coordinates + charges, size natoms */
168
169     cl_mem      f;                   /**< float3 buffer with force output array, size natoms         */
170     size_t      f_elem_size;         /**< Size in bytes for one element of f buffer      */
171
172     cl_mem      e_lj;                /**< LJ energy output, size 1                       */
173     cl_mem      e_el;                /**< Electrostatics energy input, size 1            */
174
175     cl_mem      fshift;              /**< float3 buffer with shift forces                */
176     size_t      fshift_elem_size;    /**< Size in bytes for one element of fshift buffer */
177
178     int         ntypes;              /**< number of atom types                           */
179     cl_mem      atom_types;          /**< int buffer with atom type indices, size natoms */
180     cl_mem      lj_comb;             /**< float2 buffer with sqrt(c6),sqrt(c12), size natoms */
181
182     cl_mem      shift_vec;           /**< float3 buffer with shifts values               */
183     size_t      shift_vec_elem_size; /**< Size in bytes for one element of shift_vec buffer */
184
185     cl_bool     bShiftVecUploaded;   /**< true if the shift vector has been uploaded  */
186 } cl_atomdata_t;
187
188 /*! \internal
189  * \brief Parameters required for the OpenCL nonbonded calculations.
190  */
191 typedef struct cl_nbparam
192 {
193
194     int             eeltype;           /**< type of electrostatics, takes values from #eelOcl */
195     int             vdwtype;           /**< type of VdW impl., takes values from #evdwOcl     */
196
197     float           epsfac;            /**< charge multiplication factor                      */
198     float           c_rf;              /**< Reaction-field/plain cutoff electrostatics const. */
199     float           two_k_rf;          /**< Reaction-field electrostatics constant            */
200     float           ewald_beta;        /**< Ewald/PME parameter                               */
201     float           sh_ewald;          /**< Ewald/PME correction term substracted from the direct-space potential */
202     float           sh_lj_ewald;       /**< LJ-Ewald/PME correction term added to the correction potential        */
203     float           ewaldcoeff_lj;     /**< LJ-Ewald/PME coefficient                          */
204
205     float           rcoulomb_sq;       /**< Coulomb cut-off squared                           */
206
207     float           rvdw_sq;           /**< VdW cut-off squared                               */
208     float           rvdw_switch;       /**< VdW switched cut-off                              */
209     float           rlistOuter_sq;     /**< Full, outer pair-list cut-off squared             */
210     float           rlistInner_sq;     /**< Inner, dynamic pruned pair-list cut-off squared   */
211     bool            useDynamicPruning; /**< True if we use dynamic pair-list pruning          */
212
213     shift_consts_t  dispersion_shift;  /**< VdW shift dispersion constants           */
214     shift_consts_t  repulsion_shift;   /**< VdW shift repulsion constants            */
215     switch_consts_t vdw_switch;        /**< VdW switch constants                     */
216
217     /* LJ non-bonded parameters - accessed through texture memory */
218     cl_mem                  nbfp_climg2d;      /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
219     cl_mem                  nbfp_comb_climg2d; /**< nonbonded parameter table per atom type, 2*ntype elements                          */
220
221     /* Ewald Coulomb force table data - accessed through texture memory */
222     float                  coulomb_tab_scale;   /**< table scale/spacing                        */
223     cl_mem                 coulomb_tab_climg2d; /**< pointer to the table in the device memory  */
224 } cl_nbparam_t;
225
226 /*! \internal
227  * \brief Data structure shared between the OpenCL device code and OpenCL host code
228  *
229  * Must not contain OpenCL objects (buffers)
230  * TODO: review, improve */
231 typedef struct cl_nbparam_params
232 {
233
234     int             eeltype;          /**< type of electrostatics, takes values from #eelCu */
235     int             vdwtype;          /**< type of VdW impl., takes values from #evdwCu     */
236
237     float           epsfac;           /**< charge multiplication factor                      */
238     float           c_rf;             /**< Reaction-field/plain cutoff electrostatics const. */
239     float           two_k_rf;         /**< Reaction-field electrostatics constant            */
240     float           ewald_beta;       /**< Ewald/PME parameter                               */
241     float           sh_ewald;         /**< Ewald/PME correction term substracted from the direct-space potential */
242     float           sh_lj_ewald;      /**< LJ-Ewald/PME correction term added to the correction potential        */
243     float           ewaldcoeff_lj;    /**< LJ-Ewald/PME coefficient                          */
244
245     float           rcoulomb_sq;      /**< Coulomb cut-off squared                           */
246
247     float           rvdw_sq;          /**< VdW cut-off squared                               */
248     float           rvdw_switch;      /**< VdW switched cut-off                              */
249     float           rlistOuter_sq;    /**< Full, outer pair-list cut-off squared             */
250     float           rlistInner_sq;    /**< Inner, dynamic pruned pair-list cut-off squared   */
251
252     shift_consts_t  dispersion_shift; /**< VdW shift dispersion constants           */
253     shift_consts_t  repulsion_shift;  /**< VdW shift repulsion constants            */
254     switch_consts_t vdw_switch;       /**< VdW switch constants                     */
255
256     /* Ewald Coulomb force table data - accessed through texture memory */
257     float                  coulomb_tab_scale;  /**< table scale/spacing                        */
258 } cl_nbparam_params_t;
259
260
261 /*! \internal
262  * \brief Pair list data.
263  */
264 typedef struct cl_plist
265 {
266     int              na_c;         /**< number of atoms per cluster                  */
267
268     int              nsci;         /**< size of sci, # of i clusters in the list     */
269     int              sci_nalloc;   /**< allocation size of sci                       */
270     cl_mem           sci;          /**< list of i-cluster ("super-clusters").
271                                         It contains elements of type nbnxn_sci_t     */
272
273     int              ncj4;         /**< total # of 4*j clusters                      */
274     int              cj4_nalloc;   /**< allocation size of cj4                       */
275     cl_mem           cj4;          /**< 4*j cluster list, contains j cluster number and
276                                         index into the i cluster list.
277                                         It contains elements of type nbnxn_cj4_t     */
278     int              nimask;       /**< # of 4*j clusters * # of warps               */
279     int              imask_nalloc; /**< allocation size of imask                     */
280     cl_mem           imask;        /**< imask for 2 warps for each 4*j cluster group */
281     cl_mem           excl;         /**< atom interaction bits
282                                         It contains elements of type nbnxn_excl_t    */
283     int              nexcl;        /**< count for excl                               */
284     int              excl_nalloc;  /**< allocation size of excl                      */
285
286     /* parameter+variables for normal and rolling pruning */
287     bool             haveFreshList;          /**< true after search, indictes that initial pruning with outer prunning is needed */
288     int              rollingPruningNumParts; /**< the number of parts/steps over which one cyle of roling pruning takes places */
289     int              rollingPruningPart;     /**< the next part to which the roling pruning needs to be applied */
290 }cl_plist_t;
291
292
293 /** \internal
294  * \brief Typedef of actual timer type.
295  */
296 typedef struct nbnxn_gpu_timers_t cl_timers_t;
297
298 /*! \internal
299  * \brief Main data structure for OpenCL nonbonded force calculations.
300  */
301 struct gmx_nbnxn_ocl_t
302 {
303     const gmx_device_info_t          *dev_info;    /**< OpenCL device information                              */
304     struct gmx_device_runtime_data_t *dev_rundata; /**< OpenCL runtime data (context, kernels)                 */
305
306     /**< Pointers to non-bonded kernel functions
307      * organized similar with nb_kfunc_xxx arrays in nbnxn_ocl.cpp */
308     ///@{
309     cl_kernel           kernel_noener_noprune_ptr[eelOclNR][evdwOclNR];
310     cl_kernel           kernel_ener_noprune_ptr[eelOclNR][evdwOclNR];
311     cl_kernel           kernel_noener_prune_ptr[eelOclNR][evdwOclNR];
312     cl_kernel           kernel_ener_prune_ptr[eelOclNR][evdwOclNR];
313     ///@}
314     cl_kernel           kernel_pruneonly[ePruneNR]; /**< prune kernels, ePruneKind defined the kernel kinds */
315
316     bool                bPrefetchLjParam;           /**< true if prefetching fg i-atom LJ parameters should be used in the kernels */
317
318     /**< auxiliary kernels implementing memset-like functions */
319     ///@{
320     cl_kernel           kernel_memset_f;
321     cl_kernel           kernel_memset_f2;
322     cl_kernel           kernel_memset_f3;
323     cl_kernel           kernel_zero_e_fshift;
324     ///@}
325
326     cl_bool             bUseTwoStreams;        /**< true if doing both local/non-local NB work on GPU          */
327     cl_bool             bNonLocalStreamActive; /**< true indicates that the nonlocal_done event was enqueued   */
328
329     cl_atomdata_t      *atdat;                 /**< atom data                                                  */
330     cl_nbparam_t       *nbparam;               /**< parameters required for the non-bonded calc.               */
331     cl_plist_t         *plist[2];              /**< pair-list data structures (local and non-local)            */
332     cl_nb_staging_t     nbst;                  /**< staging area where fshift/energies get downloaded          */
333
334     cl_mem              debug_buffer;          /**< debug buffer */
335
336     cl_command_queue    stream[2];             /**< local and non-local GPU queues                             */
337
338     /** events used for synchronization */
339     cl_event nonlocal_done;                     /**< event triggered when the non-local non-bonded kernel
340                                                    is done (and the local transfer can proceed) */
341     cl_event misc_ops_and_local_H2D_done;       /**< event triggered when the tasks issued in
342                                                    the local stream that need to precede the
343                                                    non-local force calculations are done
344                                                    (e.g. f buffer 0-ing, local x/q H2D) */
345
346     cl_bool                           bDoTime;  /**< True if event-based timing is enabled.                     */
347     cl_timers_t                      *timers;   /**< OpenCL event-based timers.                                 */
348     struct gmx_wallclock_gpu_nbnxn_t *timings;  /**< Timing data. TODO: deprecate this and query timers for accumulated data instead */
349 };
350
351 #ifdef __cplusplus
352 }
353 #endif
354
355 #endif  /* NBNXN_OPENCL_TYPES_H */