c155d5bdba6bc02ec97f6cd9ebb5219868f54994
[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,2018, 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_gpu_types_common.h"
53 #include "gromacs/mdlib/nbnxn_pairlist.h"
54 #include "gromacs/mdtypes/interaction_const.h"
55 #include "gromacs/timing/gpu_timing.h"
56
57
58 /*! \brief Macro definining default for the prune kernel's j4 processing concurrency.
59  *
60  *  The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
61  */
62 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
63 #define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY 4
64 #endif
65 /*! \brief Default for the prune kernel's j4 processing concurrency.
66  *
67  *  Initialized using the #GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro which allows compile-time override.
68  */
69 const int c_cudaPruneKernelJ4Concurrency = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY;
70
71 /* TODO: consider moving this to kernel_utils */
72 /* Convenience defines */
73 /*! \brief number of clusters per supercluster. */
74 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
75 /*! \brief cluster size = number of atoms per cluster. */
76 static const int c_clSize          = c_nbnxnGpuClusterSize;
77
78 #ifdef __cplusplus
79 extern "C" {
80 #endif
81
82 /*! \brief Electrostatic CUDA kernel flavors.
83  *
84  *  Types of electrostatics implementations available in the CUDA non-bonded
85  *  force kernels. These represent both the electrostatics types implemented
86  *  by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
87  *  enums.h) as well as encode implementation details analytical/tabulated
88  *  and single or twin cut-off (for Ewald kernels).
89  *  Note that the cut-off and RF kernels have only analytical flavor and unlike
90  *  in the CPU kernels, the tabulated kernels are ATM Ewald-only.
91  *
92  *  The row-order of pointers to different electrostatic kernels defined in
93  *  nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
94  *  should match the order of enumerated types below.
95  */
96 enum eelCu {
97     eelCuCUT, eelCuRF, eelCuEWALD_TAB, eelCuEWALD_TAB_TWIN, eelCuEWALD_ANA, eelCuEWALD_ANA_TWIN, eelCuNR
98 };
99
100 /*! \brief VdW CUDA kernel flavors.
101  *
102  * The enumerates values correspond to the LJ implementations in the CUDA non-bonded
103  * kernels.
104  *
105  * The column-order of pointers to different electrostatic kernels defined in
106  * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
107  * should match the order of enumerated types below.
108  */
109 enum evdwCu {
110     evdwCuCUT, evdwCuCUTCOMBGEOM, evdwCuCUTCOMBLB, evdwCuFSWITCH, evdwCuPSWITCH, evdwCuEWALDGEOM, evdwCuEWALDLB, evdwCuNR
111 };
112
113 /* All structs prefixed with "cu_" hold data used in GPU calculations and
114  * are passed to the kernels, except cu_timers_t. */
115 /*! \cond */
116 typedef struct cu_atomdata  cu_atomdata_t;
117 typedef struct cu_nbparam   cu_nbparam_t;
118 typedef struct nb_staging   nb_staging_t;
119 /*! \endcond */
120
121
122 /** \internal
123  * \brief Staging area for temporary data downloaded from the GPU.
124  *
125  *  The energies/shift forces get downloaded here first, before getting added
126  *  to the CPU-side aggregate values.
127  */
128 struct nb_staging
129 {
130     float   *e_lj;      /**< LJ energy            */
131     float   *e_el;      /**< electrostatic energy */
132     float3  *fshift;    /**< shift forces         */
133 };
134
135 /** \internal
136  * \brief Nonbonded atom data - both inputs and outputs.
137  */
138 struct cu_atomdata
139 {
140     int      natoms;            /**< number of atoms                              */
141     int      natoms_local;      /**< number of local atoms                        */
142     int      nalloc;            /**< allocation size for the atom data (xq, f)    */
143
144     float4  *xq;                /**< atom coordinates + charges, size natoms      */
145     float3  *f;                 /**< force output array, size natoms              */
146
147     float   *e_lj;              /**< LJ energy output, size 1                     */
148     float   *e_el;              /**< Electrostatics energy input, size 1          */
149
150     float3  *fshift;            /**< shift forces                                 */
151
152     int      ntypes;            /**< number of atom types                         */
153     int     *atom_types;        /**< atom type indices, size natoms               */
154     float2  *lj_comb;           /**< sqrt(c6),sqrt(c12) size natoms               */
155
156     float3  *shift_vec;         /**< shifts                                       */
157     bool     bShiftVecUploaded; /**< true if the shift vector has been uploaded   */
158 };
159
160 /** \internal
161  * \brief Parameters required for the CUDA nonbonded calculations.
162  */
163 struct cu_nbparam
164 {
165
166     int             eeltype;              /**< type of electrostatics, takes values from #eelCu */
167     int             vdwtype;              /**< type of VdW impl., takes values from #evdwCu     */
168
169     float           epsfac;               /**< charge multiplication factor                      */
170     float           c_rf;                 /**< Reaction-field/plain cutoff electrostatics const. */
171     float           two_k_rf;             /**< Reaction-field electrostatics constant            */
172     float           ewald_beta;           /**< Ewald/PME parameter                               */
173     float           sh_ewald;             /**< Ewald/PME correction term substracted from the direct-space potential */
174     float           sh_lj_ewald;          /**< LJ-Ewald/PME correction term added to the correction potential        */
175     float           ewaldcoeff_lj;        /**< LJ-Ewald/PME coefficient                          */
176
177     float           rcoulomb_sq;          /**< Coulomb cut-off squared                           */
178
179     float           rvdw_sq;              /**< VdW cut-off squared                               */
180     float           rvdw_switch;          /**< VdW switched cut-off                              */
181     float           rlistOuter_sq;        /**< Full, outer pair-list cut-off squared             */
182     float           rlistInner_sq;        /**< Inner, dynamic pruned pair-list cut-off squared   */
183     bool            useDynamicPruning;    /**< True if we use dynamic pair-list pruning          */
184
185     shift_consts_t  dispersion_shift;     /**< VdW shift dispersion constants           */
186     shift_consts_t  repulsion_shift;      /**< VdW shift repulsion constants            */
187     switch_consts_t vdw_switch;           /**< VdW switch constants                     */
188
189     /* LJ non-bonded parameters - accessed through texture memory */
190     float               *nbfp;             /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
191     cudaTextureObject_t  nbfp_texobj;      /**< texture object bound to nbfp                                                       */
192     float               *nbfp_comb;        /**< nonbonded parameter table per atom type, 2*ntype elements                          */
193     cudaTextureObject_t  nbfp_comb_texobj; /**< texture object bound to nbfp_texobj                                                */
194
195     /* Ewald Coulomb force table data - accessed through texture memory */
196     float                coulomb_tab_scale;  /**< table scale/spacing                        */
197     float               *coulomb_tab;        /**< pointer to the table in the device memory  */
198     cudaTextureObject_t  coulomb_tab_texobj; /**< texture object bound to coulomb_tab        */
199 };
200
201 /** \internal
202  * \brief Pair list data.
203  */
204 using cu_plist_t = gpu_plist;
205
206 /** \internal
207  * \brief Typedef of actual timer type.
208  */
209 typedef struct nbnxn_gpu_timers_t cu_timers_t;
210
211 /** \internal
212  * \brief Main data structure for CUDA nonbonded force calculations.
213  */
214 struct gmx_nbnxn_cuda_t
215 {
216     const gmx_device_info_t  *dev_info;       /**< CUDA device information                              */
217     bool                      bUseTwoStreams; /**< true if doing both local/non-local NB work on GPU    */
218     cu_atomdata_t            *atdat;          /**< atom data                                            */
219     cu_nbparam_t             *nbparam;        /**< parameters required for the non-bonded calc.         */
220     cu_plist_t               *plist[2];       /**< pair-list data structures (local and non-local)      */
221     nb_staging_t              nbst;           /**< staging area where fshift/energies get downloaded    */
222
223     cudaStream_t              stream[2];      /**< local and non-local GPU streams                      */
224
225     /** events used for synchronization */
226     cudaEvent_t    nonlocal_done;               /**< event triggered when the non-local non-bonded kernel
227                                                    is done (and the local transfer can proceed)           */
228     cudaEvent_t    misc_ops_and_local_H2D_done; /**< event triggered when the tasks issued in
229                                                    the local stream that need to precede the
230                                                    non-local force calculations are done
231                                                    (e.g. f buffer 0-ing, local x/q H2D) */
232
233     /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
234      * concurrent streams, so we won't time if both l/nl work is done on GPUs.
235      * Timer init/uninit is still done even with timing off so only the condition
236      * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
237     bool                       bDoTime;   /**< True if event-based timing is enabled.               */
238     cu_timers_t               *timers;    /**< CUDA event-based timers.                             */
239     gmx_wallclock_gpu_nbnxn_t *timings;   /**< Timing data. TODO: deprecate this and query timers for accumulated data instead */
240 };
241
242 #ifdef __cplusplus
243 }
244 #endif
245
246 #endif  /* NBNXN_CUDA_TYPES_H */