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