2 * This file is part of the GROMACS molecular simulation package.
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.
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 force gathering in CUDA.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
46 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
47 #include "gromacs/gpu_utils/typecasts.cuh"
50 #include "pme_calculate_splines.cuh"
51 #include "pme_gpu_utils.h"
55 * An inline CUDA function: unroll the dynamic index accesses to the constant grid sizes to avoid local memory operations.
57 __device__ __forceinline__ float read_grid_size(const float* realGridSizeFP, const int dimIndex)
61 case XX: return realGridSizeFP[XX];
62 case YY: return realGridSizeFP[YY];
63 case ZZ: return realGridSizeFP[ZZ];
69 /*! \brief Reduce the partial force contributions.
71 * \tparam[in] order The PME order (must be 4).
72 * \tparam[in] atomDataSize The number of partial force contributions for each atom (currently
73 * order^2 == 16) \tparam[in] blockSize The CUDA block size \param[out] sm_forces Shared
74 * memory array with the output forces (number of elements is number of atoms per block) \param[in]
75 * atomIndexLocal Local atom index \param[in] splineIndex Spline index \param[in]
76 * lineIndex Line index (same as threadLocalId) \param[in] realGridSizeFP Local grid
77 * size constant \param[in] fx Input force partial component X \param[in] fy Input
78 * force partial component Y \param[in] fz Input force partial component Z
80 template<const int order, const int atomDataSize, const int blockSize>
81 __device__ __forceinline__ void reduce_atom_forces(float3* __restrict__ sm_forces,
82 const int atomIndexLocal,
83 const int splineIndex,
85 const float* realGridSizeFP,
90 if (!(order & (order - 1))) // Only for orders of power of 2
92 const unsigned int activeMask = c_fullWarpMask;
94 // A tricky shuffle reduction inspired by reduce_force_j_warp_shfl
95 // TODO: find out if this is the best in terms of transactions count
96 static_assert(order == 4, "Only order of 4 is implemented");
97 static_assert(atomDataSize <= warp_size,
98 "TODO: rework for atomDataSize > warp_size (order 8 or larger)");
99 const int width = atomDataSize;
101 fx += __shfl_down_sync(activeMask, fx, 1, width);
102 fy += __shfl_up_sync(activeMask, fy, 1, width);
103 fz += __shfl_down_sync(activeMask, fz, 1, width);
110 fx += __shfl_down_sync(activeMask, fx, 2, width);
111 fz += __shfl_up_sync(activeMask, fz, 2, width);
118 // By now fx contains intermediate quad sums of all 3 components:
119 // splineIndex 0 1 2 and 3 4 5 6 and 7 8...
120 // sum of... fx0 to fx3 fy0 to fy3 fz0 to fz3 fx4 to fx7 fy4 to fy7 fz4 to fz7 etc.
122 // We have to just further reduce those groups of 4
123 for (int delta = 4; delta < atomDataSize; delta <<= 1)
125 fx += __shfl_down_sync(activeMask, fx, delta, width);
128 const int dimIndex = splineIndex;
131 const float n = read_grid_size(realGridSizeFP, dimIndex);
132 *((float*)(&sm_forces[atomIndexLocal]) + dimIndex) = fx * n;
137 // We use blockSize shared memory elements to read fx, or fy, or fz, and then reduce them to
138 // fit into smemPerDim elements which are stored separately (first 2 dimensions only)
139 const int smemPerDim = warp_size;
140 const int smemReserved = (DIM)*smemPerDim;
141 __shared__ float sm_forceReduction[smemReserved + blockSize];
142 __shared__ float* sm_forceTemp[DIM];
144 const int numWarps = blockSize / smemPerDim;
145 const int minStride =
146 max(1, atomDataSize / numWarps); // order 4: 128 threads => 4, 256 threads => 2, etc
149 for (int dimIndex = 0; dimIndex < DIM; dimIndex++)
151 int elementIndex = smemReserved + lineIndex;
152 // Store input force contributions
153 sm_forceReduction[elementIndex] = (dimIndex == XX) ? fx : (dimIndex == YY) ? fy : fz;
154 // sync here because two warps write data that the first one consumes below
156 // Reduce to fit into smemPerDim (warp size)
158 for (int redStride = atomDataSize / 2; redStride > minStride; redStride >>= 1)
160 if (splineIndex < redStride)
162 sm_forceReduction[elementIndex] += sm_forceReduction[elementIndex + redStride];
166 // Last iteration - packing everything to be nearby, storing convenience pointer
167 sm_forceTemp[dimIndex] = sm_forceReduction + dimIndex * smemPerDim;
168 int redStride = minStride;
169 if (splineIndex < redStride)
171 const int packedIndex = atomIndexLocal * redStride + splineIndex;
172 sm_forceTemp[dimIndex][packedIndex] =
173 sm_forceReduction[elementIndex] + sm_forceReduction[elementIndex + redStride];
178 assert((blockSize / warp_size) >= DIM);
179 // assert (atomsPerBlock <= warp_size);
181 const int warpIndex = lineIndex / warp_size;
182 const int dimIndex = warpIndex;
184 // First 3 warps can now process 1 dimension each
187 int sourceIndex = lineIndex % warp_size;
189 for (int redStride = minStride / 2; redStride > 1; redStride >>= 1)
191 if (!(splineIndex & redStride))
193 sm_forceTemp[dimIndex][sourceIndex] += sm_forceTemp[dimIndex][sourceIndex + redStride];
199 const float n = read_grid_size(realGridSizeFP, dimIndex);
200 const int atomIndex = sourceIndex / minStride;
202 if (sourceIndex == minStride * atomIndex)
204 *((float*)(&sm_forces[atomIndex]) + dimIndex) =
205 (sm_forceTemp[dimIndex][sourceIndex] + sm_forceTemp[dimIndex][sourceIndex + 1]) * n;
212 * A CUDA kernel which gathers the atom forces from the grid.
213 * The grid is assumed to be wrapped in dimension Z.
215 * \tparam[in] order The PME order (must be 4 currently).
216 * \tparam[in] overwriteForces True: the forces are written to the output buffer;
217 * False: the forces are added non-atomically to the output buffer (e.g. to the bonded forces).
218 * \tparam[in] wrapX Tells if the grid is wrapped in the X dimension.
219 * \tparam[in] wrapY Tells if the grid is wrapped in the Y dimension.
220 * \tparam[in] readGlobal Tells if we should read spline values from global memory
221 * \tparam[in] useOrderThreads Tells if we should use order threads per atom (order*order used if false)
222 * \param[in] kernelParams All the PME GPU data.
224 template<const int order, const bool overwriteForces, const bool wrapX, const bool wrapY, const bool readGlobal, const bool useOrderThreads>
225 __launch_bounds__(c_gatherMaxThreadsPerBlock, c_gatherMinBlocksPerMP) __global__
226 void pme_gather_kernel(const PmeGpuCudaKernelParams kernelParams)
228 /* Global memory pointers */
229 const float* __restrict__ gm_coefficients = kernelParams.atoms.d_coefficients;
230 const float* __restrict__ gm_grid = kernelParams.grid.d_realGrid;
231 float* __restrict__ gm_forces = kernelParams.atoms.d_forces;
233 /* Global memory pointers for readGlobal */
234 const float* __restrict__ gm_theta = kernelParams.atoms.d_theta;
235 const float* __restrict__ gm_dtheta = kernelParams.atoms.d_dtheta;
236 const int* __restrict__ gm_gridlineIndices = kernelParams.atoms.d_gridlineIndices;
242 const int atomsPerBlock =
243 useOrderThreads ? (c_gatherMaxThreadsPerBlock / c_pmeSpreadGatherThreadsPerAtom4ThPerAtom)
244 : (c_gatherMaxThreadsPerBlock / c_pmeSpreadGatherThreadsPerAtom);
245 const int blockIndex = blockIdx.y * gridDim.x + blockIdx.x;
247 /* Number of data components and threads for a single atom */
248 const int atomDataSize = useOrderThreads ? c_pmeSpreadGatherThreadsPerAtom4ThPerAtom
249 : c_pmeSpreadGatherThreadsPerAtom;
250 const int atomsPerWarp = useOrderThreads ? c_pmeSpreadGatherAtomsPerWarp4ThPerAtom
251 : c_pmeSpreadGatherAtomsPerWarp;
253 const int blockSize = atomsPerBlock * atomDataSize;
254 assert(blockSize == blockDim.x * blockDim.y * blockDim.z);
256 /* These are the atom indices - for the shared and global memory */
257 const int atomIndexLocal = threadIdx.z;
258 const int atomIndexOffset = blockIndex * atomsPerBlock;
259 const int atomIndexGlobal = atomIndexOffset + atomIndexLocal;
261 /* Early return for fully empty blocks at the end
262 * (should only happen for billions of input atoms)
264 if (atomIndexOffset >= kernelParams.atoms.nAtoms)
268 // 4 warps per block, 8 atoms per warp *3 *4
269 const int splineParamsSize = atomsPerBlock * DIM * order;
270 const int gridlineIndicesSize = atomsPerBlock * DIM;
271 __shared__ int sm_gridlineIndices[gridlineIndicesSize];
272 __shared__ float sm_theta[splineParamsSize];
273 __shared__ float sm_dtheta[splineParamsSize];
275 /* Spline Z coordinates */
276 const int ithz = threadIdx.x;
278 /* These are the spline contribution indices in shared memory */
279 const int splineIndex = threadIdx.y * blockDim.x + threadIdx.x;
280 const int lineIndex = (threadIdx.z * (blockDim.x * blockDim.y))
281 + splineIndex; /* And to all the block's particles */
283 const int threadLocalId =
284 (threadIdx.z * (blockDim.x * blockDim.y)) + blockDim.x * threadIdx.y + threadIdx.x;
285 const int threadLocalIdMax = blockDim.x * blockDim.y * blockDim.z;
290 const int localGridlineIndicesIndex = threadLocalId;
291 const int globalGridlineIndicesIndex = blockIndex * gridlineIndicesSize + localGridlineIndicesIndex;
292 const int globalCheckIndices = pme_gpu_check_atom_data_index(
293 globalGridlineIndicesIndex, kernelParams.atoms.nAtoms * DIM);
294 if ((localGridlineIndicesIndex < gridlineIndicesSize) & globalCheckIndices)
296 sm_gridlineIndices[localGridlineIndicesIndex] = gm_gridlineIndices[globalGridlineIndicesIndex];
297 assert(sm_gridlineIndices[localGridlineIndicesIndex] >= 0);
299 /* The loop needed for order threads per atom to make sure we load all data values, as each thread must load multiple values
300 with order*order threads per atom, it is only required for each thread to load one data value */
303 const int iMax = useOrderThreads ? 3 : 1;
305 for (int i = iMin; i < iMax; i++)
307 int localSplineParamsIndex =
309 + i * threadLocalIdMax; /* i will always be zero for order*order threads per atom */
310 int globalSplineParamsIndex = blockIndex * splineParamsSize + localSplineParamsIndex;
311 int globalCheckSplineParams = pme_gpu_check_atom_data_index(
312 globalSplineParamsIndex, kernelParams.atoms.nAtoms * DIM * order);
313 if ((localSplineParamsIndex < splineParamsSize) && globalCheckSplineParams)
315 sm_theta[localSplineParamsIndex] = gm_theta[globalSplineParamsIndex];
316 sm_dtheta[localSplineParamsIndex] = gm_dtheta[globalSplineParamsIndex];
317 assert(isfinite(sm_theta[localSplineParamsIndex]));
318 assert(isfinite(sm_dtheta[localSplineParamsIndex]));
325 const float3* __restrict__ gm_coordinates = asFloat3(kernelParams.atoms.d_coordinates);
326 /* Recaclulate Splines */
327 if (c_useAtomDataPrefetch)
330 __shared__ float sm_coefficients[atomsPerBlock];
332 __shared__ float3 sm_coordinates[atomsPerBlock];
333 /* Staging coefficients/charges */
334 pme_gpu_stage_atom_data<float, atomsPerBlock, 1>(kernelParams, sm_coefficients, gm_coefficients);
336 /* Staging coordinates */
337 pme_gpu_stage_atom_data<float3, atomsPerBlock, 1>(kernelParams, sm_coordinates, gm_coordinates);
339 atomX = sm_coordinates[atomIndexLocal];
340 atomCharge = sm_coefficients[atomIndexLocal];
344 atomX = gm_coordinates[atomIndexGlobal];
345 atomCharge = gm_coefficients[atomIndexGlobal];
347 calculate_splines<order, atomsPerBlock, atomsPerWarp, true, false>(
348 kernelParams, atomIndexOffset, atomX, atomCharge, sm_theta, sm_dtheta, sm_gridlineIndices);
355 const int globalCheck = pme_gpu_check_atom_data_index(atomIndexGlobal, kernelParams.atoms.nAtoms);
356 const int chargeCheck = pme_gpu_check_atom_charge(gm_coefficients[atomIndexGlobal]);
358 if (chargeCheck & globalCheck)
360 const int nx = kernelParams.grid.realGridSize[XX];
361 const int ny = kernelParams.grid.realGridSize[YY];
362 const int nz = kernelParams.grid.realGridSize[ZZ];
363 const int pny = kernelParams.grid.realGridSizePadded[YY];
364 const int pnz = kernelParams.grid.realGridSizePadded[ZZ];
366 const int atomWarpIndex = atomIndexLocal % atomsPerWarp;
367 const int warpIndex = atomIndexLocal / atomsPerWarp;
369 const int splineIndexBase = getSplineParamIndexBase<order, atomsPerWarp>(warpIndex, atomWarpIndex);
370 const int splineIndexZ = getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, ZZ, ithz);
371 const float2 tdz = make_float2(sm_theta[splineIndexZ], sm_dtheta[splineIndexZ]);
373 int iz = sm_gridlineIndices[atomIndexLocal * DIM + ZZ] + ithz;
374 const int ixBase = sm_gridlineIndices[atomIndexLocal * DIM + XX];
382 const int ithyMin = useOrderThreads ? 0 : threadIdx.y;
383 const int ithyMax = useOrderThreads ? order : threadIdx.y + 1;
384 for (int ithy = ithyMin; ithy < ithyMax; ithy++)
386 const int splineIndexY = getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, YY, ithy);
387 const float2 tdy = make_float2(sm_theta[splineIndexY], sm_dtheta[splineIndexY]);
389 iy = sm_gridlineIndices[atomIndexLocal * DIM + YY] + ithy;
390 if (wrapY & (iy >= ny))
394 constOffset = iy * pnz + iz;
397 for (int ithx = 0; (ithx < order); ithx++)
399 int ix = ixBase + ithx;
400 if (wrapX & (ix >= nx))
404 const int gridIndexGlobal = ix * pny * pnz + constOffset;
405 assert(gridIndexGlobal >= 0);
406 const float gridValue = gm_grid[gridIndexGlobal];
407 assert(isfinite(gridValue));
408 const int splineIndexX =
409 getSplineParamIndex<order, atomsPerWarp>(splineIndexBase, XX, ithx);
410 const float2 tdx = make_float2(sm_theta[splineIndexX], sm_dtheta[splineIndexX]);
411 const float fxy1 = tdz.x * gridValue;
412 const float fz1 = tdz.y * gridValue;
413 fx += tdx.y * tdy.x * fxy1;
414 fy += tdx.x * tdy.y * fxy1;
415 fz += tdx.x * tdy.x * fz1;
420 // Reduction of partial force contributions
421 __shared__ float3 sm_forces[atomsPerBlock];
422 reduce_atom_forces<order, atomDataSize, blockSize>(sm_forces, atomIndexLocal, splineIndex, lineIndex,
423 kernelParams.grid.realGridSizeFP, fx, fy, fz);
426 /* Calculating the final forces with no component branching, atomsPerBlock threads */
427 const int forceIndexLocal = threadLocalId;
428 const int forceIndexGlobal = atomIndexOffset + forceIndexLocal;
429 const int calcIndexCheck = pme_gpu_check_atom_data_index(forceIndexGlobal, kernelParams.atoms.nAtoms);
430 if ((forceIndexLocal < atomsPerBlock) & calcIndexCheck)
432 const float3 atomForces = sm_forces[forceIndexLocal];
433 const float negCoefficient = -gm_coefficients[forceIndexGlobal];
435 result.x = negCoefficient * kernelParams.current.recipBox[XX][XX] * atomForces.x;
436 result.y = negCoefficient
437 * (kernelParams.current.recipBox[XX][YY] * atomForces.x
438 + kernelParams.current.recipBox[YY][YY] * atomForces.y);
439 result.z = negCoefficient
440 * (kernelParams.current.recipBox[XX][ZZ] * atomForces.x
441 + kernelParams.current.recipBox[YY][ZZ] * atomForces.y
442 + kernelParams.current.recipBox[ZZ][ZZ] * atomForces.z);
443 sm_forces[forceIndexLocal] = result;
447 assert(atomsPerBlock <= warp_size);
449 /* Writing or adding the final forces component-wise, single warp */
450 const int blockForcesSize = atomsPerBlock * DIM;
451 const int numIter = (blockForcesSize + warp_size - 1) / warp_size;
452 const int iterThreads = blockForcesSize / numIter;
453 if (threadLocalId < iterThreads)
456 for (int i = 0; i < numIter; i++)
458 int outputIndexLocal = i * iterThreads + threadLocalId;
459 int outputIndexGlobal = blockIndex * blockForcesSize + outputIndexLocal;
460 const int globalOutputCheck =
461 pme_gpu_check_atom_data_index(outputIndexGlobal, kernelParams.atoms.nAtoms * DIM);
462 if (globalOutputCheck)
464 const float outputForceComponent = ((float*)sm_forces)[outputIndexLocal];
467 gm_forces[outputIndexGlobal] = outputForceComponent;
471 gm_forces[outputIndexGlobal] += outputForceComponent;
478 //! Kernel instantiations
479 template __global__ void pme_gather_kernel<4, true, true, true, true, true>(const PmeGpuCudaKernelParams);
480 template __global__ void pme_gather_kernel<4, true, true, true, true, false>(const PmeGpuCudaKernelParams);
481 template __global__ void pme_gather_kernel<4, false, true, true, true, true>(const PmeGpuCudaKernelParams);
482 template __global__ void pme_gather_kernel<4, false, true, true, true, false>(const PmeGpuCudaKernelParams);
483 template __global__ void pme_gather_kernel<4, true, true, true, false, true>(const PmeGpuCudaKernelParams);
484 template __global__ void pme_gather_kernel<4, true, true, true, false, false>(const PmeGpuCudaKernelParams);
485 template __global__ void pme_gather_kernel<4, false, true, true, false, true>(const PmeGpuCudaKernelParams);
486 template __global__ void pme_gather_kernel<4, false, true, true, false, false>(const PmeGpuCudaKernelParams);