SYCL NBNXM offload support
[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 /** \internal
67  * \brief Parameters required for the GPU nonbonded calculations.
68  */
69 struct NBParamGpu
70 {
71
72     //! type of electrostatics
73     enum Nbnxm::ElecType elecType;
74     //! type of VdW impl.
75     enum Nbnxm::VdwType vdwType;
76
77     //! charge multiplication factor
78     float epsfac;
79     //! Reaction-field/plain cutoff electrostatics const.
80     float c_rf;
81     //! Reaction-field electrostatics constant
82     float two_k_rf;
83     //! Ewald/PME parameter
84     float ewald_beta;
85     //! Ewald/PME correction term subtracted from the direct-space potential
86     float sh_ewald;
87     //! LJ-Ewald/PME correction term added to the correction potential
88     float sh_lj_ewald;
89     //! LJ-Ewald/PME coefficient
90     float ewaldcoeff_lj;
91
92     //! Coulomb cut-off squared
93     float rcoulomb_sq;
94
95     //! VdW cut-off squared
96     float rvdw_sq;
97     //! VdW switched cut-off
98     float rvdw_switch;
99     //! Full, outer pair-list cut-off squared
100     float rlistOuter_sq;
101     //! Inner, dynamic pruned pair-list cut-off squared
102     float rlistInner_sq;
103     //! True if we use dynamic pair-list pruning
104     bool useDynamicPruning;
105
106     //! VdW shift dispersion constants
107     shift_consts_t dispersion_shift;
108     //! VdW shift repulsion constants
109     shift_consts_t repulsion_shift;
110     //! VdW switch constants
111     switch_consts_t vdw_switch;
112
113     /* LJ non-bonded parameters - accessed through texture memory */
114     //! nonbonded parameter table with C6/C12 pairs per atom type-pair, 2*ntype^2 elements
115     DeviceBuffer<float> nbfp;
116     //! texture object bound to nbfp
117     DeviceTexture nbfp_texobj;
118     //! nonbonded parameter table per atom type, 2*ntype elements
119     DeviceBuffer<float> nbfp_comb;
120     //! texture object bound to nbfp_comb
121     DeviceTexture nbfp_comb_texobj;
122
123     /* Ewald Coulomb force table data - accessed through texture memory */
124     //! table scale/spacing
125     float coulomb_tab_scale;
126     //! pointer to the table in the device memory
127     DeviceBuffer<float> coulomb_tab;
128     //! texture object bound to coulomb_tab
129     DeviceTexture coulomb_tab_texobj;
130 };
131
132 namespace Nbnxm
133 {
134
135 using gmx::AtomLocality;
136 using gmx::InteractionLocality;
137
138 /*! \internal
139  * \brief GPU region timers used for timing GPU kernels and H2D/D2H transfers.
140  *
141  * The two-sized arrays hold the local and non-local values and should always
142  * be indexed with eintLocal/eintNonlocal.
143  */
144 struct gpu_timers_t
145 {
146     /*! \internal
147      * \brief Timers for local or non-local coordinate/force transfers
148      */
149     struct XFTransfers
150     {
151         //! timer for x/q H2D transfers (l/nl, every step)
152         GpuRegionTimer nb_h2d;
153         //! timer for f D2H transfer (l/nl, every step)
154         GpuRegionTimer nb_d2h;
155     };
156
157     /*! \internal
158      * \brief Timers for local or non-local interaction related operations
159      */
160     struct Interaction
161     {
162         //! timer for pair-list H2D transfers (l/nl, every PS step)
163         GpuRegionTimer pl_h2d;
164         //! true when a pair-list transfer has been done at this step
165         bool didPairlistH2D = false;
166         //! timer for non-bonded kernels (l/nl, every step)
167         GpuRegionTimer nb_k;
168         //! timer for the 1st pass list pruning kernel (l/nl, every PS step)
169         GpuRegionTimer prune_k;
170         //! true when we timed pruning and the timings need to be accounted for
171         bool didPrune = false;
172         //! timer for rolling pruning kernels (l/nl, frequency depends on chunk size)
173         GpuRegionTimer rollingPrune_k;
174         //! true when we timed rolling pruning (at the previous step) and the timings need to be accounted for
175         bool didRollingPrune = false;
176     };
177
178     //! timer for atom data transfer (every PS step)
179     GpuRegionTimer atdat;
180     //! timers for coordinate/force transfers (every step)
181     gmx::EnumerationArray<AtomLocality, XFTransfers> xf;
182     //! timers for interaction related transfers
183     gmx::EnumerationArray<InteractionLocality, Nbnxm::gpu_timers_t::Interaction> interaction;
184 };
185
186 /*! \internal
187  * \brief GPU pair list structure */
188 struct gpu_plist
189 {
190     //! number of atoms per cluster
191     int na_c;
192
193     //! size of sci, # of i clusters in the list
194     int nsci;
195     //! allocation size of sci
196     int sci_nalloc;
197     //! list of i-cluster ("super-clusters")
198     DeviceBuffer<nbnxn_sci_t> sci;
199
200     //! total # of 4*j clusters
201     int ncj4;
202     //! allocation size of cj4
203     int cj4_nalloc;
204     //! 4*j cluster list, contains j cluster number and index into the i cluster list
205     DeviceBuffer<nbnxn_cj4_t> cj4;
206     //! # of 4*j clusters * # of warps
207     int nimask;
208     //! allocation size of imask
209     int imask_nalloc;
210     //! imask for 2 warps for each 4*j cluster group
211     DeviceBuffer<unsigned int> imask;
212     //! atom interaction bits
213     DeviceBuffer<nbnxn_excl_t> excl;
214     //! count for excl
215     int nexcl;
216     //! allocation size of excl
217     int excl_nalloc;
218
219     /* parameter+variables for normal and rolling pruning */
220     //! true after search, indicates that initial pruning with outer pruning is needed
221     bool haveFreshList;
222     //! the number of parts/steps over which one cycle of rolling pruning takes places
223     int rollingPruningNumParts;
224     //! the next part to which the rolling pruning needs to be applied
225     int rollingPruningPart;
226 };
227
228 } // namespace Nbnxm
229
230 #endif