609123f346f784ebb27d7eaae10b0808b87846df
[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,2015,2016,2017, 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/gpu_utils/cuda_arch_utils.cuh"
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/mdlib/nbnxn_consts.h"
52 #include "gromacs/mdlib/nbnxn_pairlist.h"
53 #include "gromacs/mdtypes/interaction_const.h"
54 #include "gromacs/timing/gpu_timing.h"
55
56
57 /*! \brief Macro definining default for the prune kernel's j4 processing concurrency.
58  *
59  *  The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
60  */
61 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
62 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY 4
63 #endif
64 /*! \brief Default for the prune kernel's j4 processing concurrency.
65  *
66  *  Initialized using the #GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro which allows compile-time override.
67  */
68 const int c_cudaPruneKernelJ4Concurrency = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY;
69
70 /* TODO: consider moving this to kernel_utils */
71 /* Convenience defines */
72 /*! \brief number of clusters per supercluster. */
73 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
74 /*! \brief cluster size = number of atoms per cluster. */
75 static const int c_clSize          = c_nbnxnGpuClusterSize;
76
77 /*! \brief True if the use of texture fetch in the CUDA kernels is disabled. */
78 static const bool c_disableCudaTextures = DISABLE_CUDA_TEXTURES;
79
80
81 #ifdef __cplusplus
82 extern "C" {
83 #endif
84
85 /*! \brief Electrostatic CUDA kernel flavors.
86  *
87  *  Types of electrostatics implementations available in the CUDA non-bonded
88  *  force kernels. These represent both the electrostatics types implemented
89  *  by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
90  *  enums.h) as well as encode implementation details analytical/tabulated
91  *  and single or twin cut-off (for Ewald kernels).
92  *  Note that the cut-off and RF kernels have only analytical flavor and unlike
93  *  in the CPU kernels, the tabulated kernels are ATM Ewald-only.
94  *
95  *  The row-order of pointers to different electrostatic kernels defined in
96  *  nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
97  *  should match the order of enumerated types below.
98  */
99 enum eelCu {
100     eelCuCUT, eelCuRF, eelCuEWALD_TAB, eelCuEWALD_TAB_TWIN, eelCuEWALD_ANA, eelCuEWALD_ANA_TWIN, eelCuNR
101 };
102
103 /*! \brief VdW CUDA kernel flavors.
104  *
105  * The enumerates values correspond to the LJ implementations in the CUDA non-bonded
106  * kernels.
107  *
108  * The column-order of pointers to different electrostatic kernels defined in
109  * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
110  * should match the order of enumerated types below.
111  */
112 enum evdwCu {
113     evdwCuCUT, evdwCuCUTCOMBGEOM, evdwCuCUTCOMBLB, evdwCuFSWITCH, evdwCuPSWITCH, evdwCuEWALDGEOM, evdwCuEWALDLB, evdwCuNR
114 };
115
116 /* All structs prefixed with "cu_" hold data used in GPU calculations and
117  * are passed to the kernels, except cu_timers_t. */
118 /*! \cond */
119 typedef struct cu_plist     cu_plist_t;
120 typedef struct cu_atomdata  cu_atomdata_t;
121 typedef struct cu_nbparam   cu_nbparam_t;
122 typedef struct cu_timers    cu_timers_t;
123 typedef struct nb_staging   nb_staging_t;
124 /*! \endcond */
125
126
127 /** \internal
128  * \brief Staging area for temporary data downloaded from the GPU.
129  *
130  *  The energies/shift forces get downloaded here first, before getting added
131  *  to the CPU-side aggregate values.
132  */
133 struct nb_staging
134 {
135     float   *e_lj;      /**< LJ energy            */
136     float   *e_el;      /**< electrostatic energy */
137     float3  *fshift;    /**< shift forces         */
138 };
139
140 /** \internal
141  * \brief Nonbonded atom data - both inputs and outputs.
142  */
143 struct cu_atomdata
144 {
145     int      natoms;            /**< number of atoms                              */
146     int      natoms_local;      /**< number of local atoms                        */
147     int      nalloc;            /**< allocation size for the atom data (xq, f)    */
148
149     float4  *xq;                /**< atom coordinates + charges, size natoms      */
150     float3  *f;                 /**< force output array, size natoms              */
151
152     float   *e_lj;              /**< LJ energy output, size 1                     */
153     float   *e_el;              /**< Electrostatics energy input, size 1          */
154
155     float3  *fshift;            /**< shift forces                                 */
156
157     int      ntypes;            /**< number of atom types                         */
158     int     *atom_types;        /**< atom type indices, size natoms               */
159     float2  *lj_comb;           /**< sqrt(c6),sqrt(c12) size natoms               */
160
161     float3  *shift_vec;         /**< shifts                                       */
162     bool     bShiftVecUploaded; /**< true if the shift vector has been uploaded   */
163 };
164
165 /** \internal
166  * \brief Parameters required for the CUDA nonbonded calculations.
167  */
168 struct cu_nbparam
169 {
170
171     int             eeltype;              /**< type of electrostatics, takes values from #eelCu */
172     int             vdwtype;              /**< type of VdW impl., takes values from #evdwCu     */
173
174     float           epsfac;               /**< charge multiplication factor                      */
175     float           c_rf;                 /**< Reaction-field/plain cutoff electrostatics const. */
176     float           two_k_rf;             /**< Reaction-field electrostatics constant            */
177     float           ewald_beta;           /**< Ewald/PME parameter                               */
178     float           sh_ewald;             /**< Ewald/PME correction term substracted from the direct-space potential */
179     float           sh_lj_ewald;          /**< LJ-Ewald/PME correction term added to the correction potential        */
180     float           ewaldcoeff_lj;        /**< LJ-Ewald/PME coefficient                          */
181
182     float           rcoulomb_sq;          /**< Coulomb cut-off squared                           */
183
184     float           rvdw_sq;              /**< VdW cut-off squared                               */
185     float           rvdw_switch;          /**< VdW switched cut-off                              */
186     float           rlistOuter_sq;        /**< Full, outer pair-list cut-off squared             */
187     float           rlistInner_sq;        /**< Inner, dynamic pruned pair-list cut-off squared   */
188     bool            useDynamicPruning;    /**< True if we use dynamic pair-list pruning          */
189
190     shift_consts_t  dispersion_shift;     /**< VdW shift dispersion constants           */
191     shift_consts_t  repulsion_shift;      /**< VdW shift repulsion constants            */
192     switch_consts_t vdw_switch;           /**< VdW switch constants                     */
193
194     /* LJ non-bonded parameters - accessed through texture memory */
195     float               *nbfp;             /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
196     cudaTextureObject_t  nbfp_texobj;      /**< texture object bound to nbfp                                                       */
197     float               *nbfp_comb;        /**< nonbonded parameter table per atom type, 2*ntype elements                          */
198     cudaTextureObject_t  nbfp_comb_texobj; /**< texture object bound to nbfp_texobj                                                */
199
200     /* Ewald Coulomb force table data - accessed through texture memory */
201     float                coulomb_tab_scale;  /**< table scale/spacing                        */
202     float               *coulomb_tab;        /**< pointer to the table in the device memory  */
203     cudaTextureObject_t  coulomb_tab_texobj; /**< texture object bound to coulomb_tab        */
204 };
205
206 /** \internal
207  * \brief Pair list data.
208  */
209 struct cu_plist
210 {
211     int              na_c;         /**< number of atoms per cluster                  */
212
213     int              nsci;         /**< size of sci, # of i clusters in the list     */
214     int              sci_nalloc;   /**< allocation size of sci                       */
215     nbnxn_sci_t     *sci;          /**< list of i-cluster ("super-clusters")         */
216
217     int              ncj4;         /**< total # of 4*j clusters                      */
218     int              cj4_nalloc;   /**< allocation size of cj4                       */
219     nbnxn_cj4_t     *cj4;          /**< 4*j cluster list, contains j cluster number
220                                         and index into the i cluster list            */
221     int              nimask;       /**< # of 4*j clusters * # of warps               */
222     int              imask_nalloc; /**< allocation size of imask                     */
223     unsigned int    *imask;        /**< imask for 2 warps for each 4*j cluster group */
224     nbnxn_excl_t    *excl;         /**< atom interaction bits                        */
225     int              nexcl;        /**< count for excl                               */
226     int              excl_nalloc;  /**< allocation size of excl                      */
227
228     /* parameter+variables for normal and rolling pruning */
229     bool             haveFreshList;          /**< true after search, indictes that initial pruning with outer prunning is needed */
230     int              rollingPruningNumParts; /**< the number of parts/steps over which one cyle of roling pruning takes places */
231     int              rollingPruningPart;     /**< the next part to which the roling pruning needs to be applied */
232 };
233
234 /** \internal
235  * \brief CUDA events used for timing GPU kernels and H2D/D2H transfers.
236  *
237  * The two-sized arrays hold the local and non-local values and should always
238  * be indexed with eintLocal/eintNonlocal.
239  */
240 struct cu_timers
241 {
242     cudaEvent_t start_atdat;             /**< start event for atom data transfer (every PS step)             */
243     cudaEvent_t stop_atdat;              /**< stop event for atom data transfer (every PS step)              */
244     cudaEvent_t start_nb_h2d[2];         /**< start events for x/q H2D transfers (l/nl, every step)          */
245     cudaEvent_t stop_nb_h2d[2];          /**< stop events for x/q H2D transfers (l/nl, every step)           */
246     cudaEvent_t start_nb_d2h[2];         /**< start events for f D2H transfer (l/nl, every step)             */
247     cudaEvent_t stop_nb_d2h[2];          /**< stop events for f D2H transfer (l/nl, every step)              */
248     cudaEvent_t start_pl_h2d[2];         /**< start events for pair-list H2D transfers (l/nl, every PS step) */
249     cudaEvent_t stop_pl_h2d[2];          /**< start events for pair-list H2D transfers (l/nl, every PS step) */
250     bool        didPairlistH2D[2];       /**< true when a pair-list transfer has been done at this step      */
251     cudaEvent_t start_nb_k[2];           /**< start event for non-bonded kernels (l/nl, every step)          */
252     cudaEvent_t stop_nb_k[2];            /**< stop event non-bonded kernels (l/nl, every step)               */
253     cudaEvent_t start_prune_k[2];        /**< start event for the 1st pass list pruning kernel (l/nl, every PS step)   */
254     cudaEvent_t stop_prune_k[2];         /**< stop event for the 1st pass list pruning kernel (l/nl, every PS step)   */
255     bool        didPrune[2];             /**< true when we timed pruning and the timings need to be accounted for */
256     cudaEvent_t start_rollingPrune_k[2]; /**< start event for rolling pruning kernels (l/nl, frequency depends on chunk size)   */
257     cudaEvent_t stop_rollingPrune_k[2];  /**< stop event for rolling pruning kernels (l/nl, frequency depends on chunk size)   */
258     bool        didRollingPrune[2];      /**< true when we timed rolling pruning (at the previous step) and the timings need to be accounted for */
259 };
260
261 /** \internal
262  * \brief Main data structure for CUDA nonbonded force calculations.
263  */
264 struct gmx_nbnxn_cuda_t
265 {
266     const gmx_device_info_t  *dev_info;       /**< CUDA device information                              */
267     bool                      bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU    */
268     cu_atomdata_t            *atdat;          /**< atom data                                            */
269     cu_nbparam_t             *nbparam;        /**< parameters required for the non-bonded calc.         */
270     cu_plist_t               *plist[2];       /**< pair-list data structures (local and non-local)      */
271     nb_staging_t              nbst;           /**< staging area where fshift/energies get downloaded    */
272
273     cudaStream_t              stream[2];      /**< local and non-local GPU streams                      */
274
275     /** events used for synchronization */
276     cudaEvent_t    nonlocal_done;               /**< event triggered when the non-local non-bonded kernel
277                                                    is done (and the local transfer can proceed)           */
278     cudaEvent_t    misc_ops_and_local_H2D_done; /**< event triggered when the tasks issued in
279                                                    the local stream that need to precede the
280                                                    non-local force calculations are done
281                                                    (e.g. f buffer 0-ing, local x/q H2D) */
282
283     /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
284      * concurrent streams, so we won't time if both l/nl work is done on GPUs.
285      * Timer init/uninit is still done even with timing off so only the condition
286      * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
287     bool                 bDoTime;   /**< True if event-based timing is enabled.               */
288     cu_timers_t         *timers;    /**< CUDA event-based timers.                             */
289     gmx_wallclock_gpu_t *timings;   /**< Timing data.                                         */
290 };
291
292 #ifdef __cplusplus
293 }
294 #endif
295
296 #endif  /* NBNXN_CUDA_TYPES_H */