Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / ewald / pme_solve.clh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 Implements PME OpenCL solving kernel.
38  * When including this and other PME OpenCL kernel files, plenty of common
39  * constants/macros are expected to be defined.
40  * For details, please see how pme_program.cl is compiled in pme_gpu_program_impl_ocl.cpp.
41  *
42  * This file's solving kernel specifically expects the following definitions for its flavors:
43  *
44  * - gridOrdering must be defined to either XYZ or YZX
45  * and corresponds to the dimension order of the grid (GridOrdering enum in CUDA kernels);
46  * - computeEnergyAndVirial must evaluate to true or false, and expresses
47  * whether the reduction is performed.
48  *
49  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
50  */
51
52 #include "pme_gpu_types.h"
53 #include "gromacs/gpu_utils/vectype_ops.clh"
54
55 /*! \brief
56  * PME complex grid solver kernel function.
57  * Please see the file description for additional defines which this kernel expects.
58  *
59  * \param[in]     kernelParams         Input PME GPU data in constant memory.
60  * \param[in]     gm_splineModuli      B-Spline moduli.
61  * \param[out]    gm_virialAndEnergy   Reduced virial and enrgy (only with computeEnergyAndVirial ==
62  * true) \param[in,out] gm_grid              Fourier grid to transform.
63  */
64 __attribute__((work_group_size_hint(c_solveMaxWorkGroupSize, 1, 1)))
65 __kernel void CUSTOMIZED_KERNEL_NAME(pme_solve_kernel)(const struct PmeOpenCLKernelParams kernelParams,
66                                                        __global const float* __restrict__ gm_splineModuli,
67                                                        __global float* __restrict__ gm_virialAndEnergy,
68                                                        __global float2* __restrict__ gm_grid)
69 {
70     /* This kernel supports 2 different grid dimension orderings: YZX and XYZ */
71     int majorDim, middleDim, minorDim;
72     if (gridOrdering == YZX)
73     {
74         majorDim  = YY;
75         middleDim = ZZ;
76         minorDim  = XX;
77     }
78     if (gridOrdering == XYZ)
79     {
80         majorDim  = XX;
81         middleDim = YY;
82         minorDim  = ZZ;
83     }
84
85     __global const float* __restrict__ gm_splineValueMajor =
86             gm_splineModuli + kernelParams.grid.splineValuesOffset[majorDim];
87     __global const float* __restrict__ gm_splineValueMiddle =
88             gm_splineModuli + kernelParams.grid.splineValuesOffset[middleDim];
89     __global const float* __restrict__ gm_splineValueMinor =
90             gm_splineModuli + kernelParams.grid.splineValuesOffset[minorDim];
91
92     /* Various grid sizes and indices */
93     const int localOffsetMinor = 0, localOffsetMajor = 0, localOffsetMiddle = 0; // unused
94     const int localSizeMinor   = kernelParams.grid.complexGridSizePadded[minorDim];
95     const int localSizeMiddle  = kernelParams.grid.complexGridSizePadded[middleDim];
96     const int localCountMiddle = kernelParams.grid.complexGridSize[middleDim];
97     const int localCountMinor  = kernelParams.grid.complexGridSize[minorDim];
98     const int nMajor           = kernelParams.grid.realGridSize[majorDim];
99     const int nMiddle          = kernelParams.grid.realGridSize[middleDim];
100     const int nMinor           = kernelParams.grid.realGridSize[minorDim];
101     const int maxkMajor        = (nMajor + 1) / 2;  // X or Y
102     const int maxkMiddle       = (nMiddle + 1) / 2; // Y OR Z => only check for !YZX
103     const int maxkMinor        = (nMinor + 1) / 2;  // Z or X => only check for YZX
104
105     /* Each thread works on one cell of the Fourier space complex 3D grid (gm_grid).
106      * Each block handles up to c_solveMaxWorkGroupSize cells -
107      * depending on the grid contiguous dimension size,
108      * that can range from a part of a single gridline to several complete gridlines.
109      */
110     const int threadLocalId     = get_local_id(XX);
111     const int gridLineSize      = localCountMinor;
112     const int gridLineIndex     = threadLocalId / gridLineSize;
113     const int gridLineCellIndex = threadLocalId - gridLineSize * gridLineIndex;
114     const int gridLinesPerBlock = max((int)(get_local_size(XX)) / gridLineSize, 1);
115     const int activeWarps       = (get_local_size(XX) / warp_size);
116     const int indexMinor        = get_group_id(XX) * get_local_size(XX) + gridLineCellIndex;
117     const int indexMiddle       = get_group_id(YY) * gridLinesPerBlock + gridLineIndex;
118     const int indexMajor        = get_group_id(ZZ);
119
120     /* Optional outputs */
121     float energy = 0.0f;
122     float virxx  = 0.0f;
123     float virxy  = 0.0f;
124     float virxz  = 0.0f;
125     float viryy  = 0.0f;
126     float viryz  = 0.0f;
127     float virzz  = 0.0f;
128
129     assert(indexMajor < kernelParams.grid.complexGridSize[majorDim]);
130     if ((indexMiddle < localCountMiddle) & (indexMinor < localCountMinor)
131         & (gridLineIndex < gridLinesPerBlock))
132     {
133         /* The offset should be equal to the global thread index for coalesced access */
134         const int gridIndex = (indexMajor * localSizeMiddle + indexMiddle) * localSizeMinor + indexMinor;
135         __global float2* __restrict__ gm_gridCell = gm_grid + gridIndex;
136
137         const int kMajor = indexMajor + localOffsetMajor;
138         /* Checking either X in XYZ, or Y in YZX cases */
139         const float mMajor = (kMajor < maxkMajor) ? kMajor : (kMajor - nMajor);
140
141         const int kMiddle = indexMiddle + localOffsetMiddle;
142         float     mMiddle = kMiddle;
143         /* Checking Y in XYZ case */
144         if (gridOrdering == XYZ)
145         {
146             mMiddle = (kMiddle < maxkMiddle) ? kMiddle : (kMiddle - nMiddle);
147         }
148         const int kMinor = localOffsetMinor + indexMinor;
149         float     mMinor = kMinor;
150         /* Checking X in YZX case */
151         if (gridOrdering == YZX)
152         {
153             mMinor = (kMinor < maxkMinor) ? kMinor : (kMinor - nMinor);
154         }
155         /* We should skip the k-space point (0,0,0) */
156         const bool notZeroPoint = (kMinor > 0) | (kMajor > 0) | (kMiddle > 0);
157
158         float mX, mY, mZ;
159         if (gridOrdering == YZX)
160         {
161             mX = mMinor;
162             mY = mMajor;
163             mZ = mMiddle;
164         }
165         if (gridOrdering == XYZ)
166         {
167             mX = mMajor;
168             mY = mMiddle;
169             mZ = mMinor;
170         }
171
172         /* 0.5 correction factor for the first and last components of a Z dimension */
173         float corner_fac = 1.0f;
174         if (gridOrdering == YZX)
175         {
176             if ((kMiddle == 0) | (kMiddle == maxkMiddle))
177             {
178                 corner_fac = 0.5f;
179             }
180         }
181         if (gridOrdering == XYZ)
182         {
183             if ((kMinor == 0) | (kMinor == maxkMinor))
184             {
185                 corner_fac = 0.5f;
186             }
187         }
188
189         if (notZeroPoint)
190         {
191             const float mhxk = mX * kernelParams.current.recipBox[XX][XX];
192             const float mhyk = mX * kernelParams.current.recipBox[XX][YY]
193                                + mY * kernelParams.current.recipBox[YY][YY];
194             const float mhzk = mX * kernelParams.current.recipBox[XX][ZZ]
195                                + mY * kernelParams.current.recipBox[YY][ZZ]
196                                + mZ * kernelParams.current.recipBox[ZZ][ZZ];
197
198             const float m2k = mhxk * mhxk + mhyk * mhyk + mhzk * mhzk;
199             assert(m2k != 0.0f);
200             const float denom = m2k * M_PI_F * kernelParams.current.boxVolume * gm_splineValueMajor[kMajor]
201                                 * gm_splineValueMiddle[kMiddle] * gm_splineValueMinor[kMinor];
202             assert(isfinite(denom));
203             assert(denom != 0.0f);
204             const float tmp1   = exp(-kernelParams.grid.ewaldFactor * m2k);
205             const float etermk = kernelParams.constants.elFactor * tmp1 / denom;
206
207             float2       gridValue    = *gm_gridCell;
208             const float2 oldGridValue = gridValue;
209
210             gridValue.x *= etermk;
211             gridValue.y *= etermk;
212             *gm_gridCell = gridValue;
213
214             if (computeEnergyAndVirial)
215             {
216                 const float tmp1k =
217                         2.0f * (gridValue.x * oldGridValue.x + gridValue.y * oldGridValue.y);
218
219                 const float vfactor = (kernelParams.grid.ewaldFactor + 1.0f / m2k) * 2.0f;
220                 const float ets2    = corner_fac * tmp1k;
221                 energy              = ets2;
222
223                 const float ets2vf = ets2 * vfactor;
224
225                 virxx = ets2vf * mhxk * mhxk - ets2;
226                 virxy = ets2vf * mhxk * mhyk;
227                 virxz = ets2vf * mhxk * mhzk;
228                 viryy = ets2vf * mhyk * mhyk - ets2;
229                 viryz = ets2vf * mhyk * mhzk;
230                 virzz = ets2vf * mhzk * mhzk - ets2;
231             }
232         }
233     }
234
235     // This is only for reduction below. OpenCL 1.2: all local memory must be declared on kernel scope.
236     __local float sm_virialAndEnergy[c_virialAndEnergyCount * warp_size];
237
238     /* Optional energy/virial reduction */
239     if (computeEnergyAndVirial)
240     {
241         // TODO: implement AMD intrinsics reduction, like with shuffles in CUDA version. #2514
242
243         /* Shared memory reduction with atomics.
244          * Each component is first reduced into warp_size positions in the shared memory;
245          * Then first few warps reduce everything further and add to the global memory.
246          * This can likely be improved, but is anyway faster than the previous straightforward reduction,
247          * which was using too much shared memory (for storing all 7 floats on each thread).
248          */
249
250         const int  lane      = threadLocalId & (warp_size - 1);
251         const int  warpIndex = threadLocalId / warp_size;
252         const bool firstWarp = (warpIndex == 0);
253         if (firstWarp)
254         {
255             sm_virialAndEnergy[0 * warp_size + lane] = virxx;
256             sm_virialAndEnergy[1 * warp_size + lane] = viryy;
257             sm_virialAndEnergy[2 * warp_size + lane] = virzz;
258             sm_virialAndEnergy[3 * warp_size + lane] = virxy;
259             sm_virialAndEnergy[4 * warp_size + lane] = virxz;
260             sm_virialAndEnergy[5 * warp_size + lane] = viryz;
261             sm_virialAndEnergy[6 * warp_size + lane] = energy;
262         }
263         barrier(CLK_LOCAL_MEM_FENCE);
264         if (!firstWarp)
265         {
266             atomicAdd_l_f(sm_virialAndEnergy + 0 * warp_size + lane, virxx);
267             atomicAdd_l_f(sm_virialAndEnergy + 1 * warp_size + lane, viryy);
268             atomicAdd_l_f(sm_virialAndEnergy + 2 * warp_size + lane, virzz);
269             atomicAdd_l_f(sm_virialAndEnergy + 3 * warp_size + lane, virxy);
270             atomicAdd_l_f(sm_virialAndEnergy + 4 * warp_size + lane, virxz);
271             atomicAdd_l_f(sm_virialAndEnergy + 5 * warp_size + lane, viryz);
272             atomicAdd_l_f(sm_virialAndEnergy + 6 * warp_size + lane, energy);
273         }
274         barrier(CLK_LOCAL_MEM_FENCE);
275
276         const int numIter = (c_virialAndEnergyCount + activeWarps - 1) / activeWarps;
277         for (int i = 0; i < numIter; i++)
278         {
279             const int componentIndex = i * activeWarps + warpIndex;
280             if (componentIndex < c_virialAndEnergyCount)
281             {
282                 const int targetIndex = componentIndex * warp_size + lane;
283 #pragma unroll
284                 for (int reductionStride = warp_size >> 1; reductionStride >= 1; reductionStride >>= 1)
285                 {
286                     if (lane < reductionStride)
287                     {
288                         sm_virialAndEnergy[targetIndex] +=
289                                 sm_virialAndEnergy[targetIndex + reductionStride];
290                     }
291 #ifdef _NVIDIA_SOURCE_
292                     /* FIXME: this execution happens within execution width aka warp, but somehow
293                      * NVIDIA OpenCL of all things fails without the memory barrier here. #2519
294                      */
295                     barrier(CLK_LOCAL_MEM_FENCE);
296 #endif
297                 }
298                 if (lane == 0)
299                 {
300                     atomicAdd_g_f(gm_virialAndEnergy + componentIndex, sm_virialAndEnergy[targetIndex]);
301                 }
302             }
303         }
304     }
305 }