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