Gather the PME GPU constants/macros in a single header
[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 "pme-gpu-constants.h"
54 #include "pme-gpu-internal.h"
55 #include "pme-gpu-types.h"
56 #include "pme-gpu-types-host.h"
57 #include "pme-timings.cuh"
58
59 class GpuParallel3dFft;
60
61 /*! \internal \brief
62  * Gets a base of the unique index to an element in a spline parameter buffer (theta/dtheta),
63  * which is laid out for GPU spread/gather kernels. The base only corresponds to the atom index within the execution block.
64  * Feed the result into getSplineParamIndex() to get a full index.
65  * TODO: it's likely that both parameters can be just replaced with a single atom index, as they are derived from it.
66  * Do that, verifying that the generated code is not bloated, and/or revise the spline indexing scheme.
67  * Removing warp dependency would also be nice (and would probably coincide with removing PME_SPREADGATHER_ATOMS_PER_WARP).
68  *
69  * \tparam    order            PME order
70  * \param[in] warpIndex        Warp index wrt the block.
71  * \param[in] atomWarpIndex    Atom index wrt the warp (from 0 to PME_SPREADGATHER_ATOMS_PER_WARP - 1).
72  *
73  * \returns Index into theta or dtheta array using GPU layout.
74  */
75 template <int order>
76 int __host__ __device__ __forceinline__ getSplineParamIndexBase(int warpIndex, int atomWarpIndex)
77 {
78     assert((atomWarpIndex >= 0) && (atomWarpIndex < PME_SPREADGATHER_ATOMS_PER_WARP));
79     const int dimIndex    = 0;
80     const int splineIndex = 0;
81     // The zeroes are here to preserve the full index formula for reference
82     return (((splineIndex + order * warpIndex) * DIM + dimIndex) * PME_SPREADGATHER_ATOMS_PER_WARP + atomWarpIndex);
83 }
84
85 /*! \internal \brief
86  * Gets a unique index to an element in a spline parameter buffer (theta/dtheta),
87  * which is laid out for GPU spread/gather kernels. The index is wrt to the execution block,
88  * in range(0, atomsPerBlock * order * DIM).
89  * This function consumes result of getSplineParamIndexBase() and adjusts it for \p dimIndex and \p splineIndex.
90  *
91  * \tparam    order            PME order
92  * \param[in] paramIndexBase   Must be result of getSplineParamIndexBase().
93  * \param[in] dimIndex         Dimension index (from 0 to 2)
94  * \param[in] splineIndex      Spline contribution index (from 0 to \p order - 1)
95  *
96  * \returns Index into theta or dtheta array using GPU layout.
97  */
98 template <int order>
99 int __host__ __device__ __forceinline__ getSplineParamIndex(int paramIndexBase, int dimIndex, int splineIndex)
100 {
101     assert((dimIndex >= XX) && (dimIndex < DIM));
102     assert((splineIndex >= 0) && (splineIndex < order));
103     return (paramIndexBase + (splineIndex * DIM + dimIndex) * PME_SPREADGATHER_ATOMS_PER_WARP);
104 }
105
106 /*! \brief \internal
107  * An inline CUDA function for checking the global atom data indices against the atom data array sizes.
108  *
109  * \param[in] atomDataIndexGlobal  The atom data index.
110  * \param[in] nAtomData            The atom data array element count.
111  * \returns                        Non-0 if index is within bounds (or PME data padding is enabled), 0 otherwise.
112  *
113  * This is called from the spline_and_spread and gather PME kernels.
114  * The goal is to isolate the global range checks, and allow avoiding them with c_usePadding enabled.
115  */
116 int __device__ __forceinline__ pme_gpu_check_atom_data_index(const int atomDataIndex, const int nAtomData)
117 {
118     return c_usePadding ? 1 : (atomDataIndex < nAtomData);
119 }
120
121 /*! \brief \internal
122  * An inline CUDA function for skipping the zero-charge atoms.
123  *
124  * \returns                        Non-0 if atom should be processed, 0 otherwise.
125  * \param[in] coefficient          The atom charge.
126  *
127  * This is called from the spline_and_spread and gather PME kernels.
128  */
129 int __device__ __forceinline__ pme_gpu_check_atom_charge(const float coefficient)
130 {
131     assert(isfinite(coefficient));
132     return c_skipNeutralAtoms ? (coefficient != 0.0f) : 1;
133 }
134
135 /*! \brief \internal
136  * Given possibly large \p blockCount, returns a compact 1D or 2D grid for kernel scheduling,
137  * to minimize number of unused blocks.
138  */
139 template <typename PmeGpu>
140 dim3 __host__ inline pmeGpuCreateGrid(const PmeGpu *pmeGpu, int blockCount)
141 {
142     // How many maximum widths in X do we need (hopefully just one)
143     const int minRowCount = (blockCount + pmeGpu->maxGridWidthX - 1) / pmeGpu->maxGridWidthX;
144     // Trying to make things even
145     const int colCount = (blockCount + minRowCount - 1) / minRowCount;
146     GMX_ASSERT((colCount * minRowCount - blockCount) >= 0, "pmeGpuCreateGrid: totally wrong");
147     GMX_ASSERT((colCount * minRowCount - blockCount) < minRowCount, "pmeGpuCreateGrid: excessive blocks");
148     return dim3(colCount, minRowCount);
149 }
150
151 /*! \brief \internal
152  * The main PME CUDA-specific host data structure, included in the PME GPU structure by the archSpecific pointer.
153  */
154 struct PmeGpuCuda
155 {
156     /*! \brief The CUDA stream where everything related to the PME happens. */
157     cudaStream_t pmeStream;
158
159     /*! \brief
160      * A handle to the GPU context.
161      * TODO: this is currently extracted from the implementation of pmeGpu->programHandle_,
162      * but should be a constructor parameter to PmeGpu, as well as PmeGpuProgram,
163      * managed by high-level code.
164      */
165     Context context;
166
167     /* Synchronization events */
168     /*! \brief Triggered after the grid has been copied to the host (after the spreading stage). */
169     cudaEvent_t syncSpreadGridD2H;
170
171     // TODO: consider moving some things below into the non-CUDA struct.
172
173     /* Settings which are set at the start of the run */
174     /*! \brief A boolean which tells whether the complex and real grids for cuFFT are different or same. Currenty true. */
175     bool performOutOfPlaceFFT;
176     /*! \brief A boolean which tells if the CUDA timing events are enabled.
177      *  False by default, can be enabled by setting the environment variable GMX_ENABLE_GPU_TIMING.
178      *  Note: will not be reliable when multiple GPU tasks are running concurrently on the same device context,
179      * as CUDA events on multiple streams are untrustworthy.
180      */
181     bool                                             useTiming;
182
183     std::vector<std::unique_ptr<GpuParallel3dFft > > fftSetup;
184
185     std::array<GpuRegionTimer, gtPME_EVENT_COUNT>    timingEvents;
186
187     std::set<size_t>                                 activeTimers; // indices into timingEvents
188
189     /* GPU arrays element counts (not the arrays sizes in bytes!).
190      * They might be larger than the actual meaningful data sizes.
191      * These are paired: the actual element count + the maximum element count that can fit in the current allocated memory.
192      * These integer pairs are mostly meaningful for the reallocateDeviceBuffer calls.
193      * As such, if DeviceBuffer is refactored into a class, they can be freely changed, too.
194      * The only exceptions are realGridSize and complexGridSize which are also used for grid clearing/copying.
195      * TODO: these should live in a clean buffered container type, and be refactored in the NB/cudautils as well.
196      */
197     /*! \brief The kernelParams.atoms.coordinates float element count (actual)*/
198     int coordinatesSize;
199     /*! \brief The kernelParams.atoms.coordinates float element count (reserved) */
200     int coordinatesSizeAlloc;
201     /*! \brief The kernelParams.atoms.forces float element count (actual) */
202     int forcesSize;
203     /*! \brief The kernelParams.atoms.forces float element count (reserved) */
204     int forcesSizeAlloc;
205     /*! \brief The kernelParams.atoms.gridlineIndices int element count (actual) */
206     int gridlineIndicesSize;
207     /*! \brief The kernelParams.atoms.gridlineIndices int element count (reserved) */
208     int gridlineIndicesSizeAlloc;
209     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (actual) */
210     int splineDataSize;
211     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (reserved) */
212     int splineDataSizeAlloc;
213     /*! \brief The kernelParams.atoms.coefficients float element count (actual) */
214     int coefficientsSize;
215     /*! \brief The kernelParams.atoms.coefficients float element count (reserved) */
216     int coefficientsSizeAlloc;
217     /*! \brief The kernelParams.grid.splineValuesArray float element count (actual) */
218     int splineValuesSize;
219     /*! \brief The kernelParams.grid.splineValuesArray float element count (reserved) */
220     int splineValuesSizeAlloc;
221     /*! \brief The kernelParams.grid.realGrid float element count (actual) */
222     int realGridSize;
223     /*! \brief The kernelParams.grid.realGrid float element count (reserved) */
224     int realGridSizeAlloc;
225     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (actual) */
226     int complexGridSize;
227     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (reserved) */
228     int complexGridSizeAlloc;
229 };
230
231
232 /*! \brief \internal
233  * A single structure encompassing all the PME data used in CUDA kernels.
234  * This inherits from PmeGpuKernelParamsBase and adds a couple cudaTextureObject_t handles,
235  * which we would like to avoid in plain C++.
236  */
237 struct PmeGpuCudaKernelParams : PmeGpuKernelParamsBase
238 {
239     /* These are CUDA texture objects, related to the grid size. */
240     /*! \brief CUDA texture object for accessing grid.d_fractShiftsTable */
241     cudaTextureObject_t fractShiftsTableTexture;
242     /*! \brief CUDA texture object for accessing grid.d_gridlineIndicesTable */
243     cudaTextureObject_t gridlineIndicesTableTexture;
244 };
245
246 #endif