Separate PME GPU host-only and host/device data structures
[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,2018, 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 #include <set>
52
53 #include "gromacs/gpu_utils/cuda_arch_utils.cuh" // for warp_size
54
55 #include "pme-gpu-internal.h"                    // for the general PME GPU behaviour defines
56 #include "pme-gpu-types.h"
57 #include "pme-gpu-types-host.h"
58 #include "pme-timings.cuh"
59
60 class GpuParallel3dFft;
61
62 /* Some defines for PME behaviour follow */
63
64 /*
65     Here is a current memory layout for the theta/dtheta B-spline float parameter arrays.
66     This is the data in global memory used both by spreading and gathering kernels (with same scheduling).
67     This example has PME order 4 and 2 particles per warp/data chunk.
68     Each particle has 16 threads assigned to it, each thread works on 4 non-sequential global grid contributions.
69
70     ----------------------------------------------------------------------------
71     particles 0, 1                                        | particles 2, 3     | ...
72     ----------------------------------------------------------------------------
73     order index 0           | index 1 | index 2 | index 3 | order index 0 .....
74     ----------------------------------------------------------------------------
75     tx0 tx1 ty0 ty1 tz0 tz1 | ..........
76     ----------------------------------------------------------------------------
77
78     Each data chunk for a single warp is 24 floats. This goes both for theta and dtheta.
79     24 = 2 particles per warp * order 4 * 3 dimensions. 48 floats (1.5 warp size) per warp in total.
80     I have also tried intertwining theta and theta in a single array (they are used in pairs in gathering stage anyway)
81     and it didn't seem to make a performance difference.
82
83     The spline indexing is isolated in the 2 inline functions below:
84     getSplineParamIndexBase() return a base shared memory index corresponding to the atom in the block;
85     getSplineParamIndex() consumes its results and adds offsets for dimension and spline value index.
86
87     The corresponding defines follow.
88  */
89
90 /*! \brief
91  * The number of GPU threads used for computing spread/gather contributions of a single atom as function of the PME order.
92  * The assumption is currently that any thread processes only a single atom's contributions.
93  */
94 #define PME_SPREADGATHER_THREADS_PER_ATOM (order * order)
95
96 /*! \brief
97  * The number of atoms processed by a single warp in spread/gather.
98  * This macro depends on the templated order parameter (2 atoms per warp for order 4).
99  * It is mostly used for spline data layout tweaked for coalesced access.
100  */
101 #define PME_SPREADGATHER_ATOMS_PER_WARP (warp_size / PME_SPREADGATHER_THREADS_PER_ATOM)
102
103 /*! \brief
104  * Atom data alignment (in terms of number of atoms).
105  * If the GPU atom data buffers are padded (c_usePadding == true),
106  * Then the numbers of atoms which would fit in the padded GPU buffers has to be divisible by this.
107  * The literal number (16) expresses maximum spread/gather block width in warps.
108  * Accordingly, spread and gather block widths in warps should be divisors of this
109  * (e.g. in the pme-spread.cu: constexpr int c_spreadMaxThreadsPerBlock = 8 * warp_size;).
110  * There are debug asserts for this divisibility.
111  */
112 #define PME_ATOM_DATA_ALIGNMENT (16 * PME_SPREADGATHER_ATOMS_PER_WARP)
113
114 /*! \internal \brief
115  * Gets a base of the unique index to an element in a spline parameter buffer (theta/dtheta),
116  * which is laid out for GPU spread/gather kernels. The base only corresponds to the atom index within the execution block.
117  * Feed the result into getSplineParamIndex() to get a full index.
118  * TODO: it's likely that both parameters can be just replaced with a single atom index, as they are derived from it.
119  * Do that, verifying that the generated code is not bloated, and/or revise the spline indexing scheme.
120  * Removing warp dependency would also be nice (and would probably coincide with removing PME_SPREADGATHER_ATOMS_PER_WARP).
121  *
122  * \tparam    order            PME order
123  * \param[in] warpIndex        Warp index wrt the block.
124  * \param[in] atomWarpIndex    Atom index wrt the warp (from 0 to PME_SPREADGATHER_ATOMS_PER_WARP - 1).
125  *
126  * \returns Index into theta or dtheta array using GPU layout.
127  */
128 template <int order>
129 int __host__ __device__ __forceinline__ getSplineParamIndexBase(int warpIndex, int atomWarpIndex)
130 {
131     assert((atomWarpIndex >= 0) && (atomWarpIndex < PME_SPREADGATHER_ATOMS_PER_WARP));
132     const int dimIndex    = 0;
133     const int splineIndex = 0;
134     // The zeroes are here to preserve the full index formula for reference
135     return (((splineIndex + order * warpIndex) * DIM + dimIndex) * PME_SPREADGATHER_ATOMS_PER_WARP + atomWarpIndex);
136 }
137
138 /*! \internal \brief
139  * Gets a unique index to an element in a spline parameter buffer (theta/dtheta),
140  * which is laid out for GPU spread/gather kernels. The index is wrt to the execution block,
141  * in range(0, atomsPerBlock * order * DIM).
142  * This function consumes result of getSplineParamIndexBase() and adjusts it for \p dimIndex and \p splineIndex.
143  *
144  * \tparam    order            PME order
145  * \param[in] paramIndexBase   Must be result of getSplineParamIndexBase().
146  * \param[in] dimIndex         Dimension index (from 0 to 2)
147  * \param[in] splineIndex      Spline contribution index (from 0 to \p order - 1)
148  *
149  * \returns Index into theta or dtheta array using GPU layout.
150  */
151 template <int order>
152 int __host__ __device__ __forceinline__ getSplineParamIndex(int paramIndexBase, int dimIndex, int splineIndex)
153 {
154     assert((dimIndex >= XX) && (dimIndex < DIM));
155     assert((splineIndex >= 0) && (splineIndex < order));
156     return (paramIndexBase + (splineIndex * DIM + dimIndex) * PME_SPREADGATHER_ATOMS_PER_WARP);
157 }
158
159 /*! \brief \internal
160  * An inline CUDA function for checking the global atom data indices against the atom data array sizes.
161  *
162  * \param[in] atomDataIndexGlobal  The atom data index.
163  * \param[in] nAtomData            The atom data array element count.
164  * \returns                        Non-0 if index is within bounds (or PME data padding is enabled), 0 otherwise.
165  *
166  * This is called from the spline_and_spread and gather PME kernels.
167  * The goal is to isolate the global range checks, and allow avoiding them with c_usePadding enabled.
168  */
169 int __device__ __forceinline__ pme_gpu_check_atom_data_index(const int atomDataIndex, const int nAtomData)
170 {
171     return c_usePadding ? 1 : (atomDataIndex < nAtomData);
172 }
173
174 /*! \brief \internal
175  * An inline CUDA function for skipping the zero-charge atoms.
176  *
177  * \returns                        Non-0 if atom should be processed, 0 otherwise.
178  * \param[in] coefficient          The atom charge.
179  *
180  * This is called from the spline_and_spread and gather PME kernels.
181  */
182 int __device__ __forceinline__ pme_gpu_check_atom_charge(const float coefficient)
183 {
184     assert(isfinite(coefficient));
185     return c_skipNeutralAtoms ? (coefficient != 0.0f) : 1;
186 }
187
188 /*! \brief \internal
189  * Given possibly large \p blockCount, returns a compact 1D or 2D grid for kernel scheduling,
190  * to minimize number of unused blocks.
191  */
192 template <typename PmeGpu>
193 dim3 __host__ inline pmeGpuCreateGrid(const PmeGpu *pmeGpu, int blockCount)
194 {
195     // How many maximum widths in X do we need (hopefully just one)
196     const int minRowCount = (blockCount + pmeGpu->maxGridWidthX - 1) / pmeGpu->maxGridWidthX;
197     // Trying to make things even
198     const int colCount = (blockCount + minRowCount - 1) / minRowCount;
199     GMX_ASSERT((colCount * minRowCount - blockCount) >= 0, "pmeGpuCreateGrid: totally wrong");
200     GMX_ASSERT((colCount * minRowCount - blockCount) < minRowCount, "pmeGpuCreateGrid: excessive blocks");
201     return dim3(colCount, minRowCount);
202 }
203
204 /*! \brief \internal
205  * The main PME CUDA-specific host data structure, included in the PME GPU structure by the archSpecific pointer.
206  */
207 struct PmeGpuCuda
208 {
209     /*! \brief The CUDA stream where everything related to the PME happens. */
210     cudaStream_t pmeStream;
211
212     /* Synchronization events */
213     /*! \brief Triggered after the grid has been copied to the host (after the spreading stage). */
214     cudaEvent_t syncSpreadGridD2H;
215
216     // TODO: consider moving some things below into the non-CUDA struct.
217
218     /* Settings which are set at the start of the run */
219     /*! \brief A boolean which tells whether the complex and real grids for cuFFT are different or same. Currenty true. */
220     bool performOutOfPlaceFFT;
221     /*! \brief A boolean which tells if the CUDA timing events are enabled.
222      *  False by default, can be enabled by setting the environment variable GMX_ENABLE_GPU_TIMING.
223      *  Note: will not be reliable when multiple GPU tasks are running concurrently on the same device context,
224      * as CUDA events on multiple streams are untrustworthy.
225      */
226     bool                                             useTiming;
227
228     std::vector<std::unique_ptr<GpuParallel3dFft > > fftSetup;
229
230     std::array<GpuRegionTimer, gtPME_EVENT_COUNT>    timingEvents;
231
232     std::set<size_t>                                 activeTimers; // indices into timingEvents
233
234     /* GPU arrays element counts (not the arrays sizes in bytes!).
235      * They might be larger than the actual meaningful data sizes.
236      * These are paired: the actual element count + the maximum element count that can fit in the current allocated memory.
237      * These integer pairs are mostly meaningful for the reallocateDeviceBuffer calls.
238      * As such, if DeviceBuffer is refactored into a class, they can be freely changed, too.
239      * The only exceptions are realGridSize and complexGridSize which are also used for grid clearing/copying.
240      * TODO: these should live in a clean buffered container type, and be refactored in the NB/cudautils as well.
241      */
242     /*! \brief The kernelParams.atoms.coordinates float element count (actual)*/
243     int coordinatesSize;
244     /*! \brief The kernelParams.atoms.coordinates float element count (reserved) */
245     int coordinatesSizeAlloc;
246     /*! \brief The kernelParams.atoms.forces float element count (actual) */
247     int forcesSize;
248     /*! \brief The kernelParams.atoms.forces float element count (reserved) */
249     int forcesSizeAlloc;
250     /*! \brief The kernelParams.atoms.gridlineIndices int element count (actual) */
251     int gridlineIndicesSize;
252     /*! \brief The kernelParams.atoms.gridlineIndices int element count (reserved) */
253     int gridlineIndicesSizeAlloc;
254     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (actual) */
255     int splineDataSize;
256     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (reserved) */
257     int splineDataSizeAlloc;
258     /*! \brief The kernelParams.atoms.coefficients float element count (actual) */
259     int coefficientsSize;
260     /*! \brief The kernelParams.atoms.coefficients float element count (reserved) */
261     int coefficientsSizeAlloc;
262     /*! \brief The kernelParams.grid.splineValuesArray float element count (actual) */
263     int splineValuesSize;
264     /*! \brief The kernelParams.grid.splineValuesArray float element count (reserved) */
265     int splineValuesSizeAlloc;
266     /*! \brief The kernelParams.grid.realGrid float element count (actual) */
267     int realGridSize;
268     /*! \brief The kernelParams.grid.realGrid float element count (reserved) */
269     int realGridSizeAlloc;
270     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (actual) */
271     int complexGridSize;
272     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (reserved) */
273     int complexGridSizeAlloc;
274 };
275
276
277 /*! \brief \internal
278  * A single structure encompassing all the PME data used in CUDA kernels.
279  * This inherits from PmeGpuKernelParamsBase and adds a couple cudaTextureObject_t handles,
280  * which we would like to avoid in plain C++.
281  */
282 struct PmeGpuCudaKernelParams : PmeGpuKernelParamsBase
283 {
284     /* These are CUDA texture objects, related to the grid size. */
285     /*! \brief CUDA texture object for accessing grid.d_fractShiftsTable */
286     cudaTextureObject_t fractShiftsTableTexture;
287     /*! \brief CUDA texture object for accessing grid.d_gridlineIndicesTable */
288     cudaTextureObject_t gridlineIndicesTableTexture;
289 };
290
291 #endif