e576a5604503a969973f8448429284844e0e3436
[alexxy/gromacs.git] / src / gromacs / ewald / pme-gpu-utils.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018, 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 #ifndef GMX_EWALD_PME_GPU_UTILS_H
36 #define GMX_EWALD_PME_GPU_UTILS_H
37
38 /*! \internal \file
39  * \brief This file defines small PME GPU inline host/device functions.
40  * Note that OpenCL device-side functions can't use C++ features, so they are
41  * located in a similar file pme-gpu-utils.clh.
42  * Be sure to keep the logic in sync in both files when changing it!
43  *
44  * \author Aleksei Iupinov <a.yupinov@gmail.com>
45  * \ingroup module_ewald
46  */
47
48 #include "config.h"
49
50 #include <cassert>
51
52 #include "pme-gpu-constants.h"
53
54 //! A macro for inline GPU functions.
55 #if GMX_GPU == GMX_GPU_CUDA
56 #define INLINE_EVERYWHERE __host__ __device__ __forceinline__
57 #else
58 #define INLINE_EVERYWHERE inline
59 #endif
60
61 /*! \internal \brief
62  * Gets a base of the unique index to an element in a spline parameter buffer (theta/dtheta),
63  * which is laid out for GPU spread/gather kernels. The base only corresponds to the atom index within the execution block.
64  * Feed the result into getSplineParamIndex() to get a full index.
65  * TODO: it's likely that both parameters can be just replaced with a single atom index, as they are derived from it.
66  * Do that, verifying that the generated code is not bloated, and/or revise the spline indexing scheme.
67  * Removing warp dependency would also be nice (and would probably coincide with removing c_pmeSpreadGatherAtomsPerWarp).
68  *
69  * \tparam order               PME order
70  * \tparam atomsPerWarp        Number of atoms processed by a warp
71  * \param[in] warpIndex        Warp index wrt the block.
72  * \param[in] atomWarpIndex    Atom index wrt the warp (from 0 to atomsPerWarp - 1).
73  *
74  * \returns Index into theta or dtheta array using GPU layout.
75  */
76 template <int order, int atomsPerWarp>
77 int INLINE_EVERYWHERE getSplineParamIndexBase(int warpIndex, int atomWarpIndex)
78 {
79     assert((atomWarpIndex >= 0) && (atomWarpIndex < atomsPerWarp));
80     const int dimIndex    = 0;
81     const int splineIndex = 0;
82     // The zeroes are here to preserve the full index formula for reference
83     return (((splineIndex + order * warpIndex) * DIM + dimIndex) * atomsPerWarp + atomWarpIndex);
84 }
85
86 /*! \internal \brief
87  * Gets a unique index to an element in a spline parameter buffer (theta/dtheta),
88  * which is laid out for GPU spread/gather kernels. The index is wrt to the execution block,
89  * in range(0, atomsPerBlock * order * DIM).
90  * This function consumes result of getSplineParamIndexBase() and adjusts it for \p dimIndex and \p splineIndex.
91  *
92  * \tparam order               PME order
93  * \tparam atomsPerWarp        Number of atoms processed by a warp
94  * \param[in] paramIndexBase   Must be result of getSplineParamIndexBase().
95  * \param[in] dimIndex         Dimension index (from 0 to 2)
96  * \param[in] splineIndex      Spline contribution index (from 0 to \p order - 1)
97  *
98  * \returns Index into theta or dtheta array using GPU layout.
99  */
100 template <int order, int atomsPerWarp>
101 int INLINE_EVERYWHERE getSplineParamIndex(int paramIndexBase, int dimIndex, int splineIndex)
102 {
103     assert((dimIndex >= XX) && (dimIndex < DIM));
104     assert((splineIndex >= 0) && (splineIndex < order));
105     return (paramIndexBase + (splineIndex * DIM + dimIndex) * atomsPerWarp);
106 }
107
108 #if GMX_GPU == GMX_GPU_CUDA
109 // CUDA device code helpers below
110
111 /*! \internal \brief
112  * An inline CUDA function for checking the global atom data indices against the atom data array sizes.
113  *
114  * \param[in] atomDataIndex        The atom data index.
115  * \param[in] nAtomData            The atom data array element count.
116  * \returns                        Non-0 if index is within bounds (or PME data padding is enabled), 0 otherwise.
117  *
118  * This is called from the spline_and_spread and gather PME kernels.
119  * The goal is to isolate the global range checks, and allow avoiding them with c_usePadding enabled.
120  */
121 int __device__ __forceinline__ pme_gpu_check_atom_data_index(const int atomDataIndex, const int nAtomData)
122 {
123     return c_usePadding ? 1 : (atomDataIndex < nAtomData);
124 }
125
126 /*! \internal \brief
127  * An inline CUDA function for skipping the zero-charge atoms.
128  *
129  * \returns                        Non-0 if atom should be processed, 0 otherwise.
130  * \param[in] coefficient          The atom charge.
131  *
132  * This is called from the spline_and_spread and gather PME kernels.
133  */
134 int __device__ __forceinline__ pme_gpu_check_atom_charge(const float coefficient)
135 {
136     assert(isfinite(coefficient));
137     return c_skipNeutralAtoms ? (coefficient != 0.0f) : 1;
138 }
139
140 #endif
141
142 #endif