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