Implement OpenCL support
[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, 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  *  \ingroup module_mdlib
42  */
43
44 #ifndef NBNXN_OPENCL_TYPES_H
45 #define NBNXN_OPENCL_TYPES_H
46
47 #ifdef __APPLE__
48 #    include <OpenCL/opencl.h>
49 #else
50 #    include <CL/opencl.h>
51 #endif
52
53 #include "gromacs/legacyheaders/types/interaction_const.h"
54 #include "gromacs/mdlib/nbnxn_pairlist.h"
55 #include "gromacs/utility/real.h"
56
57 /* kernel does #include "gromacs/math/utilities.h" */
58 /* Move the actual useful stuff here: */
59
60 //! Define 1/sqrt(pi)
61 #define M_FLOAT_1_SQRTPI 0.564189583547756f
62
63
64 #ifdef __cplusplus
65 extern "C" {
66 #endif
67
68 /*! \brief Electrostatic OpenCL kernel flavors.
69  *
70  *  Types of electrostatics implementations available in the OpenCL non-bonded
71  *  force kernels. These represent both the electrostatics types implemented
72  *  by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
73  *  enums.h) as well as encode implementation details analytical/tabulated
74  *  and single or twin cut-off (for Ewald kernels).
75  *  Note that the cut-off and RF kernels have only analytical flavor and unlike
76  *  in the CPU kernels, the tabulated kernels are ATM Ewald-only.
77  *
78  *  The row-order of pointers to different electrostatic kernels defined in
79  *  nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
80  *  should match the order of enumerated types below.
81  */
82 enum eelOcl {
83     eelOclCUT, eelOclRF, eelOclEWALD_TAB, eelOclEWALD_TAB_TWIN, eelOclEWALD_ANA, eelOclEWALD_ANA_TWIN, eelOclNR
84 };
85
86 /*! \brief VdW OpenCL kernel flavors.
87  *
88  * The enumerates values correspond to the LJ implementations in the OpenCL non-bonded
89  * kernels.
90  *
91  * The column-order of pointers to different electrostatic kernels defined in
92  * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
93  * should match the order of enumerated types below.
94  */
95 enum evdwOcl {
96     evdwOclCUT, evdwOclFSWITCH, evdwOclPSWITCH, evdwOclEWALDGEOM, evdwOclEWALDLB, evdwOclNR
97 };
98
99 /*! \internal
100  * \brief Staging area for temporary data downloaded from the GPU.
101  *
102  *  The energies/shift forces get downloaded here first, before getting added
103  *  to the CPU-side aggregate values.
104  */
105 typedef struct cl_nb_staging
106 {
107     float    *e_lj;           /**< LJ energy                       */
108     float    *e_el;           /**< electrostatic energy            */
109     float   (*fshift)[3];     /**< float3 buffer with shift forces */
110 } cl_nb_staging_t;
111
112 /*! \internal
113  * \brief Nonbonded atom data - both inputs and outputs.
114  */
115 typedef struct cl_atomdata
116 {
117     int         natoms;              /**< number of atoms                              */
118     int         natoms_local;        /**< number of local atoms                        */
119     int         nalloc;              /**< allocation size for the atom data (xq, f)    */
120
121     cl_mem      xq;                  /**< float4 buffer with atom coordinates + charges, size natoms */
122
123     cl_mem      f;                   /**< float3 buffer with force output array, size natoms         */
124     size_t      f_elem_size;         /**< Size in bytes for one element of f buffer      */
125
126     cl_mem      e_lj;                /**< LJ energy output, size 1                       */
127     cl_mem      e_el;                /**< Electrostatics energy input, size 1            */
128
129     cl_mem      fshift;              /**< float3 buffer with shift forces                */
130     size_t      fshift_elem_size;    /**< Size in bytes for one element of fshift buffer */
131
132     int         ntypes;              /**< number of atom types                           */
133     cl_mem      atom_types;          /**< int buffer with atom type indices, size natoms */
134
135     cl_mem      shift_vec;           /**< float3 buffer with shifts values               */
136     size_t      shift_vec_elem_size; /**< Size in bytes for one element of shift_vec buffer */
137
138     cl_bool     bShiftVecUploaded;   /**< true if the shift vector has been uploaded  */
139 } cl_atomdata_t;
140
141 /*! \internal
142  * \brief Parameters required for the OpenCL nonbonded calculations.
143  */
144 typedef struct cl_nbparam
145 {
146
147     int             eeltype;          /**< type of electrostatics, takes values from #eelOcl */
148     int             vdwtype;          /**< type of VdW impl., takes values from #evdwOcl     */
149
150     float           epsfac;           /**< charge multiplication factor                      */
151     float           c_rf;             /**< Reaction-field/plain cutoff electrostatics const. */
152     float           two_k_rf;         /**< Reaction-field electrostatics constant            */
153     float           ewald_beta;       /**< Ewald/PME parameter                               */
154     float           sh_ewald;         /**< Ewald/PME correction term substracted from the direct-space potential */
155     float           sh_lj_ewald;      /**< LJ-Ewald/PME correction term added to the correction potential        */
156     float           ewaldcoeff_lj;    /**< LJ-Ewald/PME coefficient                          */
157
158     float           rcoulomb_sq;      /**< Coulomb cut-off squared                           */
159
160     float           rvdw_sq;          /**< VdW cut-off squared                               */
161     float           rvdw_switch;      /**< VdW switched cut-off                              */
162     float           rlist_sq;         /**< pair-list cut-off squared                         */
163
164     shift_consts_t  dispersion_shift; /**< VdW shift dispersion constants           */
165     shift_consts_t  repulsion_shift;  /**< VdW shift repulsion constants            */
166     switch_consts_t vdw_switch;       /**< VdW switch constants                     */
167
168     /* LJ non-bonded parameters - accessed through texture memory */
169     cl_mem                  nbfp_climg2d;      /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
170     cl_mem                  nbfp_comb_climg2d; /**< nonbonded parameter table per atom type, 2*ntype elements                          */
171
172     /* Ewald Coulomb force table data - accessed through texture memory */
173     int                    coulomb_tab_size;    /**< table size (s.t. it fits in texture cache) */
174     float                  coulomb_tab_scale;   /**< table scale/spacing                        */
175     cl_mem                 coulomb_tab_climg2d; /**< pointer to the table in the device memory  */
176 } cl_nbparam_t;
177
178 /*! \internal
179  * \brief Data structure shared between the OpenCL device code and OpenCL host code
180  *
181  * Must not contain OpenCL objects (buffers)
182  * TODO: review, improve */
183 typedef struct cl_nbparam_params
184 {
185
186     int             eeltype;          /**< type of electrostatics, takes values from #eelCu */
187     int             vdwtype;          /**< type of VdW impl., takes values from #evdwCu     */
188
189     float           epsfac;           /**< charge multiplication factor                      */
190     float           c_rf;             /**< Reaction-field/plain cutoff electrostatics const. */
191     float           two_k_rf;         /**< Reaction-field electrostatics constant            */
192     float           ewald_beta;       /**< Ewald/PME parameter                               */
193     float           sh_ewald;         /**< Ewald/PME correction term substracted from the direct-space potential */
194     float           sh_lj_ewald;      /**< LJ-Ewald/PME correction term added to the correction potential        */
195     float           ewaldcoeff_lj;    /**< LJ-Ewald/PME coefficient                          */
196
197     float           rcoulomb_sq;      /**< Coulomb cut-off squared                           */
198
199     float           rvdw_sq;          /**< VdW cut-off squared                               */
200     float           rvdw_switch;      /**< VdW switched cut-off                              */
201     float           rlist_sq;         /**< pair-list cut-off squared                         */
202
203     shift_consts_t  dispersion_shift; /**< VdW shift dispersion constants           */
204     shift_consts_t  repulsion_shift;  /**< VdW shift repulsion constants            */
205     switch_consts_t vdw_switch;       /**< VdW switch constants                     */
206
207     /* Ewald Coulomb force table data - accessed through texture memory */
208     int                    coulomb_tab_size;   /**< table size (s.t. it fits in texture cache) */
209     float                  coulomb_tab_scale;  /**< table scale/spacing                        */
210 } cl_nbparam_params_t;
211
212
213 /*! \internal
214  * \brief Pair list data.
215  */
216 typedef struct cl_plist
217 {
218     int              na_c;        /**< number of atoms per cluster                  */
219
220     int              nsci;        /**< size of sci, # of i clusters in the list     */
221     int              sci_nalloc;  /**< allocation size of sci                       */
222     cl_mem           sci;         /**< list of i-cluster ("super-clusters").
223                                        It contains elements of type nbnxn_sci_t     */
224
225     int              ncj4;        /**< total # of 4*j clusters                      */
226     int              cj4_nalloc;  /**< allocation size of cj4                       */
227     cl_mem           cj4;         /**< 4*j cluster list, contains j cluster number and
228                                        index into the i cluster list.
229                                        It contains elements of type nbnxn_cj4_t     */
230     cl_mem           excl;        /**< atom interaction bits
231                                        It contains elements of type nbnxn_excl_t    */
232     int              nexcl;       /**< count for excl                               */
233     int              excl_nalloc; /**< allocation size of excl                      */
234
235     cl_bool          bDoPrune;    /**< true if pair-list pruning needs to be
236                                        done during the  current step                */
237 }cl_plist_t;
238
239
240 /*! \internal
241  * \brief OpenCL events used for timing GPU kernels and H2D/D2H transfers.
242  *
243  * The two-sized arrays hold the local and non-local values and should always
244  * be indexed with eintLocal/eintNonlocal.
245  */
246 typedef struct cl_timers
247 {
248     cl_event atdat;             /**< event for atom data transfer (every PS step)                 */
249
250     cl_event nb_h2d[2];         /**< events for x/q H2D transfers (l/nl, every step)              */
251
252     cl_event nb_d2h_f[2];       /**< events for f D2H transfer (l/nl, every step)                 */
253     cl_event nb_d2h_fshift[2];  /**< events for fshift D2H transfer (l/nl, every step)            */
254     cl_event nb_d2h_e_el[2];    /**< events for e_el D2H transfer (l/nl, every step)              */
255     cl_event nb_d2h_e_lj[2];    /**< events for e_lj D2H transfer (l/nl, every step)              */
256
257     cl_event pl_h2d_sci[2];     /**< events for pair-list sci H2D transfers (l/nl, every PS step) */
258     cl_event pl_h2d_cj4[2];     /**< events for pair-list cj4 H2D transfers (l/nl, every PS step) */
259     cl_event pl_h2d_excl[2];    /**< events for pair-list excl H2D transfers (l/nl, every PS step)*/
260
261     cl_event nb_k[2];           /**< event for non-bonded kernels (l/nl, every step)              */
262 } cl_timers_t;
263
264 /*! \internal
265  * \brief Main data structure for OpenCL nonbonded force calculations.
266  */
267 struct gmx_nbnxn_ocl_t
268 {
269     struct gmx_device_info_t *dev_info;        /**< OpenCL device information                                  */
270
271     /**< Pointers to non-bonded kernel functions
272      * organized similar with nb_kfunc_xxx arrays in nbnxn_ocl.cpp */
273     ///@{
274     cl_kernel           kernel_noener_noprune_ptr[eelOclNR][evdwOclNR];
275     cl_kernel           kernel_ener_noprune_ptr[eelOclNR][evdwOclNR];
276     cl_kernel           kernel_noener_prune_ptr[eelOclNR][evdwOclNR];
277     cl_kernel           kernel_ener_prune_ptr[eelOclNR][evdwOclNR];
278     ///@}
279
280     /**< auxiliary kernels implementing memset-like functions */
281     ///@{
282     cl_kernel           kernel_memset_f;
283     cl_kernel           kernel_memset_f2;
284     cl_kernel           kernel_memset_f3;
285     cl_kernel           kernel_zero_e_fshift;
286     ///@}
287
288     cl_bool             bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU          */
289
290     cl_atomdata_t      *atdat;          /**< atom data                                                  */
291     cl_nbparam_t       *nbparam;        /**< parameters required for the non-bonded calc.               */
292     cl_plist_t         *plist[2];       /**< pair-list data structures (local and non-local)            */
293     cl_nb_staging_t     nbst;           /**< staging area where fshift/energies get downloaded          */
294
295     cl_mem              debug_buffer;   /**< debug buffer */
296
297     cl_command_queue    stream[2];      /**< local and non-local GPU queues                             */
298
299     /** events used for synchronization */
300     cl_event    nonlocal_done;           /**< event triggered when the non-local non-bonded kernel
301                                               is done (and the local transfer can proceed)               */
302     cl_event    misc_ops_done;           /**< event triggered when the operations that precede the
303                                               main force calculations are done (e.g. buffer 0-ing)       */
304
305     cl_bool                     bDoTime; /**< True if event-based timing is enabled.                     */
306     cl_timers_t                *timers;  /**< OpenCL event-based timers.                                 */
307     struct gmx_wallclock_gpu_t *timings; /**< Timing data.                                               */
308 };
309
310 #ifdef __cplusplus
311 }
312 #endif
313
314 #endif  /* NBNXN_OPENCL_TYPES_H */