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