Make nbnxm headers more self-contained
[alexxy/gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_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) 2013-2019,2020, 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_nbnxm
44  */
45
46 #ifndef NBNXM_CUDA_TYPES_H
47 #define NBNXM_CUDA_TYPES_H
48
49 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/gpu_utils/devicebuffer.h"
52 #include "gromacs/gpu_utils/gputraits.cuh"
53 #include "gromacs/mdtypes/interaction_const.h"
54 #include "gromacs/nbnxm/gpu_types_common.h"
55 #include "gromacs/nbnxm/nbnxm.h"
56 #include "gromacs/nbnxm/pairlist.h"
57 #include "gromacs/timing/gpu_timing.h"
58 #include "gromacs/utility/enumerationhelpers.h"
59
60 /*! \brief Macro definining default for the prune kernel's j4 processing concurrency.
61  *
62  *  The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override.
63  */
64 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
65 #    define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY 4
66 #endif
67 /*! \brief Default for the prune kernel's j4 processing concurrency.
68  *
69  *  Initialized using the #GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro which allows compile-time override.
70  */
71 const int c_cudaPruneKernelJ4Concurrency = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY;
72
73 /* TODO: consider moving this to kernel_utils */
74 /* Convenience defines */
75 /*! \brief number of clusters per supercluster. */
76 static const int c_numClPerSupercl = c_nbnxnGpuNumClusterPerSupercluster;
77 /*! \brief cluster size = number of atoms per cluster. */
78 static const int c_clSize = c_nbnxnGpuClusterSize;
79
80 /*! \brief Electrostatic CUDA kernel flavors.
81  *
82  *  Types of electrostatics implementations available in the CUDA non-bonded
83  *  force kernels. These represent both the electrostatics types implemented
84  *  by the kernels (cut-off, RF, and Ewald - a subset of what's defined in
85  *  enums.h) as well as encode implementation details analytical/tabulated
86  *  and single or twin cut-off (for Ewald kernels).
87  *  Note that the cut-off and RF kernels have only analytical flavor and unlike
88  *  in the CPU kernels, the tabulated kernels are ATM Ewald-only.
89  *
90  *  The row-order of pointers to different electrostatic kernels defined in
91  *  nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
92  *  should match the order of enumerated types below.
93  */
94 enum eelCu
95 {
96     eelCuCUT,
97     eelCuRF,
98     eelCuEWALD_TAB,
99     eelCuEWALD_TAB_TWIN,
100     eelCuEWALD_ANA,
101     eelCuEWALD_ANA_TWIN,
102     eelCuNR
103 };
104
105 /*! \brief VdW CUDA kernel flavors.
106  *
107  * The enumerates values correspond to the LJ implementations in the CUDA non-bonded
108  * kernels.
109  *
110  * The column-order of pointers to different electrostatic kernels defined in
111  * nbnxn_cuda.cu by the nb_*_kfunc_ptr function pointer table
112  * should match the order of enumerated types below.
113  */
114 enum evdwCu
115 {
116     evdwCuCUT,
117     evdwCuCUTCOMBGEOM,
118     evdwCuCUTCOMBLB,
119     evdwCuFSWITCH,
120     evdwCuPSWITCH,
121     evdwCuEWALDGEOM,
122     evdwCuEWALDLB,
123     evdwCuNR
124 };
125
126 /* All structs prefixed with "cu_" hold data used in GPU calculations and
127  * are passed to the kernels, except cu_timers_t. */
128 /*! \cond */
129 typedef struct cu_atomdata cu_atomdata_t;
130 typedef struct cu_nbparam  cu_nbparam_t;
131 typedef struct nb_staging  nb_staging_t;
132 /*! \endcond */
133
134
135 /** \internal
136  * \brief Staging area for temporary data downloaded from the GPU.
137  *
138  *  The energies/shift forces get downloaded here first, before getting added
139  *  to the CPU-side aggregate values.
140  */
141 struct nb_staging
142 {
143     float*  e_lj;   /**< LJ energy            */
144     float*  e_el;   /**< electrostatic energy */
145     float3* fshift; /**< shift forces         */
146 };
147
148 /** \internal
149  * \brief Nonbonded atom data - both inputs and outputs.
150  */
151 struct cu_atomdata
152 {
153     int natoms;       /**< number of atoms                              */
154     int natoms_local; /**< number of local atoms                        */
155     int nalloc;       /**< allocation size for the atom data (xq, f)    */
156
157     float4* xq; /**< atom coordinates + charges, size natoms      */
158     float3* f;  /**< force output array, size natoms              */
159
160     float* e_lj; /**< LJ energy output, size 1                     */
161     float* e_el; /**< Electrostatics energy input, size 1          */
162
163     float3* fshift; /**< shift forces                                 */
164
165     int     ntypes;     /**< number of atom types                         */
166     int*    atom_types; /**< atom type indices, size natoms               */
167     float2* lj_comb;    /**< sqrt(c6),sqrt(c12) size natoms               */
168
169     float3* shift_vec;         /**< shifts                                       */
170     bool    bShiftVecUploaded; /**< true if the shift vector has been uploaded   */
171 };
172
173 /** \internal
174  * \brief Parameters required for the CUDA nonbonded calculations.
175  */
176 struct cu_nbparam
177 {
178
179     int eeltype; /**< type of electrostatics, takes values from #eelCu */
180     int vdwtype; /**< type of VdW impl., takes values from #evdwCu     */
181
182     float epsfac;      /**< charge multiplication factor                      */
183     float c_rf;        /**< Reaction-field/plain cutoff electrostatics const. */
184     float two_k_rf;    /**< Reaction-field electrostatics constant            */
185     float ewald_beta;  /**< Ewald/PME parameter                               */
186     float sh_ewald;    /**< Ewald/PME correction term substracted from the direct-space potential */
187     float sh_lj_ewald; /**< LJ-Ewald/PME correction term added to the correction potential        */
188     float ewaldcoeff_lj; /**< LJ-Ewald/PME coefficient                          */
189
190     float rcoulomb_sq; /**< Coulomb cut-off squared                           */
191
192     float rvdw_sq;           /**< VdW cut-off squared                               */
193     float rvdw_switch;       /**< VdW switched cut-off                              */
194     float rlistOuter_sq;     /**< Full, outer pair-list cut-off squared             */
195     float rlistInner_sq;     /**< Inner, dynamic pruned pair-list cut-off squared   */
196     bool  useDynamicPruning; /**< True if we use dynamic pair-list pruning          */
197
198     shift_consts_t  dispersion_shift; /**< VdW shift dispersion constants           */
199     shift_consts_t  repulsion_shift;  /**< VdW shift repulsion constants            */
200     switch_consts_t vdw_switch;       /**< VdW switch constants                     */
201
202     /* LJ non-bonded parameters - accessed through texture memory */
203     float*              nbfp; /**< nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements */
204     cudaTextureObject_t nbfp_texobj; /**< texture object bound to nbfp */
205     float*              nbfp_comb; /**< nonbonded parameter table per atom type, 2*ntype elements */
206     cudaTextureObject_t nbfp_comb_texobj; /**< texture object bound to nbfp_texobj */
207
208     /* Ewald Coulomb force table data - accessed through texture memory */
209     float               coulomb_tab_scale;  /**< table scale/spacing                        */
210     float*              coulomb_tab;        /**< pointer to the table in the device memory  */
211     cudaTextureObject_t coulomb_tab_texobj; /**< texture object bound to coulomb_tab        */
212 };
213
214 /** \internal
215  * \brief Pair list data.
216  */
217 using cu_plist_t = Nbnxm::gpu_plist;
218
219 /** \internal
220  * \brief Typedef of actual timer type.
221  */
222 typedef struct Nbnxm::gpu_timers_t cu_timers_t;
223
224 class GpuEventSynchronizer;
225
226 /*! \internal
227  * \brief Main data structure for CUDA nonbonded force calculations.
228  */
229 struct gmx_nbnxm_gpu_t
230 {
231     /*! \brief CUDA device information */
232     const gmx_device_info_t* dev_info;
233     /*! \brief true if doing both local/non-local NB work on GPU */
234     bool bUseTwoStreams;
235     /*! \brief atom data */
236     cu_atomdata_t* atdat;
237     /*! \brief f buf ops cell index mapping */
238     int* cell;
239     /*! \brief number of indices in cell buffer */
240     int ncell;
241     /*! \brief number of indices allocated in cell buffer */
242     int ncell_alloc;
243     /*! \brief array of atom indices */
244     int* atomIndices;
245     /*! \brief size of atom indices */
246     int atomIndicesSize;
247     /*! \brief size of atom indices allocated in device buffer */
248     int atomIndicesSize_alloc;
249     /*! \brief x buf ops num of atoms */
250     int* cxy_na;
251     /*! \brief number of elements in cxy_na */
252     int ncxy_na;
253     /*! \brief number of elements allocated allocated in device buffer */
254     int ncxy_na_alloc;
255     /*! \brief x buf ops cell index mapping */
256     int* cxy_ind;
257     /*! \brief number of elements in cxy_ind */
258     int ncxy_ind;
259     /*! \brief number of elements allocated allocated in device buffer */
260     int ncxy_ind_alloc;
261     /*! \brief parameters required for the non-bonded calc. */
262     cu_nbparam_t* nbparam;
263     /*! \brief pair-list data structures (local and non-local) */
264     gmx::EnumerationArray<Nbnxm::InteractionLocality, cu_plist_t*> plist;
265     /*! \brief staging area where fshift/energies get downloaded */
266     nb_staging_t nbst;
267     /*! \brief local and non-local GPU streams */
268     gmx::EnumerationArray<Nbnxm::InteractionLocality, cudaStream_t> stream;
269
270     /*! \brief Events used for synchronization */
271     /*! \{ */
272     /*! \brief Event triggered when the non-local non-bonded
273      * kernel is done (and the local transfer can proceed) */
274     cudaEvent_t nonlocal_done;
275     /*! \brief Event triggered when the tasks issued in the local
276      * stream that need to precede the non-local force or buffer
277      * operation calculations are done (e.g. f buffer 0-ing, local
278      * x/q H2D, buffer op initialization in local stream that is
279      * required also by nonlocal stream ) */
280     cudaEvent_t misc_ops_and_local_H2D_done;
281     /*! \} */
282
283     /*! \brief True if there is work for the current domain in the
284      * respective locality.
285      *
286      * This includes local/nonlocal GPU work, either bonded or
287      * nonbonded, scheduled to be executed in the current
288      * domain. As long as bonded work is not split up into
289      * local/nonlocal, if there is bonded GPU work, both flags
290      * will be true. */
291     gmx::EnumerationArray<Nbnxm::InteractionLocality, bool> haveWork;
292
293     /*! \brief Pointer to event synchronizer triggered when the local
294      * GPU buffer ops / reduction is complete
295      *
296      * \note That the synchronizer is managed outside of this module
297      * in StatePropagatorDataGpu.
298      */
299     GpuEventSynchronizer* localFReductionDone;
300
301     /*! \brief Event triggered when non-local coordinate buffer
302      * has been copied from device to host. */
303     GpuEventSynchronizer* xNonLocalCopyD2HDone;
304
305     /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
306      * concurrent streams, so we won't time if both l/nl work is done on GPUs.
307      * Timer init/uninit is still done even with timing off so only the condition
308      * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
309     /*! \brief True if event-based timing is enabled. */
310     bool bDoTime;
311     /*! \brief CUDA event-based timers. */
312     cu_timers_t* timers;
313     /*! \brief Timing data. TODO: deprecate this and query timers for accumulated data instead */
314     gmx_wallclock_gpu_nbnxn_t* timings;
315 };
316
317 #endif /* NBNXN_CUDA_TYPES_H */