More const correctness
[gromacs.git] / src / gromacs / ewald / pme-spread.cu
blob65748941467d26aab305b9a6a5657a0790423091
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013-2016,2017,2018, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
38 /*! \internal \file
39  *  \brief Implements PME GPU spline calculation and charge spreading in CUDA.
40  *  TODO: consider always pre-sorting particles (as in DD case).
41  *
42  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
43  */
45 #include "gmxpre.h"
47 #include <cassert>
49 #include "gromacs/ewald/pme.h"
50 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
51 #include "gromacs/gpu_utils/cudautils.cuh"
52 #include "gromacs/utility/exceptions.h"
53 #include "gromacs/utility/gmxassert.h"
55 #include "pme.cuh"
56 #include "pme-grid.h"
57 #include "pme-timings.cuh"
60  * This define affects the spline calculation behaviour in the kernel.
61  * 0: a single GPU thread handles a single dimension of a single particle (calculating and storing (order) spline values and derivatives).
62  * 1: (order) threads do redundant work on this same task, each one stores only a single theta and single dtheta into global arrays.
63  * The only efficiency difference is less global store operations, countered by more redundant spline computation.
64  *
65  * TODO: estimate if this should be a boolean parameter (and add it to the unit test if so).
66  */
67 #define PME_GPU_PARALLEL_SPLINE 0
70 //! Spreading max block width in warps picked among powers of 2 (2, 4, 8, 16) for max. occupancy and min. runtime in most cases
71 constexpr int c_spreadMaxWarpsPerBlock = 8;
72 /* TODO: it has been observed that the kernel can be faster with smaller block sizes (2 or 4 warps)
73  * only on some GPUs (660Ti) with large enough grid (>= 48^3) due to load/store units being overloaded
74  * (ldst_fu_utilization metric maxed out in nvprof). Runtime block size choice might be nice to have.
75  * This has been tried on architectures up to Maxwell (GTX 750) and it would be good to revisit this.
76  */
77 //! Spreading max block size in threads
78 constexpr int c_spreadMaxThreadsPerBlock = c_spreadMaxWarpsPerBlock * warp_size;
81 /*! \brief
82  * General purpose function for loading atom-related data from global to shared memory.
83  *
84  * \tparam[in] T                 Data type (float/int/...)
85  * \tparam[in] atomsPerBlock     Number of atoms processed by a block - should be accounted for in the size of the shared memory array.
86  * \tparam[in] dataCountPerAtom  Number of data elements per single atom (e.g. DIM for an rvec coordinates array).
87  * \param[in]  kernelParams      Input PME CUDA data in constant memory.
88  * \param[out] sm_destination    Shared memory array for output.
89  * \param[in]  gm_source         Global memory array for input.
90  */
91 template<typename T,
92          const int atomsPerBlock,
93          const int dataCountPerAtom>
94 __device__  __forceinline__
95 void pme_gpu_stage_atom_data(const PmeGpuCudaKernelParams       kernelParams,
96                              T * __restrict__                   sm_destination,
97                              const T * __restrict__             gm_source)
99     static_assert(c_usePadding, "With padding disabled, index checking should be fixed to account for spline theta/dtheta per-warp alignment");
100     const int blockIndex       = blockIdx.y * gridDim.x + blockIdx.x;
101     const int threadLocalIndex = ((threadIdx.z * blockDim.y + threadIdx.y) * blockDim.x) + threadIdx.x;
102     const int localIndex       = threadLocalIndex;
103     const int globalIndexBase  = blockIndex * atomsPerBlock * dataCountPerAtom;
104     const int globalIndex      = globalIndexBase + localIndex;
105     const int globalCheck      = pme_gpu_check_atom_data_index(globalIndex, kernelParams.atoms.nAtoms * dataCountPerAtom);
106     if ((localIndex < atomsPerBlock * dataCountPerAtom) & globalCheck)
107     {
108         assert(isfinite(float(gm_source[globalIndex])));
109         sm_destination[localIndex] = gm_source[globalIndex];
110     }
113 /*! \brief
114  * PME GPU spline parameter and gridline indices calculation.
115  * This corresponds to the CPU functions calc_interpolation_idx() and make_bsplines().
116  * First stage of the whole kernel.
118  * \tparam[in] order                PME interpolation order.
119  * \tparam[in] atomsPerBlock        Number of atoms processed by a block - should be accounted for
120  *                                  in the sizes of the shared memory arrays.
121  * \param[in]  kernelParams         Input PME CUDA data in constant memory.
122  * \param[in]  atomIndexOffset      Starting atom index for the execution block w.r.t. global memory.
123  * \param[in]  sm_coordinates       Atom coordinates in the shared memory.
124  * \param[in]  sm_coefficients      Atom charges/coefficients in the shared memory.
125  * \param[out] sm_theta             Atom spline values in the shared memory.
126  * \param[out] sm_gridlineIndices   Atom gridline indices in the shared memory.
127  */
128 template <const int order,
129           const int atomsPerBlock>
130 __device__ __forceinline__ void calculate_splines(const PmeGpuCudaKernelParams           kernelParams,
131                                                   const int                              atomIndexOffset,
132                                                   const float3 * __restrict__            sm_coordinates,
133                                                   const float * __restrict__             sm_coefficients,
134                                                   float * __restrict__                   sm_theta,
135                                                   int * __restrict__                     sm_gridlineIndices)
137     /* Global memory pointers for output */
138     float * __restrict__ gm_theta           = kernelParams.atoms.d_theta;
139     float * __restrict__ gm_dtheta          = kernelParams.atoms.d_dtheta;
140     int * __restrict__   gm_gridlineIndices = kernelParams.atoms.d_gridlineIndices;
142     /* Fractional coordinates */
143     __shared__ float sm_fractCoords[atomsPerBlock * DIM];
145     /* Thread index w.r.t. block */
146     const int threadLocalId = (threadIdx.z * (blockDim.x * blockDim.y))
147         + (threadIdx.y * blockDim.x) + threadIdx.x;
148     /* Warp index w.r.t. block - could probably be obtained easier? */
149     const int warpIndex = threadLocalId / warp_size;
150     /* Thread index w.r.t. warp */
151     const int threadWarpIndex = threadLocalId % warp_size;
152     /* Atom index w.r.t. warp - alternating 0 1 0 1 .. */
153     const int atomWarpIndex = threadWarpIndex % PME_SPREADGATHER_ATOMS_PER_WARP;
154     /* Atom index w.r.t. block/shared memory */
155     const int atomIndexLocal = warpIndex * PME_SPREADGATHER_ATOMS_PER_WARP + atomWarpIndex;
157     /* Atom index w.r.t. global memory */
158     const int atomIndexGlobal = atomIndexOffset + atomIndexLocal;
159     /* Spline contribution index in one dimension */
160     const int orderIndex = threadWarpIndex / (PME_SPREADGATHER_ATOMS_PER_WARP * DIM);
161     /* Dimension index */
162     const int dimIndex = (threadWarpIndex - orderIndex * (PME_SPREADGATHER_ATOMS_PER_WARP * DIM)) / PME_SPREADGATHER_ATOMS_PER_WARP;
164     /* Multi-purpose index of rvec/ivec atom data */
165     const int sharedMemoryIndex = atomIndexLocal * DIM + dimIndex;
167     /* Spline parameters need a working buffer.
168      * With PME_GPU_PARALLEL_SPLINE == 0 it is just a local array of (order) values for some of the threads, which is fine;
169      * With PME_GPU_PARALLEL_SPLINE == 1 (order) times more threads are involved, so the shared memory is used to avoid overhead.
170      * The buffer's size, striding and indexing are adjusted accordingly.
171      * The buffer is accessed with SPLINE_DATA_PTR and SPLINE_DATA macros.
172      */
173 #if PME_GPU_PARALLEL_SPLINE
174     const int        splineDataStride  = atomsPerBlock * DIM;
175     const int        splineDataIndex   = sharedMemoryIndex;
176     __shared__ float sm_splineData[splineDataStride * order];
177     float           *splineDataPtr = sm_splineData;
178 #else
179     const int        splineDataStride = 1;
180     const int        splineDataIndex  = 0;
181     float            splineData[splineDataStride * order];
182     float           *splineDataPtr = splineData;
183 #endif
185 #define SPLINE_DATA_PTR(i) (splineDataPtr + ((i) * splineDataStride + splineDataIndex))
186 #define SPLINE_DATA(i) (*SPLINE_DATA_PTR(i))
188     const int localCheck  = (dimIndex < DIM) && (orderIndex < (PME_GPU_PARALLEL_SPLINE ? order : 1));
189     const int globalCheck = pme_gpu_check_atom_data_index(atomIndexGlobal, kernelParams.atoms.nAtoms);
191     if (localCheck && globalCheck)
192     {
193         /* Indices interpolation */
195         if (orderIndex == 0)
196         {
197             int           tableIndex, tInt;
198             float         n, t;
199             const float3  x = sm_coordinates[atomIndexLocal];
200             /* Accessing fields in fshOffset/nXYZ/recipbox/... with dimIndex offset
201              * puts them into local memory(!) instead of accessing the constant memory directly.
202              * That's the reason for the switch, to unroll explicitly.
203              * The commented parts correspond to the 0 components of the recipbox.
204              */
205             switch (dimIndex)
206             {
207                 case XX:
208                     tableIndex  = kernelParams.grid.tablesOffsets[XX];
209                     n           = kernelParams.grid.realGridSizeFP[XX];
210                     t           = x.x * kernelParams.current.recipBox[dimIndex][XX] + x.y * kernelParams.current.recipBox[dimIndex][YY] + x.z * kernelParams.current.recipBox[dimIndex][ZZ];
211                     break;
213                 case YY:
214                     tableIndex  = kernelParams.grid.tablesOffsets[YY];
215                     n           = kernelParams.grid.realGridSizeFP[YY];
216                     t           = /*x.x * kernelParams.current.recipbox[dimIndex][XX] + */ x.y * kernelParams.current.recipBox[dimIndex][YY] + x.z * kernelParams.current.recipBox[dimIndex][ZZ];
217                     break;
219                 case ZZ:
220                     tableIndex  = kernelParams.grid.tablesOffsets[ZZ];
221                     n           = kernelParams.grid.realGridSizeFP[ZZ];
222                     t           = /*x.x * kernelParams.current.recipbox[dimIndex][XX] + x.y * kernelParams.current.recipbox[dimIndex][YY] + */ x.z * kernelParams.current.recipBox[dimIndex][ZZ];
223                     break;
224             }
225             const float shift = c_pmeMaxUnitcellShift;
226             /* Fractional coordinates along box vectors, adding a positive shift to ensure t is positive for triclinic boxes */
227             t    = (t + shift) * n;
228             tInt = (int)t;
229             sm_fractCoords[sharedMemoryIndex] = t - tInt;
230             tableIndex                       += tInt;
231             assert(tInt >= 0);
232             assert(tInt < c_pmeNeighborUnitcellCount * n);
234             // TODO have shared table for both parameters to share the fetch, as index is always same?
235             // TODO compare texture/LDG performance
236             sm_fractCoords[sharedMemoryIndex] +=
237                 fetchFromParamLookupTable(kernelParams.grid.d_fractShiftsTable,
238                                           kernelParams.fractShiftsTableTexture,
239                                           tableIndex);
240             sm_gridlineIndices[sharedMemoryIndex] =
241                 fetchFromParamLookupTable(kernelParams.grid.d_gridlineIndicesTable,
242                                           kernelParams.gridlineIndicesTableTexture,
243                                           tableIndex);
244             gm_gridlineIndices[atomIndexOffset * DIM + sharedMemoryIndex] = sm_gridlineIndices[sharedMemoryIndex];
245         }
247         /* B-spline calculation */
249         const int chargeCheck = pme_gpu_check_atom_charge(sm_coefficients[atomIndexLocal]);
250         if (chargeCheck)
251         {
252             float       div;
253             int         o = orderIndex; // This is an index that is set once for PME_GPU_PARALLEL_SPLINE == 1
255             const float dr = sm_fractCoords[sharedMemoryIndex];
256             assert(isfinite(dr));
258             /* dr is relative offset from lower cell limit */
259             *SPLINE_DATA_PTR(order - 1) = 0.0f;
260             *SPLINE_DATA_PTR(1)         = dr;
261             *SPLINE_DATA_PTR(0)         = 1.0f - dr;
263 #pragma unroll
264             for (int k = 3; k < order; k++)
265             {
266                 div                     = 1.0f / (k - 1.0f);
267                 *SPLINE_DATA_PTR(k - 1) = div * dr * SPLINE_DATA(k - 2);
268 #pragma unroll
269                 for (int l = 1; l < (k - 1); l++)
270                 {
271                     *SPLINE_DATA_PTR(k - l - 1) = div * ((dr + l) * SPLINE_DATA(k - l - 2) + (k - l - dr) * SPLINE_DATA(k - l - 1));
272                 }
273                 *SPLINE_DATA_PTR(0) = div * (1.0f - dr) * SPLINE_DATA(0);
274             }
276             const int thetaIndexBase        = getSplineParamIndexBase<order>(warpIndex, atomWarpIndex);
277             const int thetaGlobalOffsetBase = atomIndexOffset * DIM * order;
279             /* Differentiation and storing the spline derivatives (dtheta) */
280 #if !PME_GPU_PARALLEL_SPLINE
281             // With PME_GPU_PARALLEL_SPLINE == 1, o is already set to orderIndex;
282             // With PME_GPU_PARALLEL_SPLINE == 0, we loop o over range(order).
283 #pragma unroll
284             for (o = 0; o < order; o++)
285 #endif
286             {
287                 const int   thetaIndex       = getSplineParamIndex<order>(thetaIndexBase, dimIndex, o);
288                 const int   thetaGlobalIndex = thetaGlobalOffsetBase + thetaIndex;
290                 const float dtheta = ((o > 0) ? SPLINE_DATA(o - 1) : 0.0f) - SPLINE_DATA(o);
291                 assert(isfinite(dtheta));
292                 gm_dtheta[thetaGlobalIndex] = dtheta;
293             }
295             div  = 1.0f / (order - 1.0f);
296             *SPLINE_DATA_PTR(order - 1) = div * dr * SPLINE_DATA(order - 2);
297 #pragma unroll
298             for (int k = 1; k < (order - 1); k++)
299             {
300                 *SPLINE_DATA_PTR(order - k - 1) = div * ((dr + k) * SPLINE_DATA(order - k - 2) + (order - k - dr) * SPLINE_DATA(order - k - 1));
301             }
302             *SPLINE_DATA_PTR(0) = div * (1.0f - dr) * SPLINE_DATA(0);
304             /* Storing the spline values (theta) */
305 #if !PME_GPU_PARALLEL_SPLINE
306             // See comment for the loop above
307 #pragma unroll
308             for (o = 0; o < order; o++)
309 #endif
310             {
311                 const int thetaIndex       = getSplineParamIndex<order>(thetaIndexBase, dimIndex, o);
312                 const int thetaGlobalIndex = thetaGlobalOffsetBase + thetaIndex;
314                 sm_theta[thetaIndex]       = SPLINE_DATA(o);
315                 assert(isfinite(sm_theta[thetaIndex]));
316                 gm_theta[thetaGlobalIndex] = SPLINE_DATA(o);
317             }
318         }
319     }
322 /*! \brief
323  * Charge spreading onto the grid.
324  * This corresponds to the CPU function spread_coefficients_bsplines_thread().
325  * Second stage of the whole kernel.
327  * \tparam[in] order                PME interpolation order.
328  * \tparam[in] wrapX                A boolean which tells if the grid overlap in dimension X should be wrapped.
329  * \tparam[in] wrapY                A boolean which tells if the grid overlap in dimension Y should be wrapped.
330  * \param[in]  kernelParams         Input PME CUDA data in constant memory.
331  * \param[in]  atomIndexOffset      Starting atom index for the execution block w.r.t. global memory.
332  * \param[in]  sm_coefficients      Atom coefficents/charges in the shared memory.
333  * \param[in]  sm_gridlineIndices   Atom gridline indices in the shared memory.
334  * \param[in]  sm_theta             Atom spline values in the shared memory.
335  */
336 template <
337     const int order, const bool wrapX, const bool wrapY>
338 __device__ __forceinline__ void spread_charges(const PmeGpuCudaKernelParams           kernelParams,
339                                                int                                    atomIndexOffset,
340                                                const float * __restrict__             sm_coefficients,
341                                                const int * __restrict__               sm_gridlineIndices,
342                                                const float * __restrict__             sm_theta)
344     /* Global memory pointer to the output grid */
345     float * __restrict__ gm_grid = kernelParams.grid.d_realGrid;
348     const int nx  = kernelParams.grid.realGridSize[XX];
349     const int ny  = kernelParams.grid.realGridSize[YY];
350     const int nz  = kernelParams.grid.realGridSize[ZZ];
351     const int pny = kernelParams.grid.realGridSizePadded[YY];
352     const int pnz = kernelParams.grid.realGridSizePadded[ZZ];
354     const int offx = 0, offy = 0, offz = 0; // unused for now
356     const int atomIndexLocal  = threadIdx.z;
357     const int atomIndexGlobal = atomIndexOffset + atomIndexLocal;
359     const int globalCheck = pme_gpu_check_atom_data_index(atomIndexGlobal, kernelParams.atoms.nAtoms);
360     const int chargeCheck = pme_gpu_check_atom_charge(sm_coefficients[atomIndexLocal]);
361     if (chargeCheck & globalCheck)
362     {
363         // Spline Y/Z coordinates
364         const int ithy   = threadIdx.y;
365         const int ithz   = threadIdx.x;
366         const int ixBase = sm_gridlineIndices[atomIndexLocal * DIM + XX] - offx;
367         int       iy     = sm_gridlineIndices[atomIndexLocal * DIM + YY] - offy + ithy;
368         if (wrapY & (iy >= ny))
369         {
370             iy -= ny;
371         }
372         int iz  = sm_gridlineIndices[atomIndexLocal * DIM + ZZ] - offz + ithz;
373         if (iz >= nz)
374         {
375             iz -= nz;
376         }
378         /* Atom index w.r.t. warp - alternating 0 1 0 1 .. */
379         const int    atomWarpIndex   = atomIndexLocal % PME_SPREADGATHER_ATOMS_PER_WARP;
380         /* Warp index w.r.t. block - could probably be obtained easier? */
381         const int    warpIndex       = atomIndexLocal / PME_SPREADGATHER_ATOMS_PER_WARP;
383         const int    splineIndexBase = getSplineParamIndexBase<order>(warpIndex, atomWarpIndex);
384         const int    splineIndexZ    = getSplineParamIndex<order>(splineIndexBase, ZZ, ithz);
385         const float  thetaZ          = sm_theta[splineIndexZ];
386         const int    splineIndexY    = getSplineParamIndex<order>(splineIndexBase, YY, ithy);
387         const float  thetaY          = sm_theta[splineIndexY];
388         const float  constVal        = thetaZ * thetaY * sm_coefficients[atomIndexLocal];
389         assert(isfinite(constVal));
390         const int    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             const int   splineIndexX    = getSplineParamIndex<order>(splineIndexBase, XX, ithx);
402             const float thetaX          = sm_theta[splineIndexX];
403             assert(isfinite(thetaX));
404             assert(isfinite(gm_grid[gridIndexGlobal]));
405             atomicAdd(gm_grid + gridIndexGlobal, thetaX * constVal);
406         }
407     }
410 /*! \brief
411  * A spline computation and charge spreading kernel function.
413  * \tparam[in] order                PME interpolation order.
414  * \tparam[in] computeSplines       A boolean which tells if the spline parameter and
415  *                                  gridline indices' computation should be performed.
416  * \tparam[in] spreadCharges        A boolean which tells if the charge spreading should be performed.
417  * \tparam[in] wrapX                A boolean which tells if the grid overlap in dimension X should be wrapped.
418  * \tparam[in] wrapY                A boolean which tells if the grid overlap in dimension Y should be wrapped.
419  * \param[in]  kernelParams         Input PME CUDA data in constant memory.
420  */
421 template <
422     const int order,
423     const bool computeSplines,
424     const bool spreadCharges,
425     const bool wrapX,
426     const bool wrapY
427     >
428 __launch_bounds__(c_spreadMaxThreadsPerBlock)
429 __global__ void pme_spline_and_spread_kernel(const PmeGpuCudaKernelParams kernelParams)
431     const int        atomsPerBlock = c_spreadMaxThreadsPerBlock / PME_SPREADGATHER_THREADS_PER_ATOM;
432     // Gridline indices, ivec
433     __shared__ int   sm_gridlineIndices[atomsPerBlock * DIM];
434     // Charges
435     __shared__ float sm_coefficients[atomsPerBlock];
436     // Spline values
437     __shared__ float sm_theta[atomsPerBlock * DIM * order];
439     const int        blockIndex      = blockIdx.y * gridDim.x + blockIdx.x;
440     const int        atomIndexOffset = blockIndex * atomsPerBlock;
442     /* Early return for fully empty blocks at the end
443      * (should only happen on Fermi or billions of input atoms)
444      */
445     if (atomIndexOffset >= kernelParams.atoms.nAtoms)
446     {
447         return;
448     }
450     /* Staging coefficients/charges for both spline and spread */
451     pme_gpu_stage_atom_data<float, atomsPerBlock, 1>(kernelParams, sm_coefficients, kernelParams.atoms.d_coefficients);
453     if (computeSplines)
454     {
455         /* Staging coordinates */
456         __shared__ float sm_coordinates[DIM * atomsPerBlock];
457         pme_gpu_stage_atom_data<float, atomsPerBlock, DIM>(kernelParams, sm_coordinates, kernelParams.atoms.d_coordinates);
459         __syncthreads();
460         calculate_splines<order, atomsPerBlock>(kernelParams, atomIndexOffset, (const float3 *)sm_coordinates,
461                                                 sm_coefficients, sm_theta, sm_gridlineIndices);
462         gmx_syncwarp();
463     }
464     else
465     {
466         /* Staging the data for spread
467          * (the data is assumed to be in GPU global memory with proper layout already,
468          * as in after running the spline kernel)
469          */
470         /* Spline data - only thetas (dthetas will only be needed in gather) */
471         pme_gpu_stage_atom_data<float, atomsPerBlock, DIM * order>(kernelParams, sm_theta, kernelParams.atoms.d_theta);
472         /* Gridline indices */
473         pme_gpu_stage_atom_data<int, atomsPerBlock, DIM>(kernelParams, sm_gridlineIndices, kernelParams.atoms.d_gridlineIndices);
475         __syncthreads();
476     }
478     /* Spreading */
479     if (spreadCharges)
480     {
481         spread_charges<order, wrapX, wrapY>(kernelParams, atomIndexOffset, sm_coefficients, sm_gridlineIndices, sm_theta);
482     }
485 void pme_gpu_spread(const PmeGpu    *pmeGpu,
486                     int gmx_unused   gridIndex,
487                     real            *h_grid,
488                     bool             computeSplines,
489                     bool             spreadCharges)
491     GMX_ASSERT(computeSplines || spreadCharges, "PME spline/spread kernel has invalid input (nothing to do)");
492     cudaStream_t  stream          = pmeGpu->archSpecific->pmeStream;
493     const auto   *kernelParamsPtr = pmeGpu->kernelParams.get();
494     GMX_ASSERT(kernelParamsPtr->atoms.nAtoms > 0, "No atom data in PME GPU spread");
496     const int order         = pmeGpu->common->pme_order;
497     const int atomsPerBlock = c_spreadMaxThreadsPerBlock / PME_SPREADGATHER_THREADS_PER_ATOM;
498     // TODO: pick smaller block size in runtime if needed
499     // (e.g. on 660 Ti where 50% occupancy is ~25% faster than 100% occupancy with RNAse (~17.8k atoms))
500     // If doing so, change atomsPerBlock in the kernels as well.
501     // TODO: test varying block sizes on modern arch-s as well
502     // TODO: also consider using cudaFuncSetCacheConfig() for preferring shared memory on older architectures
503     //(for spline data mostly, together with varying PME_GPU_PARALLEL_SPLINE define)
504     GMX_ASSERT(!c_usePadding || !(PME_ATOM_DATA_ALIGNMENT % atomsPerBlock), "inconsistent atom data padding vs. spreading block size");
506     const int blockCount = pmeGpu->nAtomsPadded / atomsPerBlock;
507     auto      dimGrid    = pmeGpuCreateGrid(pmeGpu, blockCount);
508     dim3 dimBlock(order, order, atomsPerBlock);
510     // These should later check for PME decomposition
511     const bool wrapX = true;
512     const bool wrapY = true;
513     GMX_UNUSED_VALUE(wrapX);
514     GMX_UNUSED_VALUE(wrapY);
515     switch (order)
516     {
517         case 4:
518         {
519             // TODO: cleaner unroll with some template trick?
520             if (computeSplines)
521             {
522                 if (spreadCharges)
523                 {
524                     pme_gpu_start_timing(pmeGpu, gtPME_SPLINEANDSPREAD);
525                     pme_spline_and_spread_kernel<4, true, true, wrapX, wrapY> <<< dimGrid, dimBlock, 0, stream>>> (*kernelParamsPtr);
526                     CU_LAUNCH_ERR("pme_spline_and_spread_kernel");
527                     pme_gpu_stop_timing(pmeGpu, gtPME_SPLINEANDSPREAD);
528                 }
529                 else
530                 {
531                     pme_gpu_start_timing(pmeGpu, gtPME_SPLINE);
532                     pme_spline_and_spread_kernel<4, true, false, wrapX, wrapY> <<< dimGrid, dimBlock, 0, stream>>> (*kernelParamsPtr);
533                     CU_LAUNCH_ERR("pme_spline_and_spread_kernel");
534                     pme_gpu_stop_timing(pmeGpu, gtPME_SPLINE);
535                 }
536             }
537             else
538             {
539                 pme_gpu_start_timing(pmeGpu, gtPME_SPREAD);
540                 pme_spline_and_spread_kernel<4, false, true, wrapX, wrapY> <<< dimGrid, dimBlock, 0, stream>>> (*kernelParamsPtr);
541                 CU_LAUNCH_ERR("pme_spline_and_spread_kernel");
542                 pme_gpu_stop_timing(pmeGpu, gtPME_SPREAD);
543             }
544         }
545         break;
547         default:
548             GMX_THROW(gmx::NotImplementedError("The code for pme_order != 4 was not tested!"));
549     }
551     const bool copyBackGrid = spreadCharges && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_FFT(pmeGpu));
552     if (copyBackGrid)
553     {
554         pme_gpu_copy_output_spread_grid(pmeGpu, h_grid);
555     }
556     const bool copyBackAtomData = computeSplines && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_gather(pmeGpu));
557     if (copyBackAtomData)
558     {
559         pme_gpu_copy_output_spread_atom_data(pmeGpu);
560     }