Simplify PME GPU constants
[gromacs.git] / src / gromacs / ewald / pme_gather.cu
blob2edec547ef7d17e39365cf5d0ce07f99202b1173
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020, 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  */
36 /*! \internal \file
37  *  \brief Implements PME force gathering in CUDA.
38  *
39  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
40  */
42 #include "gmxpre.h"
44 #include <cassert>
46 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
47 #include "gromacs/gpu_utils/typecasts.cuh"
49 #include "pme.cuh"
50 #include "pme_gpu_calculate_splines.cuh"
51 #include "pme_grid.h"
53 /*! \brief
54  * An inline CUDA function: unroll the dynamic index accesses to the constant grid sizes to avoid local memory operations.
55  */
56 __device__ __forceinline__ float read_grid_size(const float* realGridSizeFP, const int dimIndex)
58     switch (dimIndex)
59     {
60         case XX: return realGridSizeFP[XX];
61         case YY: return realGridSizeFP[YY];
62         case ZZ: return realGridSizeFP[ZZ];
63     }
64     assert(false);
65     return 0.0f;
68 /*! \brief Reduce the partial force contributions.
69  *
70  * \tparam[in] order              The PME order (must be 4).
71  * \tparam[in] atomDataSize       The number of partial force contributions for each atom (currently
72  *                                order^2 == 16)
73  * \tparam[in] blockSize          The CUDA block size
74  *
75  * \param[out] sm_forces          Shared memory array with the output forces (number of elements
76  *                                is number of atoms per block)
77  * \param[in]  atomIndexLocal     Local atom index
78  * \param[in]  splineIndex        Spline index
79  * \param[in]  lineIndex          Line index (same as threadLocalId)
80  * \param[in]  realGridSizeFP     Local grid size constant
81  * \param[in]  fx                 Input force partial component X
82  * \param[in]  fy                 Input force partial component Y
83  * \param[in]  fz                 Input force partial component Z
84  */
85 template<int order, int atomDataSize, int blockSize>
86 __device__ __forceinline__ void reduce_atom_forces(float3* __restrict__ sm_forces,
87                                                    const int    atomIndexLocal,
88                                                    const int    splineIndex,
89                                                    const int    lineIndex,
90                                                    const float* realGridSizeFP,
91                                                    float&       fx,
92                                                    float&       fy,
93                                                    float&       fz)
95     if (!(order & (order - 1))) // Only for orders of power of 2
96     {
97         const unsigned int activeMask = c_fullWarpMask;
99         // A tricky shuffle reduction inspired by reduce_force_j_warp_shfl
100         // TODO: find out if this is the best in terms of transactions count
101         static_assert(order == 4, "Only order of 4 is implemented");
102         static_assert(atomDataSize <= warp_size,
103                       "TODO: rework for atomDataSize > warp_size (order 8 or larger)");
104         const int width = atomDataSize;
106         fx += __shfl_down_sync(activeMask, fx, 1, width);
107         fy += __shfl_up_sync(activeMask, fy, 1, width);
108         fz += __shfl_down_sync(activeMask, fz, 1, width);
110         if (splineIndex & 1)
111         {
112             fx = fy;
113         }
115         fx += __shfl_down_sync(activeMask, fx, 2, width);
116         fz += __shfl_up_sync(activeMask, fz, 2, width);
118         if (splineIndex & 2)
119         {
120             fx = fz;
121         }
123         // By now fx contains intermediate quad sums of all 3 components:
124         // splineIndex    0            1            2 and 3      4            5            6 and 7 8...
125         // sum of...      fx0 to fx3   fy0 to fy3   fz0 to fz3   fx4 to fx7   fy4 to fy7   fz4 to fz7 etc.
127         // We have to just further reduce those groups of 4
128         for (int delta = 4; delta < atomDataSize; delta <<= 1)
129         {
130             fx += __shfl_down_sync(activeMask, fx, delta, width);
131         }
133         const int dimIndex = splineIndex;
134         if (dimIndex < DIM)
135         {
136             const float n = read_grid_size(realGridSizeFP, dimIndex);
137             *((float*)(&sm_forces[atomIndexLocal]) + dimIndex) = fx * n;
138         }
139     }
140     else
141     {
142         // We use blockSize shared memory elements to read fx, or fy, or fz, and then reduce them to
143         // fit into smemPerDim elements which are stored separately (first 2 dimensions only)
144         const int         smemPerDim   = warp_size;
145         const int         smemReserved = (DIM)*smemPerDim;
146         __shared__ float  sm_forceReduction[smemReserved + blockSize];
147         __shared__ float* sm_forceTemp[DIM];
149         const int numWarps = blockSize / smemPerDim;
150         const int minStride =
151                 max(1, atomDataSize / numWarps); // order 4: 128 threads => 4, 256 threads => 2, etc
153 #pragma unroll
154         for (int dimIndex = 0; dimIndex < DIM; dimIndex++)
155         {
156             int elementIndex = smemReserved + lineIndex;
157             // Store input force contributions
158             sm_forceReduction[elementIndex] = (dimIndex == XX) ? fx : (dimIndex == YY) ? fy : fz;
159             // sync here because two warps write data that the first one consumes below
160             __syncthreads();
161             // Reduce to fit into smemPerDim (warp size)
162 #pragma unroll
163             for (int redStride = atomDataSize / 2; redStride > minStride; redStride >>= 1)
164             {
165                 if (splineIndex < redStride)
166                 {
167                     sm_forceReduction[elementIndex] += sm_forceReduction[elementIndex + redStride];
168                 }
169             }
170             __syncthreads();
171             // Last iteration - packing everything to be nearby, storing convenience pointer
172             sm_forceTemp[dimIndex] = sm_forceReduction + dimIndex * smemPerDim;
173             int redStride          = minStride;
174             if (splineIndex < redStride)
175             {
176                 const int packedIndex = atomIndexLocal * redStride + splineIndex;
177                 sm_forceTemp[dimIndex][packedIndex] =
178                         sm_forceReduction[elementIndex] + sm_forceReduction[elementIndex + redStride];
179             }
180             __syncthreads();
181         }
183         assert((blockSize / warp_size) >= DIM);
184         // assert (atomsPerBlock <= warp_size);
186         const int warpIndex = lineIndex / warp_size;
187         const int dimIndex  = warpIndex;
189         // First 3 warps can now process 1 dimension each
190         if (dimIndex < DIM)
191         {
192             int sourceIndex = lineIndex % warp_size;
193 #pragma unroll
194             for (int redStride = minStride / 2; redStride > 1; redStride >>= 1)
195             {
196                 if (!(splineIndex & redStride))
197                 {
198                     sm_forceTemp[dimIndex][sourceIndex] += sm_forceTemp[dimIndex][sourceIndex + redStride];
199                 }
200             }
202             __syncwarp();
204             const float n         = read_grid_size(realGridSizeFP, dimIndex);
205             const int   atomIndex = sourceIndex / minStride;
207             if (sourceIndex == minStride * atomIndex)
208             {
209                 *((float*)(&sm_forces[atomIndex]) + dimIndex) =
210                         (sm_forceTemp[dimIndex][sourceIndex] + sm_forceTemp[dimIndex][sourceIndex + 1]) * n;
211             }
212         }
213     }
216 /*! \brief
217  * A CUDA kernel which gathers the atom forces from the grid.
218  * The grid is assumed to be wrapped in dimension Z.
220  * \tparam[in] order                The PME order (must be 4 currently).
221  * \tparam[in] wrapX                Tells if the grid is wrapped in the X dimension.
222  * \tparam[in] wrapY                Tells if the grid is wrapped in the Y dimension.
223  * \tparam[in] readGlobal           Tells if we should read spline values from global memory
224  * \tparam[in] threadsPerAtom       How many threads work on each atom
226  * \param[in]  kernelParams         All the PME GPU data.
227  */
228 template<int order, bool wrapX, bool wrapY, bool readGlobal, ThreadsPerAtom threadsPerAtom>
229 __launch_bounds__(c_gatherMaxThreadsPerBlock, c_gatherMinBlocksPerMP) __global__
230         void pme_gather_kernel(const PmeGpuCudaKernelParams kernelParams)
232     /* Global memory pointers */
233     const float* __restrict__ gm_coefficients = kernelParams.atoms.d_coefficients;
234     const float* __restrict__ gm_grid         = kernelParams.grid.d_realGrid;
235     float* __restrict__ gm_forces             = kernelParams.atoms.d_forces;
237     /* Global memory pointers for readGlobal */
238     const float* __restrict__ gm_theta         = kernelParams.atoms.d_theta;
239     const float* __restrict__ gm_dtheta        = kernelParams.atoms.d_dtheta;
240     const int* __restrict__ gm_gridlineIndices = kernelParams.atoms.d_gridlineIndices;
242     float3 atomX;
243     float  atomCharge;
245     const int blockIndex = blockIdx.y * gridDim.x + blockIdx.x;
247     /* Number of data components and threads for a single atom */
248     const int threadsPerAtomValue = (threadsPerAtom == ThreadsPerAtom::Order) ? order : order * order;
249     const int atomDataSize        = threadsPerAtomValue;
250     const int atomsPerBlock       = c_gatherMaxThreadsPerBlock / atomDataSize;
251     // Number of atoms processed by a single warp in spread and gather
252     const int atomsPerWarp = warp_size / atomDataSize;
254     const int blockSize = atomsPerBlock * atomDataSize;
255     assert(blockSize == blockDim.x * blockDim.y * blockDim.z);
257     /* These are the atom indices - for the shared and global memory */
258     const int atomIndexLocal  = threadIdx.z;
259     const int atomIndexOffset = blockIndex * atomsPerBlock;
260     const int atomIndexGlobal = atomIndexOffset + atomIndexLocal;
262     /* Early return for fully empty blocks at the end
263      * (should only happen for billions of input atoms)
264      */
265     if (atomIndexOffset >= kernelParams.atoms.nAtoms)
266     {
267         return;
268     }
269     // 4 warps per block, 8 atoms per warp *3 *4
270     const int        splineParamsSize    = atomsPerBlock * DIM * order;
271     const int        gridlineIndicesSize = atomsPerBlock * DIM;
272     __shared__ int   sm_gridlineIndices[gridlineIndicesSize];
273     __shared__ float sm_theta[splineParamsSize];
274     __shared__ float sm_dtheta[splineParamsSize];
276     /* Spline Z coordinates */
277     const int ithz = threadIdx.x;
279     /* These are the spline contribution indices in shared memory */
280     const int splineIndex = threadIdx.y * blockDim.x + threadIdx.x;
281     const int lineIndex   = (threadIdx.z * (blockDim.x * blockDim.y))
282                           + splineIndex; /* And to all the block's particles */
284     const int threadLocalId =
285             (threadIdx.z * (blockDim.x * blockDim.y)) + blockDim.x * threadIdx.y + threadIdx.x;
286     const int threadLocalIdMax = blockDim.x * blockDim.y * blockDim.z;
288     if (readGlobal)
289     {
290         /* Read splines */
291         const int localGridlineIndicesIndex = threadLocalId;
292         const int globalGridlineIndicesIndex = blockIndex * gridlineIndicesSize + localGridlineIndicesIndex;
293         if (localGridlineIndicesIndex < gridlineIndicesSize)
294         {
295             sm_gridlineIndices[localGridlineIndicesIndex] = gm_gridlineIndices[globalGridlineIndicesIndex];
296             assert(sm_gridlineIndices[localGridlineIndicesIndex] >= 0);
297         }
298         /* The loop needed for order threads per atom to make sure we load all data values, as each thread must load multiple values
299            with order*order threads per atom, it is only required for each thread to load one data value */
301         const int iMin = 0;
302         const int iMax = (threadsPerAtom == ThreadsPerAtom::Order) ? 3 : 1;
304         for (int i = iMin; i < iMax; i++)
305         {
306             int localSplineParamsIndex =
307                     threadLocalId
308                     + i * threadLocalIdMax; /* i will always be zero for order*order threads per atom */
309             int globalSplineParamsIndex = blockIndex * splineParamsSize + localSplineParamsIndex;
310             if (localSplineParamsIndex < splineParamsSize)
311             {
312                 sm_theta[localSplineParamsIndex]  = gm_theta[globalSplineParamsIndex];
313                 sm_dtheta[localSplineParamsIndex] = gm_dtheta[globalSplineParamsIndex];
314                 assert(isfinite(sm_theta[localSplineParamsIndex]));
315                 assert(isfinite(sm_dtheta[localSplineParamsIndex]));
316             }
317         }
318         __syncthreads();
319     }
320     else
321     {
322         const float3* __restrict__ gm_coordinates = asFloat3(kernelParams.atoms.d_coordinates);
323         /* Recaclulate  Splines  */
324         if (c_useAtomDataPrefetch)
325         {
326             // charges
327             __shared__ float sm_coefficients[atomsPerBlock];
328             // Coordinates
329             __shared__ float3 sm_coordinates[atomsPerBlock];
330             /* Staging coefficients/charges */
331             pme_gpu_stage_atom_data<float, atomsPerBlock, 1>(sm_coefficients, gm_coefficients);
333             /* Staging coordinates */
334             pme_gpu_stage_atom_data<float3, atomsPerBlock, 1>(sm_coordinates, gm_coordinates);
335             __syncthreads();
336             atomX      = sm_coordinates[atomIndexLocal];
337             atomCharge = sm_coefficients[atomIndexLocal];
338         }
339         else
340         {
341             atomX      = gm_coordinates[atomIndexGlobal];
342             atomCharge = gm_coefficients[atomIndexGlobal];
343         }
344         calculate_splines<order, atomsPerBlock, atomsPerWarp, true, false>(
345                 kernelParams, atomIndexOffset, atomX, atomCharge, sm_theta, sm_dtheta, sm_gridlineIndices);
346         __syncwarp();
347     }
348     float fx = 0.0f;
349     float fy = 0.0f;
350     float fz = 0.0f;
352     const int chargeCheck = pme_gpu_check_atom_charge(gm_coefficients[atomIndexGlobal]);
354     if (chargeCheck)
355     {
356         const int nx  = kernelParams.grid.realGridSize[XX];
357         const int ny  = kernelParams.grid.realGridSize[YY];
358         const int nz  = kernelParams.grid.realGridSize[ZZ];
359         const int pny = kernelParams.grid.realGridSizePadded[YY];
360         const int pnz = kernelParams.grid.realGridSizePadded[ZZ];
362         const int atomWarpIndex = atomIndexLocal % atomsPerWarp;
363         const int warpIndex     = atomIndexLocal / atomsPerWarp;
365         const int splineIndexBase = getSplineParamIndexBase<order, atomsPerWarp>(warpIndex, atomWarpIndex);
366         const int splineIndexZ = getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, ZZ, ithz);
367         const float2 tdz       = make_float2(sm_theta[splineIndexZ], sm_dtheta[splineIndexZ]);
369         int       iz     = sm_gridlineIndices[atomIndexLocal * DIM + ZZ] + ithz;
370         const int ixBase = sm_gridlineIndices[atomIndexLocal * DIM + XX];
372         if (iz >= nz)
373         {
374             iz -= nz;
375         }
376         int constOffset, iy;
378         const int ithyMin = (threadsPerAtom == ThreadsPerAtom::Order) ? 0 : threadIdx.y;
379         const int ithyMax = (threadsPerAtom == ThreadsPerAtom::Order) ? order : threadIdx.y + 1;
380         for (int ithy = ithyMin; ithy < ithyMax; ithy++)
381         {
382             const int splineIndexY = getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, YY, ithy);
383             const float2 tdy       = make_float2(sm_theta[splineIndexY], sm_dtheta[splineIndexY]);
385             iy = sm_gridlineIndices[atomIndexLocal * DIM + YY] + ithy;
386             if (wrapY & (iy >= ny))
387             {
388                 iy -= ny;
389             }
390             constOffset = iy * pnz + iz;
392 #pragma unroll
393             for (int ithx = 0; (ithx < order); ithx++)
394             {
395                 int ix = ixBase + ithx;
396                 if (wrapX & (ix >= nx))
397                 {
398                     ix -= nx;
399                 }
400                 const int gridIndexGlobal = ix * pny * pnz + constOffset;
401                 assert(gridIndexGlobal >= 0);
402                 const float gridValue = gm_grid[gridIndexGlobal];
403                 assert(isfinite(gridValue));
404                 const int splineIndexX =
405                         getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, XX, ithx);
406                 const float2 tdx  = make_float2(sm_theta[splineIndexX], sm_dtheta[splineIndexX]);
407                 const float  fxy1 = tdz.x * gridValue;
408                 const float  fz1  = tdz.y * gridValue;
409                 fx += tdx.y * tdy.x * fxy1;
410                 fy += tdx.x * tdy.y * fxy1;
411                 fz += tdx.x * tdy.x * fz1;
412             }
413         }
414     }
416     // Reduction of partial force contributions
417     __shared__ float3 sm_forces[atomsPerBlock];
418     reduce_atom_forces<order, atomDataSize, blockSize>(sm_forces, atomIndexLocal, splineIndex, lineIndex,
419                                                        kernelParams.grid.realGridSizeFP, fx, fy, fz);
420     __syncthreads();
422     /* Calculating the final forces with no component branching, atomsPerBlock threads */
423     const int forceIndexLocal  = threadLocalId;
424     const int forceIndexGlobal = atomIndexOffset + forceIndexLocal;
425     if (forceIndexLocal < atomsPerBlock)
426     {
427         const float3 atomForces     = sm_forces[forceIndexLocal];
428         const float  negCoefficient = -gm_coefficients[forceIndexGlobal];
429         float3       result;
430         result.x = negCoefficient * kernelParams.current.recipBox[XX][XX] * atomForces.x;
431         result.y = negCoefficient
432                    * (kernelParams.current.recipBox[XX][YY] * atomForces.x
433                       + kernelParams.current.recipBox[YY][YY] * atomForces.y);
434         result.z = negCoefficient
435                    * (kernelParams.current.recipBox[XX][ZZ] * atomForces.x
436                       + kernelParams.current.recipBox[YY][ZZ] * atomForces.y
437                       + kernelParams.current.recipBox[ZZ][ZZ] * atomForces.z);
438         sm_forces[forceIndexLocal] = result;
439     }
441     __syncwarp();
442     assert(atomsPerBlock <= warp_size);
444     /* Writing or adding the final forces component-wise, single warp */
445     const int blockForcesSize = atomsPerBlock * DIM;
446     const int numIter         = (blockForcesSize + warp_size - 1) / warp_size;
447     const int iterThreads     = blockForcesSize / numIter;
448     if (threadLocalId < iterThreads)
449     {
450 #pragma unroll
451         for (int i = 0; i < numIter; i++)
452         {
453             int         outputIndexLocal     = i * iterThreads + threadLocalId;
454             int         outputIndexGlobal    = blockIndex * blockForcesSize + outputIndexLocal;
455             const float outputForceComponent = ((float*)sm_forces)[outputIndexLocal];
456             gm_forces[outputIndexGlobal]     = outputForceComponent;
457         }
458     }
461 //! Kernel instantiations
462 // clang-format off
463 template __global__ void pme_gather_kernel<4, true, true, true,  ThreadsPerAtom::Order>       (const PmeGpuCudaKernelParams);
464 template __global__ void pme_gather_kernel<4, true, true, true,  ThreadsPerAtom::OrderSquared>(const PmeGpuCudaKernelParams);
465 template __global__ void pme_gather_kernel<4, true, true, false, ThreadsPerAtom::Order>       (const PmeGpuCudaKernelParams);
466 template __global__ void pme_gather_kernel<4, true, true, false, ThreadsPerAtom::OrderSquared>(const PmeGpuCudaKernelParams);
467 // clang-format on