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