3749d5748b91af38514be5d22fd3f7b14592a564
[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,2019, 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
50  * match. 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,
78               "DeviceBuffer is defined as an 8 byte stub for OpenCL C");
79 static_assert(sizeof(DeviceBuffer<int>) == 8,
80               "DeviceBuffer is defined as an 8 byte stub for OpenCL C");
81 #else
82 #    define HIDE_FROM_OPENCL_COMPILER(x) char8
83 #endif
84
85 /* What follows is all the PME GPU function arguments,
86  * sorted into several device-side structures depending on the update rate.
87  * This is GPU agnostic (float3 replaced by float[3], etc.).
88  * The GPU-framework specifics (e.g. cudaTextureObject_t handles) are described
89  * in the larger structure PmeGpuCudaKernelParams in the pme.cuh.
90  */
91
92 /*! \internal \brief
93  * A GPU data structure for storing the constant PME data.
94  * This only has to be initialized once.
95  */
96 struct PmeGpuConstParams
97 {
98     /*! \brief Electrostatics coefficient = ONE_4PI_EPS0 / pme->epsilon_r */
99     float elFactor;
100     /*! \brief Virial and energy GPU array. Size is PME_GPU_ENERGY_AND_VIRIAL_COUNT (7) floats.
101      * The element order is virxx, viryy, virzz, virxy, virxz, viryz, energy. */
102     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_virialAndEnergy;
103 };
104
105 /*! \internal \brief
106  * A GPU data structure for storing the PME data related to the grid sizes and cut-off.
107  * This only has to be updated at every DD step.
108  */
109 struct PmeGpuGridParams
110 {
111     /*! \brief Ewald solving factor = (M_PI / pme->ewaldcoeff_q)^2 */
112     float ewaldFactor;
113
114     /* Grid sizes */
115     /*! \brief Real-space grid data dimensions. */
116     int realGridSize[DIM];
117     /*! \brief Real-space grid dimensions, only converted to floating point. */
118     float realGridSizeFP[DIM];
119     /*! \brief Real-space grid dimensions (padded). The padding as compared to realGridSize includes the (order - 1) overlap. */
120     int realGridSizePadded[DIM]; /* Is major dimension of this ever used in kernels? */
121     /*! \brief Fourier grid dimensions. This counts the complex numbers! */
122     int complexGridSize[DIM];
123     /*! \brief Fourier grid dimensions (padded). This counts the complex numbers! */
124     int complexGridSizePadded[DIM];
125
126     /*! \brief Offsets for X/Y/Z components of d_splineModuli */
127     int splineValuesOffset[DIM];
128     /*! \brief Offsets for X/Y/Z components of d_fractShiftsTable and d_gridlineIndicesTable */
129     int tablesOffsets[DIM];
130
131     /* Grid arrays */
132     /*! \brief Real space grid. */
133     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_realGrid;
134     /*! \brief Complex grid - used in FFT/solve. If inplace cu/clFFT is used, then it is the same handle as realGrid. */
135     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_fourierGrid;
136
137     /*! \brief Grid spline values as in pme->bsp_mod
138      * (laid out sequentially (XXX....XYYY......YZZZ.....Z))
139      */
140     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_splineModuli;
141     /*! \brief Fractional shifts lookup table as in pme->fshx/fshy/fshz, laid out sequentially (XXX....XYYY......YZZZ.....Z) */
142     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_fractShiftsTable;
143     /*! \brief Gridline indices lookup table
144      * (modulo lookup table as in pme->nnx/nny/nnz, laid out sequentially (XXX....XYYY......YZZZ.....Z)) */
145     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<int>) d_gridlineIndicesTable;
146 };
147
148 /*! \internal \brief
149  * A GPU data structure for storing the PME data of the atoms, local to this process' domain
150  * partition. This only has to be updated every DD step.
151  */
152 struct PmeGpuAtomParams
153 {
154     /*! \brief Number of local atoms */
155     int nAtoms;
156     /*! \brief Global GPU memory array handle with input rvec atom coordinates.
157      * The coordinates themselves change and need to be copied to the GPU for every PME computation,
158      * but reallocation happens only at DD.
159      */
160     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_coordinates;
161     /*! \brief Global GPU memory array handle with input atom charges.
162      * The charges only need to be reallocated and copied to the GPU at DD step.
163      */
164     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_coefficients;
165     /*! \brief Global GPU memory array handle with input/output rvec atom forces.
166      * The forces change and need to be copied from (and possibly to) the GPU for every PME
167      * computation, but reallocation happens only at DD.
168      */
169     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_forces;
170     /*! \brief Global GPU memory array handle with ivec atom gridline indices.
171      * Computed on GPU in the spline calculation part.
172      */
173     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<int>) d_gridlineIndices;
174     /* B-spline parameters are computed entirely on GPU for every PME computation, not copied.
175      * Unless we want to try something like GPU spread + CPU gather?
176      */
177     /*! \brief Global GPU memory array handle with B-spline values */
178     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_theta;
179     /*! \brief Global GPU memory array handle with B-spline derivative values */
180     HIDE_FROM_OPENCL_COMPILER(DeviceBuffer<float>) d_dtheta;
181 };
182
183 /*! \internal \brief
184  * A GPU data structure for storing the PME data which might change for each new PME computation.
185  */
186 struct PmeGpuDynamicParams
187 {
188     /* The box parameters. The box only changes size with pressure coupling enabled. */
189     /*! \brief
190      * Reciprocal (inverted unit cell) box.
191      *
192      * The box is transposed as compared to the CPU pme->recipbox.
193      * Basically, spread uses matrix columns (while solve and gather use rows).
194      * This storage format might be not the most optimal since the box is always triangular so there are zeroes.
195      */
196     float recipBox[DIM][DIM];
197     /*! \brief The unit cell volume for solving. */
198     float boxVolume;
199 };
200
201 /*! \internal \brief
202  * A single structure encompassing almost all the PME data used in GPU kernels on device.
203  * This is inherited by the GPU framework-specific structure
204  * (PmeGpuCudaKernelParams in pme.cuh).
205  * This way, most code preparing the kernel parameters can be GPU-agnostic by casting
206  * the kernel parameter data pointer to PmeGpuKernelParamsBase.
207  */
208 struct PmeGpuKernelParamsBase
209 {
210     /*! \brief Constant data that is set once. */
211     struct PmeGpuConstParams constants;
212     /*! \brief Data dependent on the grid size/cutoff. */
213     struct PmeGpuGridParams grid;
214     /*! \brief Data dependent on the DD and local atoms. */
215     struct PmeGpuAtomParams atoms;
216     /*! \brief Data that possibly changes for every new PME computation.
217      * This should be kept up-to-date by calling pme_gpu_prepare_computation(...)
218      * before launching spreading.
219      */
220     struct PmeGpuDynamicParams current;
221 };
222
223 #endif