Random cosmetic fixes
[alexxy/gromacs.git] / src / gromacs / ewald / pme_gpu_calculate_splines.cuh
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020,2021, 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 helper routines for PME gather and spline routines.
38  *
39  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
40  */
41
42 #include "gmxpre.h"
43
44 #include <cassert>
45
46 #include "gromacs/gpu_utils/cuda_kernel_utils.cuh"
47 #include "gromacs/gpu_utils/vectype_ops.cuh"
48
49 #include "pme.cuh"
50 #include "pme_grid.h"
51
52 /*! \internal \brief
53  * Gets a base of the unique index to an element in a spline parameter buffer (theta/dtheta),
54  * which is laid out for GPU spread/gather kernels. The base only corresponds to the atom index within the execution block.
55  * Feed the result into getSplineParamIndex() to get a full index.
56  * TODO: it's likely that both parameters can be just replaced with a single atom index, as they are derived from it.
57  * Do that, verifying that the generated code is not bloated, and/or revise the spline indexing scheme.
58  * Removing warp dependency would also be nice (and would probably coincide with removing c_pmeSpreadGatherAtomsPerWarp).
59  *
60  * \tparam order               PME order
61  * \tparam atomsPerWarp        Number of atoms processed by a warp
62  * \param[in] warpIndex        Warp index wrt the block.
63  * \param[in] atomWarpIndex    Atom index wrt the warp (from 0 to atomsPerWarp - 1).
64  *
65  * \returns Index into theta or dtheta array using GPU layout.
66  */
67 template<int order, int atomsPerWarp>
68 int __device__ __forceinline__ getSplineParamIndexBase(int warpIndex, int atomWarpIndex)
69 {
70     assert((atomWarpIndex >= 0) && (atomWarpIndex < atomsPerWarp));
71     const int dimIndex    = 0;
72     const int splineIndex = 0;
73     // The zeroes are here to preserve the full index formula for reference
74     return (((splineIndex + order * warpIndex) * DIM + dimIndex) * atomsPerWarp + atomWarpIndex);
75 }
76
77 /*! \internal \brief
78  * Gets a unique index to an element in a spline parameter buffer (theta/dtheta),
79  * which is laid out for GPU spread/gather kernels. The index is wrt to the execution block,
80  * in range(0, atomsPerBlock * order * DIM).
81  * This function consumes result of getSplineParamIndexBase() and adjusts it for \p dimIndex and \p splineIndex.
82  *
83  * \tparam order               PME order
84  * \tparam atomsPerWarp        Number of atoms processed by a warp
85  * \param[in] paramIndexBase   Must be result of getSplineParamIndexBase().
86  * \param[in] dimIndex         Dimension index (from 0 to 2)
87  * \param[in] splineIndex      Spline contribution index (from 0 to \p order - 1)
88  *
89  * \returns Index into theta or dtheta array using GPU layout.
90  */
91 template<int order, int atomsPerWarp>
92 int __device__ __forceinline__ getSplineParamIndex(int paramIndexBase, int dimIndex, int splineIndex)
93 {
94     assert((dimIndex >= XX) && (dimIndex < DIM));
95     assert((splineIndex >= 0) && (splineIndex < order));
96     return (paramIndexBase + (splineIndex * DIM + dimIndex) * atomsPerWarp);
97 }
98
99 /*! \internal \brief
100  * An inline CUDA function for skipping the zero-charge atoms.
101  *
102  * \returns                        Non-0 if atom should be processed, 0 otherwise.
103  * \param[in] coefficient          The atom charge.
104  *
105  * This is called from the spline_and_spread and gather PME kernels.
106  */
107 bool __device__ __forceinline__ pme_gpu_check_atom_charge(const float coefficient)
108 {
109     assert(isfinite(coefficient));
110     return c_skipNeutralAtoms ? (coefficient != 0.0F) : true;
111 }
112
113 //! Controls if the atom and charge data is prefeched into shared memory or loaded per thread from global
114 static const bool c_useAtomDataPrefetch = true;
115
116 /*! \brief Asserts if the argument is finite.
117  *
118  *  The function works for any data type, that can be casted to float. Note that there is also
119  *  a specialized implementation for float3 data type.
120  *
121  * \param[in] arg  Argument to check.
122  */
123 template<typename T>
124 __device__ inline void assertIsFinite(T arg);
125
126 template<>
127 __device__ inline void assertIsFinite(float3 gmx_unused arg)
128 {
129     assert(isfinite(static_cast<float>(arg.x)));
130     assert(isfinite(static_cast<float>(arg.y)));
131     assert(isfinite(static_cast<float>(arg.z)));
132 }
133
134 template<typename T>
135 __device__ inline void assertIsFinite(T gmx_unused arg)
136 {
137     assert(isfinite(static_cast<float>(arg)));
138 }
139
140 /*! \brief
141  * General purpose function for loading atom-related data from global to shared memory.
142  *
143  * \tparam[in] T                  Data type (float/int/...)
144  * \tparam[in] atomsPerBlock      Number of atoms processed by a block - should be accounted for in
145  * the size of the shared memory array.
146  * \tparam[in] dataCountPerAtom   Number of data elements per single atom (e.g. DIM for an rvec
147  * coordinates array).
148  * \param[out] sm_destination     Shared memory array for output.
149  * \param[in]  gm_source          Global memory array for input.
150  */
151 template<typename T, int atomsPerBlock, int dataCountPerAtom>
152 __device__ __forceinline__ void pme_gpu_stage_atom_data(T* __restrict__ sm_destination,
153                                                         const T* __restrict__ gm_source)
154 {
155     const int blockIndex       = blockIdx.y * gridDim.x + blockIdx.x;
156     const int threadLocalIndex = ((threadIdx.z * blockDim.y + threadIdx.y) * blockDim.x) + threadIdx.x;
157     const int localIndex       = threadLocalIndex;
158     const int globalIndexBase = blockIndex * atomsPerBlock * dataCountPerAtom;
159     const int globalIndex     = globalIndexBase + localIndex;
160     if (localIndex < atomsPerBlock * dataCountPerAtom)
161     {
162         assertIsFinite(gm_source[globalIndex]);
163         sm_destination[localIndex] = gm_source[globalIndex];
164     }
165 }
166
167 /*! \brief
168  * PME GPU spline parameter and gridline indices calculation.
169  * This corresponds to the CPU functions calc_interpolation_idx() and make_bsplines().
170  * First stage of the whole kernel.
171  *
172  * \tparam[in] order                PME interpolation order.
173  * \tparam[in] atomsPerBlock        Number of atoms processed by a block - should be accounted for
174  *                                  in the sizes of the shared memory arrays.
175  * \tparam[in] atomsPerWarp         Number of atoms processed by a warp
176  * \tparam[in] writeSmDtheta        Bool controlling if the theta derivative should be written to shared memory. Enables calculation of dtheta if set.
177  * \tparam[in] writeGlobal          A boolean which tells if the theta values and gridlines should be written to global memory. Enables calculation of dtheta if set.
178  * \param[in]  kernelParams         Input PME CUDA data in constant memory.
179  * \param[in]  atomIndexOffset      Starting atom index for the execution block w.r.t. global memory.
180  * \param[in]  atomX                Atom coordinate of atom processed by thread.
181  * \param[in]  atomCharge           Atom charge/coefficient of atom processed by thread.
182  * \param[out] sm_theta             Atom spline values in the shared memory.
183  * \param[out] sm_dtheta            Derivative of atom spline values in shared memory.
184  * \param[out] sm_gridlineIndices   Atom gridline indices in the shared memory.
185  */
186
187 template<int order, int atomsPerBlock, int atomsPerWarp, bool writeSmDtheta, bool writeGlobal>
188 __device__ __forceinline__ void calculate_splines(const PmeGpuCudaKernelParams kernelParams,
189                                                   const int                    atomIndexOffset,
190                                                   const float3                 atomX,
191                                                   const float                  atomCharge,
192                                                   float* __restrict__ sm_theta,
193                                                   float* __restrict__ sm_dtheta,
194                                                   int* __restrict__ sm_gridlineIndices)
195 {
196     /* Global memory pointers for output */
197     float* __restrict__ gm_theta         = kernelParams.atoms.d_theta;
198     float* __restrict__ gm_dtheta        = kernelParams.atoms.d_dtheta;
199     int* __restrict__ gm_gridlineIndices = kernelParams.atoms.d_gridlineIndices;
200
201     /* Fractional coordinates */
202     __shared__ float sm_fractCoords[atomsPerBlock * DIM];
203
204     /* Thread index w.r.t. block */
205     const int threadLocalId =
206             (threadIdx.z * (blockDim.x * blockDim.y)) + (threadIdx.y * blockDim.x) + threadIdx.x;
207     /* Warp index w.r.t. block - could probably be obtained easier? */
208     const int warpIndex = threadLocalId / warp_size;
209     /* Atom index w.r.t. warp - alternating 0 1 0 1 .. */
210     const int atomWarpIndex = threadIdx.z % atomsPerWarp;
211     /* Atom index w.r.t. block/shared memory */
212     const int atomIndexLocal = warpIndex * atomsPerWarp + atomWarpIndex;
213
214     /* Spline contribution index in one dimension */
215     const int threadLocalIdXY = (threadIdx.y * blockDim.x) + threadIdx.x;
216     const int orderIndex      = threadLocalIdXY / DIM;
217     /* Dimension index */
218     const int dimIndex = threadLocalIdXY % DIM;
219
220     /* Multi-purpose index of rvec/ivec atom data */
221     const int sharedMemoryIndex = atomIndexLocal * DIM + dimIndex;
222
223     float splineData[order];
224
225     const int localCheck = (dimIndex < DIM) && (orderIndex < 1);
226
227     /* we have 4 threads per atom, but can only use 3 here for the dimensions */
228     if (localCheck)
229     {
230         /* Indices interpolation */
231
232         if (orderIndex == 0)
233         {
234             int   tableIndex, tInt;
235             float n, t;
236             assert(atomIndexLocal < DIM * atomsPerBlock);
237             /* Accessing fields in fshOffset/nXYZ/recipbox/... with dimIndex offset
238              * puts them into local memory(!) instead of accessing the constant memory directly.
239              * That's the reason for the switch, to unroll explicitly.
240              * The commented parts correspond to the 0 components of the recipbox.
241              */
242             switch (dimIndex)
243             {
244                 case XX:
245                     tableIndex = kernelParams.grid.tablesOffsets[XX];
246                     n          = kernelParams.grid.realGridSizeFP[XX];
247                     t          = atomX.x * kernelParams.current.recipBox[dimIndex][XX]
248                         + atomX.y * kernelParams.current.recipBox[dimIndex][YY]
249                         + atomX.z * kernelParams.current.recipBox[dimIndex][ZZ];
250                     break;
251
252                 case YY:
253                     tableIndex = kernelParams.grid.tablesOffsets[YY];
254                     n          = kernelParams.grid.realGridSizeFP[YY];
255                     t = /*atomX.x * kernelParams.current.recipBox[dimIndex][XX] + */ atomX.y
256                                 * kernelParams.current.recipBox[dimIndex][YY]
257                         + atomX.z * kernelParams.current.recipBox[dimIndex][ZZ];
258                     break;
259
260                 case ZZ:
261                     tableIndex = kernelParams.grid.tablesOffsets[ZZ];
262                     n          = kernelParams.grid.realGridSizeFP[ZZ];
263                     t          = /*atomX.x * kernelParams.current.recipBox[dimIndex][XX] + atomX.y * kernelParams.current.recipBox[dimIndex][YY] + */ atomX
264                                 .z
265                         * kernelParams.current.recipBox[dimIndex][ZZ];
266                     break;
267             }
268             const float shift = c_pmeMaxUnitcellShift;
269             /* Fractional coordinates along box vectors, adding a positive shift to ensure t is positive for triclinic boxes */
270             t    = (t + shift) * n;
271             tInt = static_cast<int>(t);
272             assert(sharedMemoryIndex < atomsPerBlock * DIM);
273             sm_fractCoords[sharedMemoryIndex] = t - tInt;
274             tableIndex += tInt;
275             assert(tInt >= 0);
276             assert(tInt < c_pmeNeighborUnitcellCount * n);
277
278             // TODO have shared table for both parameters to share the fetch, as index is always same?
279             // TODO compare texture/LDG performance
280             sm_fractCoords[sharedMemoryIndex] += fetchFromParamLookupTable(
281                     kernelParams.grid.d_fractShiftsTable, kernelParams.fractShiftsTableTexture, tableIndex);
282             sm_gridlineIndices[sharedMemoryIndex] =
283                     fetchFromParamLookupTable(kernelParams.grid.d_gridlineIndicesTable,
284                                               kernelParams.gridlineIndicesTableTexture,
285                                               tableIndex);
286             if (writeGlobal)
287             {
288                 gm_gridlineIndices[atomIndexOffset * DIM + sharedMemoryIndex] =
289                         sm_gridlineIndices[sharedMemoryIndex];
290             }
291         }
292
293         /* B-spline calculation */
294
295         const int chargeCheck = pme_gpu_check_atom_charge(atomCharge);
296         if (chargeCheck)
297         {
298             float div;
299             int o = orderIndex; // This is an index that is set once for PME_GPU_PARALLEL_SPLINE == 1
300
301             const float dr = sm_fractCoords[sharedMemoryIndex];
302             assert(isfinite(dr));
303
304             /* dr is relative offset from lower cell limit */
305             splineData[order - 1] = 0.0F;
306             splineData[1]         = dr;
307             splineData[0]         = 1.0F - dr;
308
309 #pragma unroll
310             for (int k = 3; k < order; k++)
311             {
312                 div               = 1.0F / (k - 1.0F);
313                 splineData[k - 1] = div * dr * splineData[k - 2];
314 #pragma unroll
315                 for (int l = 1; l < (k - 1); l++)
316                 {
317                     splineData[k - l - 1] =
318                             div * ((dr + l) * splineData[k - l - 2] + (k - l - dr) * splineData[k - l - 1]);
319                 }
320                 splineData[0] = div * (1.0F - dr) * splineData[0];
321             }
322
323             const int thetaIndexBase =
324                     getSplineParamIndexBase<order, atomsPerWarp>(warpIndex, atomWarpIndex);
325             const int thetaGlobalOffsetBase = atomIndexOffset * DIM * order;
326             /* only calculate dtheta if we are saving it to shared or global memory */
327             if (writeSmDtheta || writeGlobal)
328             {
329                 /* Differentiation and storing the spline derivatives (dtheta) */
330 #pragma unroll
331                 for (o = 0; o < order; o++)
332                 {
333                     const int thetaIndex =
334                             getSplineParamIndex<order, atomsPerWarp>(thetaIndexBase, dimIndex, o);
335
336                     const float dtheta = ((o > 0) ? splineData[o - 1] : 0.0F) - splineData[o];
337                     assert(isfinite(dtheta));
338                     assert(thetaIndex < order * DIM * atomsPerBlock);
339                     if (writeSmDtheta)
340                     {
341                         sm_dtheta[thetaIndex] = dtheta;
342                     }
343                     if (writeGlobal)
344                     {
345                         const int thetaGlobalIndex  = thetaGlobalOffsetBase + thetaIndex;
346                         gm_dtheta[thetaGlobalIndex] = dtheta;
347                     }
348                 }
349             }
350
351             div                   = 1.0F / (order - 1.0F);
352             splineData[order - 1] = div * dr * splineData[order - 2];
353 #pragma unroll
354             for (int k = 1; k < (order - 1); k++)
355             {
356                 splineData[order - k - 1] = div
357                                             * ((dr + k) * splineData[order - k - 2]
358                                                + (order - k - dr) * splineData[order - k - 1]);
359             }
360             splineData[0] = div * (1.0F - dr) * splineData[0];
361
362             /* Storing the spline values (theta) */
363 #pragma unroll
364             for (o = 0; o < order; o++)
365             {
366                 const int thetaIndex =
367                         getSplineParamIndex<order, atomsPerWarp>(thetaIndexBase, dimIndex, o);
368                 assert(thetaIndex < order * DIM * atomsPerBlock);
369                 sm_theta[thetaIndex] = splineData[o];
370                 assert(isfinite(sm_theta[thetaIndex]));
371                 if (writeGlobal)
372                 {
373                     const int thetaGlobalIndex = thetaGlobalOffsetBase + thetaIndex;
374                     gm_theta[thetaGlobalIndex] = splineData[o];
375                 }
376             }
377         }
378     }
379 }