28dbe70459e4728c67b17331da1fb29e8e0b245f
[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, 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/locality.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49
50 #include "pairlist.h"
51
52 #if GMX_GPU == GMX_GPU_OPENCL
53 #    include "gromacs/gpu_utils/gpuregiontimer_ocl.h"
54 #endif
55
56 #if GMX_GPU == GMX_GPU_CUDA
57 #    include "gromacs/gpu_utils/gpuregiontimer.cuh"
58 #endif
59
60 namespace Nbnxm
61 {
62
63 using gmx::AtomLocality;
64 using gmx::InteractionLocality;
65
66 /*! \internal
67  * \brief GPU region timers used for timing GPU kernels and H2D/D2H transfers.
68  *
69  * The two-sized arrays hold the local and non-local values and should always
70  * be indexed with eintLocal/eintNonlocal.
71  */
72 struct gpu_timers_t
73 {
74     /*! \internal
75      * \brief Timers for local or non-local coordinate/force transfers
76      */
77     struct XFTransfers
78     {
79         //! timer for x/q H2D transfers (l/nl, every step)
80         GpuRegionTimer nb_h2d;
81         //! timer for f D2H transfer (l/nl, every step)
82         GpuRegionTimer nb_d2h;
83     };
84
85     /*! \internal
86      * \brief Timers for local or non-local interaction related operations
87      */
88     struct Interaction
89     {
90         //! timer for pair-list H2D transfers (l/nl, every PS step)
91         GpuRegionTimer pl_h2d;
92         //! true when a pair-list transfer has been done at this step
93         bool didPairlistH2D = false;
94         //! timer for non-bonded kernels (l/nl, every step)
95         GpuRegionTimer nb_k;
96         //! timer for the 1st pass list pruning kernel (l/nl, every PS step)
97         GpuRegionTimer prune_k;
98         //! true when we timed pruning and the timings need to be accounted for
99         bool didPrune = false;
100         //! timer for rolling pruning kernels (l/nl, frequency depends on chunk size)
101         GpuRegionTimer rollingPrune_k;
102         //! true when we timed rolling pruning (at the previous step) and the timings need to be accounted for
103         bool didRollingPrune = false;
104     };
105
106     //! timer for atom data transfer (every PS step)
107     GpuRegionTimer atdat;
108     //! timers for coordinate/force transfers (every step)
109     gmx::EnumerationArray<AtomLocality, XFTransfers> xf;
110     //! timers for interaction related transfers
111     gmx::EnumerationArray<InteractionLocality, Nbnxm::gpu_timers_t::Interaction> interaction;
112 };
113
114 /*! \internal
115  * \brief GPU pair list structure */
116 struct gpu_plist
117 {
118     //! number of atoms per cluster
119     int na_c;
120
121     //! size of sci, # of i clusters in the list
122     int nsci;
123     //! allocation size of sci
124     int sci_nalloc;
125     //! list of i-cluster ("super-clusters")
126     DeviceBuffer<nbnxn_sci_t> sci;
127
128     //! total # of 4*j clusters
129     int ncj4;
130     //! allocation size of cj4
131     int cj4_nalloc;
132     //! 4*j cluster list, contains j cluster number and index into the i cluster list
133     DeviceBuffer<nbnxn_cj4_t> cj4;
134     //! # of 4*j clusters * # of warps
135     int nimask;
136     //! allocation size of imask
137     int imask_nalloc;
138     //! imask for 2 warps for each 4*j cluster group
139     DeviceBuffer<unsigned int> imask;
140     //! atom interaction bits
141     DeviceBuffer<nbnxn_excl_t> excl;
142     //! count for excl
143     int nexcl;
144     //! allocation size of excl
145     int excl_nalloc;
146
147     /* parameter+variables for normal and rolling pruning */
148     //! true after search, indictes that initial pruning with outer prunning is needed
149     bool haveFreshList;
150     //! the number of parts/steps over which one cyle of roling pruning takes places
151     int rollingPruningNumParts;
152     //! the next part to which the roling pruning needs to be applied
153     int rollingPruningPart;
154 };
155
156 } // namespace Nbnxm
157
158 #endif