2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, 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.
35 #ifndef GMX_EWALD_PME_GPU_UTILS_CLH
36 #define GMX_EWALD_PME_GPU_UTILS_CLH
39 * \brief This file defines the small PME OpenCL inline device functions.
40 * This closely mirrors pme_gpu_utils.h (which is used in CUDA and unit tests), except with no templates.
41 * Instead of templated parameters this file expects following defines during compilation:
42 * - order - PME interpolation order;
43 * - atomsPerWarp - number of atoms processed by a warp (fixed for spread and gather kernels to be the same);
44 * - c_usePadding and c_skipNeutralAtoms - same as in pme_gpu_constants.h.
46 * \author Aleksei Iupinov <a.yupinov@gmail.com>
47 * \ingroup module_ewald
51 * Gets a base of the unique index to an element in a spline parameter buffer (theta/dtheta),
52 * which is laid out for GPU spread/gather kernels. The base only corresponds to the atom index within the execution block.
53 * Feed the result into getSplineParamIndex() to get a full index.
54 * TODO: it's likely that both parameters can be just replaced with a single atom index, as they are derived from it.
55 * Do that, verifying that the generated code is not bloated, and/or revise the spline indexing scheme.
56 * Removing warp dependency would also be nice (and would probably coincide with removing atomsPerWarp).
58 * \param[in] warpIndex Warp index wrt the block.
59 * \param[in] atomWarpIndex Atom index wrt the warp (from 0 to atomsPerWarp - 1).
61 * \returns Index into theta or dtheta array using GPU layout.
63 inline int getSplineParamIndexBase(int warpIndex, int atomWarpIndex)
65 assert((atomWarpIndex >= 0) && (atomWarpIndex < atomsPerWarp));
66 const int dimIndex = 0;
67 const int splineIndex = 0;
68 // The zeroes are here to preserve the full index formula for reference
69 return (((splineIndex + order * warpIndex) * DIM + dimIndex) * atomsPerWarp + atomWarpIndex);
73 * Gets a unique index to an element in a spline parameter buffer (theta/dtheta),
74 * which is laid out for GPU spread/gather kernels. The index is wrt to the execution block,
75 * in range(0, atomsPerBlock * order * DIM).
76 * This function consumes result of getSplineParamIndexBase() and adjusts it for \p dimIndex and \p splineIndex.
78 * \param[in] paramIndexBase Must be result of getSplineParamIndexBase().
79 * \param[in] dimIndex Dimension index (from 0 to 2)
80 * \param[in] splineIndex Spline contribution index (from 0 to \p order - 1)
81 * \param[in] order PME order
82 * \param[in] atomsPerWarp Number of atoms processed by a warp
84 * \returns Index into theta or dtheta array using GPU layout.
86 inline int getSplineParamIndex(int paramIndexBase, int dimIndex, int splineIndex)
88 assert((dimIndex >= XX) && (dimIndex < DIM));
89 assert((splineIndex >= 0) && (splineIndex < order));
90 return (paramIndexBase + (splineIndex * DIM + dimIndex) * atomsPerWarp);
94 * A function for checking the global atom data indices against the atom data array sizes.
96 * \param[in] atomDataIndexGlobal The atom data index.
97 * \param[in] nAtomData The atom data array element count.
98 * \returns Non-0 if index is within bounds (or PME data padding is enabled), 0 otherwise.
100 * This is called from the spline_and_spread and gather PME kernels.
101 * The goal is to isolate the global range checks, and allow avoiding them with c_usePadding being true.
103 inline int pme_gpu_check_atom_data_index(const size_t atomDataIndex, const size_t nAtomData)
105 return c_usePadding ? 1 : (atomDataIndex < nAtomData);
109 * A function for optionally skipping neutral charges, depending on c_skipNeutralAtoms.
111 * \param[in] coefficient The atom charge/coefficient.
112 * \returns Non-0 if atom should be processed, 0 otherwise.
114 inline int pme_gpu_check_atom_charge(const float coefficient)
116 return c_skipNeutralAtoms ? (coefficient != 0.0f) : 1;