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