Remove useless extern CUDA texture reference declarations from PME
[alexxy/gromacs.git] / src / gromacs / ewald / pme.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017, 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
36 /*! \internal \file
37  * \brief This file defines the PME CUDA-specific data structure,
38  * various compile-time constants shared among the PME CUDA kernels,
39  * and also names some PME CUDA memory management routines.
40  * TODO: consider changing defines into variables where possible; have inline getters.
41  *
42  * \author Aleksei Iupinov <a.yupinov@gmail.com>
43  */
44
45 #ifndef GMX_EWALD_PME_CUH
46 #define GMX_EWALD_PME_CUH
47
48 #include <cassert>
49
50 #include <array>
51
52 #include "gromacs/gpu_utils/cuda_arch_utils.cuh" // for warp_size
53
54 #include "pme-gpu-internal.h"                    // for the general PME GPU behaviour defines
55 #include "pme-timings.cuh"
56
57 class GpuParallel3dFft;
58
59 /* Some defines for PME behaviour follow */
60
61 /*
62     Here is a current memory layout for the theta/dtheta B-spline float parameter arrays.
63     This is the data in global memory used both by spreading and gathering kernels (with same scheduling).
64     This example has PME order 4 and 2 particles per warp/data chunk.
65     Each particle has 16 threads assigned to it, each thread works on 4 non-sequential global grid contributions.
66
67     ----------------------------------------------------------------------------
68     particles 0, 1                                        | particles 2, 3     | ...
69     ----------------------------------------------------------------------------
70     order index 0           | index 1 | index 2 | index 3 | order index 0 .....
71     ----------------------------------------------------------------------------
72     tx0 tx1 ty0 ty1 tz0 tz1 | ..........
73     ----------------------------------------------------------------------------
74
75     Each data chunk for a single warp is 24 floats. This goes both for theta and dtheta.
76     24 = 2 particles per warp * order 4 * 3 dimensions. 48 floats (1.5 warp size) per warp in total.
77     I have also tried intertwining theta and theta in a single array (they are used in pairs in gathering stage anyway)
78     and it didn't seem to make a performance difference.
79
80     The corresponding defines follow.
81  */
82
83 /* This is the distance between the neighbour theta elements - would be 2 for the intertwining layout */
84 #define PME_SPLINE_THETA_STRIDE 1
85
86 /*! \brief
87  * The number of GPU threads used for computing spread/gather contributions of a single atom as function of the PME order.
88  * The assumption is currently that any thread processes only a single atom's contributions.
89  */
90 #define PME_SPREADGATHER_THREADS_PER_ATOM (order * order)
91
92 /*! \brief
93  * The number of atoms processed by a single warp in spread/gather.
94  * This macro depends on the templated order parameter (2 atoms per warp for order 4).
95  * It is mostly used for spline data layout tweaked for coalesced access.
96  */
97 #define PME_SPREADGATHER_ATOMS_PER_WARP (warp_size / PME_SPREADGATHER_THREADS_PER_ATOM)
98
99 /*! \brief
100  * Atom data alignment (in terms of number of atoms).
101  * If the GPU atom data buffers are padded (c_usePadding == true),
102  * Then the numbers of atoms which would fit in the padded GPU buffers has to be divisible by this.
103  * The literal number (16) expresses maximum spread/gather block width in warps.
104  * Accordingly, spread and gather block widths in warps should be divisors of this
105  * (e.g. in the pme-spread.cu: constexpr int c_spreadMaxThreadsPerBlock = 8 * warp_size;).
106  * There are debug asserts for this divisibility.
107  */
108 #define PME_ATOM_DATA_ALIGNMENT (16 * PME_SPREADGATHER_ATOMS_PER_WARP)
109
110 /*! \brief \internal
111  * An inline CUDA function for checking the global atom data indices against the atom data array sizes.
112  *
113  * \param[in] atomDataIndexGlobal  The atom data index.
114  * \param[in] nAtomData            The atom data array element count.
115  * \returns                        Non-0 if index is within bounds (or PME data padding is enabled), 0 otherwise.
116  *
117  * This is called from the spline_and_spread and gather PME kernels.
118  * The goal is to isolate the global range checks, and allow avoiding them with c_usePadding enabled.
119  */
120 int __device__ __forceinline__ pme_gpu_check_atom_data_index(const int atomDataIndex, const int nAtomData)
121 {
122     return c_usePadding ? 1 : (atomDataIndex < nAtomData);
123 }
124
125 /*! \brief \internal
126  * An inline CUDA function for skipping the zero-charge atoms.
127  *
128  * \returns                        Non-0 if atom should be processed, 0 otherwise.
129  * \param[in] coefficient          The atom charge.
130  *
131  * This is called from the spline_and_spread and gather PME kernels.
132  */
133 int __device__ __forceinline__ pme_gpu_check_atom_charge(const float coefficient)
134 {
135     assert(isfinite(coefficient));
136     return c_skipNeutralAtoms ? (coefficient != 0.0f) : 1;
137 }
138
139 /*! \brief \internal
140  * The main PME CUDA-specific host data structure, included in the PME GPU structure by the archSpecific pointer.
141  */
142 struct pme_gpu_cuda_t
143 {
144     /*! \brief The CUDA stream where everything related to the PME happens. */
145     cudaStream_t pmeStream;
146
147     /* Synchronization events */
148     /*! \brief Triggered after the energy/virial have been copied to the host (after the solving stage). */
149     cudaEvent_t syncEnerVirD2H;
150     /*! \brief Triggered after the output forces have been copied to the host (after the gathering stage). */
151     cudaEvent_t syncForcesD2H;
152     /*! \brief Triggered after the grid has been copied to the host (after the spreading stage). */
153     cudaEvent_t syncSpreadGridD2H;
154     /*! \brief Triggered after the atom spline data has been copied to the host (after the spline computation). */
155     cudaEvent_t syncSplineAtomDataD2H;
156     /*! \brief Triggered after the grid hes been copied to the host (after the solving stage) */
157     cudaEvent_t syncSolveGridD2H;
158
159     // TODO: consider moving some things below into the non-CUDA struct.
160
161     /* Settings which are set at the start of the run */
162     /*! \brief A boolean which tells whether the complex and real grids for cuFFT are different or same. Currenty true. */
163     bool performOutOfPlaceFFT;
164     /*! \brief A boolean which tells if the CUDA timing events are enabled.
165      * True by default, disabled by setting the environment variable GMX_DISABLE_CUDA_TIMING.
166      * FIXME: this should also be disabled if any other GPU task is running concurrently on the same device,
167      * as CUDA events on multiple streams are untrustworthy.
168      */
169     bool                                             useTiming;
170
171     std::vector<std::unique_ptr<GpuParallel3dFft > > fftSetup;
172
173     std::array<GpuRegionTimer, gtPME_EVENT_COUNT>    timingEvents;
174
175     /* GPU arrays element counts (not the arrays sizes in bytes!).
176      * They might be larger than the actual meaningful data sizes.
177      * These are paired: the actual element count + the maximum element count that can fit in the current allocated memory.
178      * These integer pairs are mostly meaningful for the cu_realloc/free_buffered calls.
179      * As such, if cu_realloc/free_buffered is refactored, they can be freely changed, too.
180      * The only exceptions are realGridSize and complexGridSize which are also used for grid clearing/copying.
181      * TODO: these should live in a clean buffered container type, and be refactored in the NB/cudautils as well.
182      */
183     /*! \brief The kernelParams.atoms.coordinates float element count (actual)*/
184     int coordinatesSize;
185     /*! \brief The kernelParams.atoms.coordinates float element count (reserved) */
186     int coordinatesSizeAlloc;
187     /*! \brief The kernelParams.atoms.forces float element count (actual) */
188     int forcesSize;
189     /*! \brief The kernelParams.atoms.forces float element count (reserved) */
190     int forcesSizeAlloc;
191     /*! \brief The kernelParams.atoms.gridlineIndices int element count (actual) */
192     int gridlineIndicesSize;
193     /*! \brief The kernelParams.atoms.gridlineIndices int element count (reserved) */
194     int gridlineIndicesSizeAlloc;
195     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (actual) */
196     int splineDataSize;
197     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (reserved) */
198     int splineDataSizeAlloc;
199     /*! \brief The kernelParams.atoms.coefficients float element count (actual) */
200     int coefficientsSize;
201     /*! \brief The kernelParams.atoms.coefficients float element count (reserved) */
202     int coefficientsSizeAlloc;
203     /*! \brief The kernelParams.grid.splineValuesArray float element count (actual) */
204     int splineValuesSize;
205     /*! \brief The kernelParams.grid.splineValuesArray float element count (reserved) */
206     int splineValuesSizeAlloc;
207     /*! \brief The kernelParams.grid.realGrid float element count (actual) */
208     int realGridSize;
209     /*! \brief The kernelParams.grid.realGrid float element count (reserved) */
210     int realGridSizeAlloc;
211     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (actual) */
212     int complexGridSize;
213     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (reserved) */
214     int complexGridSizeAlloc;
215 };
216
217
218 /*! \brief \internal
219  * A single structure encompassing all the PME data used in CUDA kernels.
220  * This inherits from pme_gpu_kernel_params_base_t and adds a couple cudaTextureObject_t handles,
221  * which we would like to avoid in plain C++.
222  */
223 struct pme_gpu_cuda_kernel_params_t : pme_gpu_kernel_params_base_t
224 {
225     /* These are CUDA texture objects, related to the grid size. */
226     /*! \brief CUDA texture object for accessing grid.d_fractShiftsTable */
227     cudaTextureObject_t fractShiftsTableTexture;
228     /*! \brief CUDA texture object for accessing grid.d_gridlineIndicesTable */
229     cudaTextureObject_t gridlineIndicesTableTexture;
230 };
231
232 /* CUDA texture reference functions which reside in respective kernel files
233  * (due to texture references having scope of a translation unit).
234  */
235 /*! Returns the reference to the gridlineIndices texture. */
236 texture<int, 1, cudaReadModeElementType>   &pme_gpu_get_gridline_texref();
237 /*! Returns the reference to the fractShifts texture. */
238 texture<float, 1, cudaReadModeElementType> &pme_gpu_get_fract_shifts_texref();
239
240 #endif