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