Use GpuTimers in CUDA and OpenCL versions of NBNXM directly
[alexxy/gromacs.git] / src / gromacs / nbnxm / gpu_types_common.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief Implements common internal types for different NBNXN GPU implementations
37  *
38  * \author Szilárd Páll <pall.szilard@gmail.com>
39  * \ingroup module_nbnxm
40  */
41
42 #ifndef GMX_MDLIB_NBNXN_GPU_COMMON_TYPES_H
43 #define GMX_MDLIB_NBNXN_GPU_COMMON_TYPES_H
44
45 #include "config.h"
46
47 #include "gromacs/mdtypes/interaction_const.h"
48 #include "gromacs/mdtypes/locality.h"
49 #include "gromacs/utility/enumerationhelpers.h"
50
51 #include "nbnxm.h"
52 #include "pairlist.h"
53
54 #if GMX_GPU_OPENCL
55 #    include "gromacs/gpu_utils/gpuregiontimer_ocl.h"
56 #endif
57
58 #if GMX_GPU_CUDA
59 #    include "gromacs/gpu_utils/gpuregiontimer.cuh"
60 #endif
61
62 #if GMX_GPU_SYCL
63 #    include "gromacs/gpu_utils/gpuregiontimer_sycl.h"
64 #endif
65
66 /*! \brief Macro definining default for the prune kernel's j4 processing concurrency.
67  *
68  *  The GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY macro allows compile-time override with the default value of 4.
69  */
70 #ifndef GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY
71 #    define GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY 4
72 #endif
73 //! Default for the prune kernel's j4 processing concurrency.
74 static constexpr int c_pruneKernelJ4Concurrency = GMX_NBNXN_PRUNE_KERNEL_J4_CONCURRENCY;
75
76 /*! \internal
77  * \brief Staging area for temporary data downloaded from the GPU.
78  *
79  * Since SYCL buffers already have host-side storage, this is a bit redundant.
80  * But it allows prefetching of the data from GPU, and brings GPU backends closer together.
81  */
82 struct NBStagingData
83 {
84     //! LJ energy
85     float* eLJ = nullptr;
86     //! electrostatic energy
87     float* eElec = nullptr;
88     //! shift forces
89     Float3* fShift = nullptr;
90 };
91
92 /** \internal
93  * \brief Nonbonded atom data - both inputs and outputs.
94  */
95 struct NBAtomData
96 {
97     //! number of atoms
98     int numAtoms;
99     //! number of local atoms
100     int numAtomsLocal;
101     //! allocation size for the atom data (xq, f)
102     int numAtomsAlloc;
103
104     //! atom coordinates + charges, size \ref numAtoms
105     DeviceBuffer<Float4> xq;
106     //! force output array, size \ref numAtoms
107     DeviceBuffer<Float3> f;
108
109     //! LJ energy output, size 1
110     DeviceBuffer<float> eLJ;
111     //! Electrostatics energy input, size 1
112     DeviceBuffer<float> eElec;
113
114     //! shift forces
115     DeviceBuffer<Float3> fShift;
116
117     //! number of atom types
118     int numTypes;
119     //! atom type indices, size \ref numAtoms
120     DeviceBuffer<int> atomTypes;
121     //! sqrt(c6),sqrt(c12) size \ref numAtoms
122     DeviceBuffer<Float2> ljComb;
123
124     //! shifts
125     DeviceBuffer<Float3> shiftVec;
126     //! true if the shift vector has been uploaded
127     bool shiftVecUploaded;
128 };
129
130 /** \internal
131  * \brief Parameters required for the GPU nonbonded calculations.
132  */
133 struct NBParamGpu
134 {
135
136     //! type of electrostatics
137     enum Nbnxm::ElecType elecType;
138     //! type of VdW impl.
139     enum Nbnxm::VdwType vdwType;
140
141     //! charge multiplication factor
142     float epsfac;
143     //! Reaction-field/plain cutoff electrostatics const.
144     float c_rf;
145     //! Reaction-field electrostatics constant
146     float two_k_rf;
147     //! Ewald/PME parameter
148     float ewald_beta;
149     //! Ewald/PME correction term subtracted from the direct-space potential
150     float sh_ewald;
151     //! LJ-Ewald/PME correction term added to the correction potential
152     float sh_lj_ewald;
153     //! LJ-Ewald/PME coefficient
154     float ewaldcoeff_lj;
155
156     //! Coulomb cut-off squared
157     float rcoulomb_sq;
158
159     //! VdW cut-off squared
160     float rvdw_sq;
161     //! VdW switched cut-off
162     float rvdw_switch;
163     //! Full, outer pair-list cut-off squared
164     float rlistOuter_sq;
165     //! Inner, dynamic pruned pair-list cut-off squared
166     float rlistInner_sq;
167     //! True if we use dynamic pair-list pruning
168     bool useDynamicPruning;
169
170     //! VdW shift dispersion constants
171     shift_consts_t dispersion_shift;
172     //! VdW shift repulsion constants
173     shift_consts_t repulsion_shift;
174     //! VdW switch constants
175     switch_consts_t vdw_switch;
176
177     /* LJ non-bonded parameters - accessed through texture memory */
178     //! nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements
179     DeviceBuffer<float> nbfp;
180     //! texture object bound to nbfp
181     DeviceTexture nbfp_texobj;
182     //! nonbonded parameter table per atom type, 2*ntype elements
183     DeviceBuffer<float> nbfp_comb;
184     //! texture object bound to nbfp_comb
185     DeviceTexture nbfp_comb_texobj;
186
187     /* Ewald Coulomb force table data - accessed through texture memory */
188     //! table scale/spacing
189     float coulomb_tab_scale;
190     //! pointer to the table in the device memory
191     DeviceBuffer<float> coulomb_tab;
192     //! texture object bound to coulomb_tab
193     DeviceTexture coulomb_tab_texobj;
194 };
195
196 namespace Nbnxm
197 {
198
199 using gmx::AtomLocality;
200 using gmx::InteractionLocality;
201
202 /*! \internal
203  * \brief GPU region timers used for timing GPU kernels and H2D/D2H transfers.
204  *
205  * The two-sized arrays hold the local and non-local values and should always
206  * be indexed with eintLocal/eintNonlocal.
207  */
208 struct GpuTimers
209 {
210     /*! \internal
211      * \brief Timers for local or non-local coordinate/force transfers
212      */
213     struct XFTransfers
214     {
215         //! timer for x/q H2D transfers (l/nl, every step)
216         GpuRegionTimer nb_h2d;
217         //! timer for f D2H transfer (l/nl, every step)
218         GpuRegionTimer nb_d2h;
219     };
220
221     /*! \internal
222      * \brief Timers for local or non-local interaction related operations
223      */
224     struct Interaction
225     {
226         //! timer for pair-list H2D transfers (l/nl, every PS step)
227         GpuRegionTimer pl_h2d;
228         //! true when a pair-list transfer has been done at this step
229         bool didPairlistH2D = false;
230         //! timer for non-bonded kernels (l/nl, every step)
231         GpuRegionTimer nb_k;
232         //! timer for the 1st pass list pruning kernel (l/nl, every PS step)
233         GpuRegionTimer prune_k;
234         //! true when we timed pruning and the timings need to be accounted for
235         bool didPrune = false;
236         //! timer for rolling pruning kernels (l/nl, frequency depends on chunk size)
237         GpuRegionTimer rollingPrune_k;
238         //! true when we timed rolling pruning (at the previous step) and the timings need to be accounted for
239         bool didRollingPrune = false;
240     };
241
242     //! timer for atom data transfer (every PS step)
243     GpuRegionTimer atdat;
244     //! timers for coordinate/force transfers (every step)
245     gmx::EnumerationArray<AtomLocality, XFTransfers> xf;
246     //! timers for interaction related transfers
247     gmx::EnumerationArray<InteractionLocality, Nbnxm::GpuTimers::Interaction> interaction;
248 };
249
250 /*! \internal
251  * \brief GPU pair list structure */
252 struct gpu_plist
253 {
254     //! number of atoms per cluster
255     int na_c;
256
257     //! size of sci, # of i clusters in the list
258     int nsci;
259     //! allocation size of sci
260     int sci_nalloc;
261     //! list of i-cluster ("super-clusters")
262     DeviceBuffer<nbnxn_sci_t> sci;
263
264     //! total # of 4*j clusters
265     int ncj4;
266     //! allocation size of cj4
267     int cj4_nalloc;
268     //! 4*j cluster list, contains j cluster number and index into the i cluster list
269     DeviceBuffer<nbnxn_cj4_t> cj4;
270     //! # of 4*j clusters * # of warps
271     int nimask;
272     //! allocation size of imask
273     int imask_nalloc;
274     //! imask for 2 warps for each 4*j cluster group
275     DeviceBuffer<unsigned int> imask;
276     //! atom interaction bits
277     DeviceBuffer<nbnxn_excl_t> excl;
278     //! count for excl
279     int nexcl;
280     //! allocation size of excl
281     int excl_nalloc;
282
283     /* parameter+variables for normal and rolling pruning */
284     //! true after search, indicates that initial pruning with outer pruning is needed
285     bool haveFreshList;
286     //! the number of parts/steps over which one cycle of rolling pruning takes places
287     int rollingPruningNumParts;
288     //! the next part to which the rolling pruning needs to be applied
289     int rollingPruningPart;
290 };
291
292 } // namespace Nbnxm
293
294 #endif