6fb6cae82726120279bb3a0a6329949fac3be110
[alexxy/gromacs.git] / src / gromacs / ewald / pme-gpu-types.h
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 Defines the PME GPU data structures
38  * (the GPU function parameters used both on host and device sides).
39  *
40  * \author Aleksei Iupinov <a.yupinov@gmail.com>
41  * \ingroup module_ewald
42  */
43
44 #ifndef GMX_EWALD_PME_GPU_TYPES_H
45 #define GMX_EWALD_PME_GPU_TYPES_H
46
47 /*
48  * In OpenCL, the structures must be laid out on the host and device exactly the same way.
49  * If something is off, one might get an error CL_INVALID_ARG_SIZE if any structure's sizes don't match.
50  * What's worse, structures might be of same size but members might be aligned differently,
51  * resulting in wrong kernel results. The structures below are aligned manually.
52  * The pattern is ordering the members of structs from smallest to largest sizeof
53  * (arrays behave the same way as sequences of separate fields),
54  * as described in "The Lost Art of C Structure Packing".
55  *
56  * However, if the need arises at some point, they can all be aligned forcefully:
57  *
58  * #define GMX_GPU_ALIGNED __attribute__ ((aligned(8)))
59  * struct GMX_GPU_ALIGNED PmeGpuConstParams
60  * struct GMX_GPU_ALIGNED PmeGpuGridParams
61  * etc...
62  *
63  * One might also try __attribute__ ((packed)), but it doesn't work with DeviceBuffer,
64  * as it appears to not be POD.
65  */
66
67
68 /*! \brief A workaround to hide DeviceBuffer template from OpenCL kernel compilation
69  * - to turn it into a dummy of the same size as host implementation of device buffer.
70  * As we only care about 64-bit, 8 bytes is fine.
71  * TODO: what we should be doing is providing separate device-side views of the same structures -
72  * then there would be no need for macro.
73  */
74 #ifndef __OPENCL_C_VERSION__
75 #include "gromacs/gpu_utils/devicebuffer.h"
76 #define HIDE_FROM_OPENCL_COMPILER(x) x
77 static_assert(sizeof(DeviceBuffer<float>) == 8, "DeviceBuffer is defined as an 8 byte stub for OpenCL C");
78 static_assert(sizeof(DeviceBuffer<int>) == 8, "DeviceBuffer is defined as an 8 byte stub for OpenCL C");
79 #else
80 #define HIDE_FROM_OPENCL_COMPILER(x) char8
81 #endif
82
83 /* What follows is all the PME GPU function arguments,
84  * sorted into several device-side structures depending on the update rate.
85  * This is GPU agnostic (float3 replaced by float[3], etc.).
86  * The GPU-framework specifics (e.g. cudaTextureObject_t handles) are described
87  * in the larger structure PmeGpuCudaKernelParams in the pme.cuh.
88  */
89
90 /*! \internal \brief
91  * A GPU data structure for storing the constant PME data.
92  * This only has to be initialized once.
93  */
94 struct PmeGpuConstParams
95 {
96     /*! \brief Electrostatics coefficient = ONE_4PI_EPS0 / pme->epsilon_r */
97     float elFactor;
98     /*! \brief Virial and energy GPU array. Size is PME_GPU_ENERGY_AND_VIRIAL_COUNT (7) floats.
99      * The element order is virxx, viryy, virzz, virxy, virxz, viryz, energy. */
100     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_virialAndEnergy;
101 };
102
103 /*! \internal \brief
104  * A GPU data structure for storing the PME data related to the grid sizes and cut-off.
105  * This only has to be updated at every DD step.
106  */
107 struct PmeGpuGridParams
108 {
109     /*! \brief Ewald solving factor = (M_PI / pme->ewaldcoeff_q)^2 */
110     float ewaldFactor;
111
112     /* Grid sizes */
113     /*! \brief Real-space grid data dimensions. */
114     int   realGridSize[DIM];
115     /*! \brief Real-space grid dimensions, only converted to floating point. */
116     float realGridSizeFP[DIM];
117     /*! \brief Real-space grid dimensions (padded). The padding as compared to realGridSize includes the (order - 1) overlap. */
118     int   realGridSizePadded[DIM]; /* Is major dimension of this ever used in kernels? */
119     /*! \brief Fourier grid dimensions. This counts the complex numbers! */
120     int   complexGridSize[DIM];
121     /*! \brief Fourier grid dimensions (padded). This counts the complex numbers! */
122     int   complexGridSizePadded[DIM];
123
124     /*! \brief Offsets for X/Y/Z components of d_splineModuli */
125     int  splineValuesOffset[DIM];
126     /*! \brief Offsets for X/Y/Z components of d_fractShiftsTable and d_gridlineIndicesTable */
127     int  tablesOffsets[DIM];
128
129     /* Grid arrays */
130     /*! \brief Real space grid. */
131     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_realGrid;
132     /*! \brief Complex grid - used in FFT/solve. If inplace cu/clFFT is used, then it is the same handle as realGrid. */
133     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_fourierGrid;
134
135     /*! \brief Grid spline values as in pme->bsp_mod
136      * (laid out sequentially (XXX....XYYY......YZZZ.....Z))
137      */
138     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_splineModuli;
139     /*! \brief Fractional shifts lookup table as in pme->fshx/fshy/fshz, laid out sequentially (XXX....XYYY......YZZZ.....Z) */
140     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_fractShiftsTable;
141     /*! \brief Gridline indices lookup table
142      * (modulo lookup table as in pme->nnx/nny/nnz, laid out sequentially (XXX....XYYY......YZZZ.....Z)) */
143     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<int>) d_gridlineIndicesTable;
144 };
145
146 /*! \internal \brief
147  * A GPU data structure for storing the PME data of the atoms, local to this process' domain partition.
148  * This only has to be updated every DD step.
149  */
150 struct PmeGpuAtomParams
151 {
152     /*! \brief Number of local atoms */
153     int    nAtoms;
154     /*! \brief Global GPU memory array handle with input rvec atom coordinates.
155      * The coordinates themselves change and need to be copied to the GPU for every PME computation,
156      * but reallocation happens only at DD.
157      */
158     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_coordinates;
159     /*! \brief Global GPU memory array handle with input atom charges.
160      * The charges only need to be reallocated and copied to the GPU at DD step.
161      */
162     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_coefficients;
163     /*! \brief Global GPU memory array handle with input/output rvec atom forces.
164      * The forces change and need to be copied from (and possibly to) the GPU for every PME computation,
165      * but reallocation happens only at DD.
166      */
167     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_forces;
168     /*! \brief Global GPU memory array handle with ivec atom gridline indices.
169      * Computed on GPU in the spline calculation part.
170      */
171     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<int>) d_gridlineIndices;
172     /* B-spline parameters are computed entirely on GPU for every PME computation, not copied.
173      * Unless we want to try something like GPU spread + CPU gather?
174      */
175     /*! \brief Global GPU memory array handle with B-spline values */
176     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_theta;
177     /*! \brief Global GPU memory array handle with B-spline derivative values */
178     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_dtheta;
179 };
180
181 /*! \internal \brief
182  * A GPU data structure for storing the PME data which might change for each new PME computation.
183  */
184 struct PmeGpuDynamicParams
185 {
186     /* The box parameters. The box only changes size with pressure coupling enabled. */
187     /*! \brief
188      * Reciprocal (inverted unit cell) box.
189      *
190      * The box is transposed as compared to the CPU pme->recipbox.
191      * Basically, spread uses matrix columns (while solve and gather use rows).
192      * This storage format might be not the most optimal since the box is always triangular so there are zeroes.
193      */
194     float  recipBox[DIM][DIM];
195     /*! \brief The unit cell volume for solving. */
196     float  boxVolume;
197 };
198
199 /*! \internal \brief
200  * A single structure encompassing almost all the PME data used in GPU kernels on device.
201  * This is inherited by the GPU framework-specific structure
202  * (PmeGpuCudaKernelParams in pme.cuh).
203  * This way, most code preparing the kernel parameters can be GPU-agnostic by casting
204  * the kernel parameter data pointer to PmeGpuKernelParamsBase.
205  */
206 struct PmeGpuKernelParamsBase
207 {
208     /*! \brief Constant data that is set once. */
209     struct PmeGpuConstParams   constants;
210     /*! \brief Data dependent on the grid size/cutoff. */
211     struct PmeGpuGridParams    grid;
212     /*! \brief Data dependent on the DD and local atoms. */
213     struct PmeGpuAtomParams    atoms;
214     /*! \brief Data that possibly changes for every new PME computation.
215      * This should be kept up-to-date by calling pme_gpu_prepare_computation(...)
216      * before launching spreading.
217      */
218     struct PmeGpuDynamicParams current;
219 };
220
221 #endif