367a59c75f712fb03edec238896b01e579fef63a
[alexxy/gromacs.git] / src / gromacs / mdlib / nbnxn_cuda / nbnxn_cuda_types.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2012, The GROMACS development team.
6  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 /*! \internal \file
39  *  \brief
40  *  Data types used internally in the nbnxn_cuda module.
41  *
42  *  \author Szilárd Páll <pall.szilard@gmail.com>
43  *  \ingroup module_mdlib
44  */
45
46 #ifndef NBNXN_CUDA_TYPES_H
47 #define NBNXN_CUDA_TYPES_H
48
49 #include "gromacs/legacyheaders/types/interaction_const.h"
50 #include "gromacs/mdlib/nbnxn_pairlist.h"
51 #include "gromacs/legacyheaders/types/nbnxn_cuda_types_ext.h"
52 #include "../../gmxlib/cuda_tools/cudautils.cuh"
53
54 /* CUDA versions from 5.0 above support texture objects. */
55 #if CUDA_VERSION >= 5000
56 #define TEXOBJ_SUPPORTED
57 #else  /* CUDA_VERSION */
58 /** This typedef allows us to define only one version of struct cu_nbparam */
59 typedef int cudaTextureObject_t;
60 #endif /* CUDA_VERSION */
61
62 #ifdef __cplusplus
63 extern "C" {
64 #endif
65
66 /*! \brief Electrostatic CUDA kernel flavors.
67  *
68  *  Types of electrostatics implementations available in the CUDA non-bonded
69  *  force kernels. These represent both the electrostatics types implemented
70  *  by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
71  *  enums.h) as well as encode implementation details analytical/tabulated
72  *  and single or twin cut-off (for Ewald kernels).
73  *  Note that the cut-off and RF kernels have only analytical flavor and unlike
74  *  in the CPU kernels, the tabulated kernels are ATM Ewald-only.
75  *
76  *  The row-order of pointers to different electrostatic kernels defined in
77  *  nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
78  *  should match the order of enumerated types below.
79  */
80 enum eelCu {
81     eelCuCUT, eelCuRF, eelCuEWALD_TAB, eelCuEWALD_TAB_TWIN, eelCuEWALD_ANA, eelCuEWALD_ANA_TWIN, eelCuNR
82 };
83
84 /*! \brief VdW CUDA kernel flavors.
85  *
86  * The enumerates values correspond to the LJ implementations in the CUDA non-bonded
87  * kernels.
88  *
89  * The column-order of pointers to different electrostatic kernels defined in
90  * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
91  * should match the order of enumerated types below.
92  */
93 enum evdwCu {
94     evdwCuCUT, evdwCuFSWITCH, evdwCuPSWITCH, evdwCuEWALDGEOM, evdwCuEWALDLB, evdwCuNR
95 };
96
97 /* All structs prefixed with "cu_" hold data used in GPU calculations and
98  * are passed to the kernels, except cu_timers_t. */
99 /*! \cond */
100 typedef struct cu_plist     cu_plist_t;
101 typedef struct cu_atomdata  cu_atomdata_t;
102 typedef struct cu_nbparam   cu_nbparam_t;
103 typedef struct cu_timers    cu_timers_t;
104 typedef struct nb_staging   nb_staging_t;
105 /*! \endcond */
106
107
108 /** \internal
109  * \brief Staging area for temporary data downloaded from the GPU.
110  *
111  *  The energies/shift forces get downloaded here first, before getting added
112  *  to the CPU-side aggregate values.
113  */
114 struct nb_staging
115 {
116     float   *e_lj;      /**< LJ energy            */
117     float   *e_el;      /**< electrostatic energy */
118     float3  *fshift;    /**< shift forces         */
119 };
120
121 /** \internal
122  * \brief Nonbonded atom data - both inputs and outputs.
123  */
124 struct cu_atomdata
125 {
126     int      natoms;            /**< number of atoms                              */
127     int      natoms_local;      /**< number of local atoms                        */
128     int      nalloc;            /**< allocation size for the atom data (xq, f)    */
129
130     float4  *xq;                /**< atom coordinates + charges, size natoms      */
131     float3  *f;                 /**< force output array, size natoms              */
132
133     float   *e_lj;              /**< LJ energy output, size 1                     */
134     float   *e_el;              /**< Electrostatics energy input, size 1          */
135
136     float3  *fshift;            /**< shift forces                                 */
137
138     int      ntypes;            /**< number of atom types                         */
139     int     *atom_types;        /**< atom type indices, size natoms               */
140
141     float3  *shift_vec;         /**< shifts                                       */
142     bool     bShiftVecUploaded; /**< true if the shift vector has been uploaded   */
143 };
144
145 /** \internal
146  * \brief Parameters required for the CUDA nonbonded calculations.
147  */
148 struct cu_nbparam
149 {
150
151     int             eeltype;          /**< type of electrostatics, takes values from #eelCu */
152     int             vdwtype;          /**< type of VdW impl., takes values from #evdwCu     */
153
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                          */
161
162     float           rcoulomb_sq;      /**< Coulomb cut-off squared                           */
163
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                         */
167
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                     */
171
172     /* LJ non-bonded parameters - accessed through texture memory */
173     float               *nbfp;             /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
174     cudaTextureObject_t  nbfp_texobj;      /**< texture object bound to nbfp                                                       */
175     float               *nbfp_comb;        /**< nonbonded parameter table per atom type, 2*ntype elements                          */
176     cudaTextureObject_t  nbfp_comb_texobj; /**< texture object bound to nbfp_texobj                                                */
177
178     /* Ewald Coulomb force table data - accessed through texture memory */
179     int                  coulomb_tab_size;   /**< table size (s.t. it fits in texture cache) */
180     float                coulomb_tab_scale;  /**< table scale/spacing                        */
181     float               *coulomb_tab;        /**< pointer to the table in the device memory  */
182     cudaTextureObject_t  coulomb_tab_texobj; /**< texture object bound to coulomb_tab        */
183 };
184
185 /** \internal
186  * \brief Pair list data.
187  */
188 struct cu_plist
189 {
190     int              na_c;        /**< number of atoms per cluster                  */
191
192     int              nsci;        /**< size of sci, # of i clusters in the list     */
193     int              sci_nalloc;  /**< allocation size of sci                       */
194     nbnxn_sci_t     *sci;         /**< list of i-cluster ("super-clusters")         */
195
196     int              ncj4;        /**< total # of 4*j clusters                      */
197     int              cj4_nalloc;  /**< allocation size of cj4                       */
198     nbnxn_cj4_t     *cj4;         /**< 4*j cluster list, contains j cluster number
199                                        and index into the i cluster list            */
200     nbnxn_excl_t    *excl;        /**< atom interaction bits                        */
201     int              nexcl;       /**< count for excl                               */
202     int              excl_nalloc; /**< allocation size of excl                      */
203
204     bool             bDoPrune;    /**< true if pair-list pruning needs to be
205                                        done during the  current step                */
206 };
207
208 /** \internal
209  * \brief CUDA events used for timing GPU kernels and H2D/D2H transfers.
210  *
211  * The two-sized arrays hold the local and non-local values and should always
212  * be indexed with eintLocal/eintNonlocal.
213  */
214 struct cu_timers
215 {
216     cudaEvent_t start_atdat;     /**< start event for atom data transfer (every PS step)             */
217     cudaEvent_t stop_atdat;      /**< stop event for atom data transfer (every PS step)              */
218     cudaEvent_t start_nb_h2d[2]; /**< start events for x/q H2D transfers (l/nl, every step)          */
219     cudaEvent_t stop_nb_h2d[2];  /**< stop events for x/q H2D transfers (l/nl, every step)           */
220     cudaEvent_t start_nb_d2h[2]; /**< start events for f D2H transfer (l/nl, every step)             */
221     cudaEvent_t stop_nb_d2h[2];  /**< stop events for f D2H transfer (l/nl, every step)              */
222     cudaEvent_t start_pl_h2d[2]; /**< start events for pair-list H2D transfers (l/nl, every PS step) */
223     cudaEvent_t stop_pl_h2d[2];  /**< start events for pair-list H2D transfers (l/nl, every PS step) */
224     cudaEvent_t start_nb_k[2];   /**< start event for non-bonded kernels (l/nl, every step)          */
225     cudaEvent_t stop_nb_k[2];    /**< stop event non-bonded kernels (l/nl, every step)               */
226 };
227
228 /** \internal
229  * \brief Main data structure for CUDA nonbonded force calculations.
230  */
231 struct nbnxn_cuda
232 {
233     cuda_dev_info_t *dev_info;       /**< CUDA device information                              */
234     bool             bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU    */
235     bool             bUseStreamSync; /**< true if the standard cudaStreamSynchronize is used
236                                           and not memory polling-based waiting                 */
237     cu_atomdata_t   *atdat;          /**< atom data                                            */
238     cu_nbparam_t    *nbparam;        /**< parameters required for the non-bonded calc.         */
239     cu_plist_t      *plist[2];       /**< pair-list data structures (local and non-local)      */
240     nb_staging_t     nbst;           /**< staging area where fshift/energies get downloaded    */
241
242     cudaStream_t     stream[2];      /**< local and non-local GPU streams                      */
243
244     /** events used for synchronization */
245     cudaEvent_t    nonlocal_done;    /**< event triggered when the non-local non-bonded kernel
246                                         is done (and the local transfer can proceed)           */
247     cudaEvent_t    misc_ops_done;    /**< event triggered when the operations that precede the
248                                           main force calculations are done (e.g. buffer 0-ing) */
249
250     /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
251      * concurrent streams, so we won't time if both l/nl work is done on GPUs.
252      * Timer init/uninit is still done even with timing off so only the condition
253      * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
254     bool             bDoTime;       /**< True if event-based timing is enabled.               */
255     cu_timers_t     *timers;        /**< CUDA event-based timers.                             */
256     wallclock_gpu_t *timings;       /**< Timing data.                                         */
257 };
258
259 #ifdef __cplusplus
260 }
261 #endif
262
263 #endif  /* NBNXN_CUDA_TYPES_H */