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