2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * \brief Implements PME GPU Fourier grid solving in CUDA.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
46 #include <math_constants.h>
48 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
53 * PME complex grid solver kernel function.
55 * \tparam[in] gridOrdering Specifies the dimension ordering of the complex grid.
56 * \tparam[in] computeEnergyAndVirial Tells if the reciprocal energy and virial should be computed.
57 * \param[in] kernelParams Input PME CUDA data in constant memory.
60 GridOrdering gridOrdering,
61 bool computeEnergyAndVirial
63 __launch_bounds__(c_solveMaxThreadsPerBlock)
64 __global__ void pme_solve_kernel(const struct PmeGpuCudaKernelParams kernelParams)
66 /* This kernel supports 2 different grid dimension orderings: YZX and XYZ */
67 int majorDim, middleDim, minorDim;
70 case GridOrdering::YZX:
76 case GridOrdering::XYZ:
86 /* Global memory pointers */
87 const float * __restrict__ gm_splineValueMajor = kernelParams.grid.d_splineModuli + kernelParams.grid.splineValuesOffset[majorDim];
88 const float * __restrict__ gm_splineValueMiddle = kernelParams.grid.d_splineModuli + kernelParams.grid.splineValuesOffset[middleDim];
89 const float * __restrict__ gm_splineValueMinor = kernelParams.grid.d_splineModuli + kernelParams.grid.splineValuesOffset[minorDim];
90 float * __restrict__ gm_virialAndEnergy = kernelParams.constants.d_virialAndEnergy;
91 float2 * __restrict__ gm_grid = (float2 *)kernelParams.grid.d_fourierGrid;
93 /* Various grid sizes and indices */
94 const int localOffsetMinor = 0, localOffsetMajor = 0, localOffsetMiddle = 0; //unused
95 const int localSizeMinor = kernelParams.grid.complexGridSizePadded[minorDim];
96 const int localSizeMiddle = kernelParams.grid.complexGridSizePadded[middleDim];
97 const int localCountMiddle = kernelParams.grid.complexGridSize[middleDim];
98 const int localCountMinor = kernelParams.grid.complexGridSize[minorDim];
99 const int nMajor = kernelParams.grid.realGridSize[majorDim];
100 const int nMiddle = kernelParams.grid.realGridSize[middleDim];
101 const int nMinor = kernelParams.grid.realGridSize[minorDim];
102 const int maxkMajor = (nMajor + 1) / 2; // X or Y
103 const int maxkMiddle = (nMiddle + 1) / 2; // Y OR Z => only check for !YZX
104 const int maxkMinor = (nMinor + 1) / 2; // Z or X => only check for YZX
106 /* Each thread works on one cell of the Fourier space complex 3D grid (gm_grid).
107 * Each block handles up to c_solveMaxThreadsPerBlock cells -
108 * depending on the grid contiguous dimension size,
109 * that can range from a part of a single gridline to several complete gridlines.
111 const int threadLocalId = threadIdx.x;
112 const int gridLineSize = localCountMinor;
113 const int gridLineIndex = threadLocalId / gridLineSize;
114 const int gridLineCellIndex = threadLocalId - gridLineSize * gridLineIndex;
115 const int gridLinesPerBlock = blockDim.x / gridLineSize;
116 const int activeWarps = (blockDim.x / warp_size);
117 const int indexMinor = blockIdx.x * blockDim.x + gridLineCellIndex;
118 const int indexMiddle = blockIdx.y * gridLinesPerBlock + gridLineIndex;
119 const int indexMajor = blockIdx.z;
121 /* Optional outputs */
130 assert(indexMajor < kernelParams.grid.complexGridSize[majorDim]);
131 if ((indexMiddle < localCountMiddle) & (indexMinor < localCountMinor) & (gridLineIndex < gridLinesPerBlock))
133 /* The offset should be equal to the global thread index for coalesced access */
134 const int gridIndex = (indexMajor * localSizeMiddle + indexMiddle) * localSizeMinor + indexMinor;
135 float2 * __restrict__ gm_gridCell = gm_grid + gridIndex;
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);
141 const int kMiddle = indexMiddle + localOffsetMiddle;
142 float mMiddle = kMiddle;
143 /* Checking Y in XYZ case */
144 if (gridOrdering == GridOrdering::XYZ)
146 mMiddle = (kMiddle < maxkMiddle) ? kMiddle : (kMiddle - nMiddle);
148 const int kMinor = localOffsetMinor + indexMinor;
149 float mMinor = kMinor;
150 /* Checking X in YZX case */
151 if (gridOrdering == GridOrdering::YZX)
153 mMinor = (kMinor < maxkMinor) ? kMinor : (kMinor - nMinor);
155 /* We should skip the k-space point (0,0,0) */
156 const bool notZeroPoint = (kMinor > 0) | (kMajor > 0) | (kMiddle > 0);
159 switch (gridOrdering)
161 case GridOrdering::YZX:
167 case GridOrdering::XYZ:
177 /* 0.5 correction factor for the first and last components of a Z dimension */
178 float corner_fac = 1.0f;
179 switch (gridOrdering)
181 case GridOrdering::YZX:
182 if ((kMiddle == 0) | (kMiddle == maxkMiddle))
188 case GridOrdering::XYZ:
189 if ((kMinor == 0) | (kMinor == maxkMinor))
201 const float mhxk = mX * kernelParams.current.recipBox[XX][XX];
202 const float mhyk = mX * kernelParams.current.recipBox[XX][YY] + mY * kernelParams.current.recipBox[YY][YY];
203 const float mhzk = mX * kernelParams.current.recipBox[XX][ZZ] + mY * kernelParams.current.recipBox[YY][ZZ] + mZ * kernelParams.current.recipBox[ZZ][ZZ];
205 const float m2k = mhxk * mhxk + mhyk * mhyk + mhzk * mhzk;
207 //TODO: use LDG/textures for gm_splineValue
208 float denom = m2k * float(CUDART_PI_F) * kernelParams.current.boxVolume * gm_splineValueMajor[kMajor] * gm_splineValueMiddle[kMiddle] * gm_splineValueMinor[kMinor];
209 assert(isfinite(denom));
210 assert(denom != 0.0f);
211 const float tmp1 = expf(-kernelParams.grid.ewaldFactor * m2k);
212 const float etermk = kernelParams.constants.elFactor * tmp1 / denom;
214 float2 gridValue = *gm_gridCell;
215 const float2 oldGridValue = gridValue;
216 gridValue.x *= etermk;
217 gridValue.y *= etermk;
218 *gm_gridCell = gridValue;
220 if (computeEnergyAndVirial)
222 const float tmp1k = 2.0f * (gridValue.x * oldGridValue.x + gridValue.y * oldGridValue.y);
224 float vfactor = (kernelParams.grid.ewaldFactor + 1.0f / m2k) * 2.0f;
225 float ets2 = corner_fac * tmp1k;
228 float ets2vf = ets2 * vfactor;
230 virxx = ets2vf * mhxk * mhxk - ets2;
231 virxy = ets2vf * mhxk * mhyk;
232 virxz = ets2vf * mhxk * mhzk;
233 viryy = ets2vf * mhyk * mhyk - ets2;
234 viryz = ets2vf * mhyk * mhzk;
235 virzz = ets2vf * mhzk * mhzk - ets2;
240 /* Optional energy/virial reduction */
241 if (computeEnergyAndVirial)
243 /* A tricky shuffle reduction inspired by reduce_force_j_warp_shfl.
244 * The idea is to reduce 7 energy/virial components into a single variable (aligned by 8).
245 * We will reduce everything into virxx.
248 /* We can only reduce warp-wise */
249 const int width = warp_size;
250 const unsigned int activeMask = c_fullWarpMask;
252 /* Making pair sums */
253 virxx += gmx_shfl_down_sync(activeMask, virxx, 1, width);
254 viryy += gmx_shfl_up_sync (activeMask, viryy, 1, width);
255 virzz += gmx_shfl_down_sync(activeMask, virzz, 1, width);
256 virxy += gmx_shfl_up_sync (activeMask, virxy, 1, width);
257 virxz += gmx_shfl_down_sync(activeMask, virxz, 1, width);
258 viryz += gmx_shfl_up_sync (activeMask, viryz, 1, width);
259 energy += gmx_shfl_down_sync(activeMask, energy, 1, width);
260 if (threadLocalId & 1)
262 virxx = viryy; // virxx now holds virxx and viryy pair sums
263 virzz = virxy; // virzz now holds virzz and virxy pair sums
264 virxz = viryz; // virxz now holds virxz and viryz pair sums
267 /* Making quad sums */
268 virxx += gmx_shfl_down_sync(activeMask, virxx, 2, width);
269 virzz += gmx_shfl_up_sync (activeMask, virzz, 2, width);
270 virxz += gmx_shfl_down_sync(activeMask, virxz, 2, width);
271 energy += gmx_shfl_up_sync (activeMask, energy, 2, width);
272 if (threadLocalId & 2)
274 virxx = virzz; // virxx now holds quad sums of virxx, virxy, virzz and virxy
275 virxz = energy; // virxz now holds quad sums of virxz, viryz, energy and unused paddings
278 /* Making octet sums */
279 virxx += gmx_shfl_down_sync(activeMask, virxx, 4, width);
280 virxz += gmx_shfl_up_sync (activeMask, virxz, 4, width);
281 if (threadLocalId & 4)
283 virxx = virxz; // virxx now holds all 7 components' octet sums + unused paddings
286 /* We only need to reduce virxx now */
288 for (int delta = 8; delta < width; delta <<= 1)
290 virxx += gmx_shfl_down_sync(activeMask, virxx, delta, width);
292 /* Now first 7 threads of each warp have the full output contributions in virxx */
294 const int componentIndex = threadLocalId & (warp_size - 1);
295 const bool validComponentIndex = (componentIndex < c_virialAndEnergyCount);
296 /* Reduce 7 outputs per warp in the shared memory */
297 const int stride = 8; // this is c_virialAndEnergyCount==7 rounded up to power of 2 for convenience, hence the assert
298 assert(c_virialAndEnergyCount == 7);
299 const int reductionBufferSize = (c_solveMaxThreadsPerBlock / warp_size) * stride;
300 __shared__ float sm_virialAndEnergy[reductionBufferSize];
302 if (validComponentIndex)
304 const int warpIndex = threadLocalId / warp_size;
305 sm_virialAndEnergy[warpIndex * stride + componentIndex] = virxx;
309 /* Reduce to the single warp size */
310 const int targetIndex = threadLocalId;
312 for (int reductionStride = reductionBufferSize >> 1; reductionStride >= warp_size; reductionStride >>= 1)
314 const int sourceIndex = targetIndex + reductionStride;
315 if ((targetIndex < reductionStride) & (sourceIndex < activeWarps * stride))
317 // TODO: the second conditional is only needed on first iteration, actually - see if compiler eliminates it!
318 sm_virialAndEnergy[targetIndex] += sm_virialAndEnergy[sourceIndex];
323 /* Now use shuffle again */
324 if (threadLocalId < warp_size)
326 float output = sm_virialAndEnergy[threadLocalId];
328 for (int delta = stride; delta < warp_size; delta <<= 1)
330 output += gmx_shfl_down_sync(activeMask, output, delta, warp_size);
333 if (validComponentIndex)
335 assert(isfinite(output));
336 atomicAdd(gm_virialAndEnergy + componentIndex, output);
342 //! Kernel instantiations
343 template __global__ void pme_solve_kernel<GridOrdering::YZX, true>(const PmeGpuCudaKernelParams);
344 template __global__ void pme_solve_kernel<GridOrdering::YZX, false>(const PmeGpuCudaKernelParams);
345 template __global__ void pme_solve_kernel<GridOrdering::XYZ, true>(const PmeGpuCudaKernelParams);
346 template __global__ void pme_solve_kernel<GridOrdering::XYZ, false>(const PmeGpuCudaKernelParams);