Avoid allocating SYCL buffer on each call to PME solve
[alexxy/gromacs.git] / src / gromacs / ewald / pme_gpu_program_impl.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \internal \file
36  * \brief
37  * Declares PmeGpuProgramImpl, which stores PME GPU (compiled) kernel handles.
38  *
39  * \author Aleksei Iupinov <a.yupinov@gmail.com>
40  * \ingroup module_ewald
41  */
42 #ifndef GMX_EWALD_PME_PME_GPU_PROGRAM_IMPL_H
43 #define GMX_EWALD_PME_PME_GPU_PROGRAM_IMPL_H
44
45 #include "config.h"
46
47 #include <memory>
48
49 #include "gromacs/gpu_utils/device_context.h"
50 #include "gromacs/utility/classhelpers.h"
51
52 class ISyclKernelFunctor;
53 class DeviceContext;
54 struct DeviceInformation;
55
56 /*! \internal
57  * \brief
58  * PME GPU persistent host program/kernel data, which should be initialized once for the whole execution.
59  *
60  * Primary purpose of this is to not recompile GPU kernels for each OpenCL unit test,
61  * while the relevant GPU context (e.g. cl_context) instance persists.
62  * In CUDA, this just assigns the kernel function pointers.
63  * This also implicitly relies on the fact that reasonable share of the kernels are always used.
64  * If there were more template parameters, even smaller share of all possible kernels would be used.
65  *
66  * \todo In future if we would need to react to either user input or
67  * auto-tuning to compile different kernels, then we might wish to
68  * revisit the number of kernels we pre-compile, and/or the management
69  * of their lifetime.
70  *
71  * This also doesn't manage cuFFT/clFFT kernels, which depend on the PME grid dimensions.
72  *
73  * TODO: pass cl_context to the constructor and not create it inside.
74  * See also Issue #2522.
75  */
76 struct PmeGpuProgramImpl
77 {
78     /*! \brief
79      * This is a handle to the GPU context, which is just a dummy in CUDA,
80      * but is created/destroyed by this class in OpenCL.
81      */
82     const DeviceContext& deviceContext_;
83
84     //! Conveniently all the PME kernels use the same single argument type
85 #if GMX_GPU_CUDA
86     using PmeKernelHandle = void (*)(const struct PmeGpuCudaKernelParams);
87 #elif GMX_GPU_OPENCL
88     using PmeKernelHandle = cl_kernel;
89 #else
90     using PmeKernelHandle = ISyclKernelFunctor*;
91 #endif
92
93     /*! \brief
94      * Maximum synchronous GPU thread group execution width.
95      * "Warp" is a CUDA term which we end up reusing in OpenCL kernels as well.
96      * For CUDA, this is a static value that comes from gromacs/gpu_utils/cuda_arch_utils.cuh;
97      * for OpenCL, we have to query it dynamically.
98      */
99     size_t warpSize_;
100
101     //@{
102     /**
103      * Spread/spline kernels are compiled only for order of 4.
104      * There are multiple versions of each kernel, paramaretized according to
105      *   Number of threads per atom. Using either order(4) or order*order (16) threads per atom is
106      * supported If the spline data is written in the spline/spread kernel and loaded in the gather
107      *   or recalculated in the gather.
108      * Spreading kernels also have hardcoded X/Y indices wrapping parameters,
109      * as a placeholder for implementing 1/2D decomposition.
110      * The kernels are templated separately for spreading on one grid (one or
111      * two sets of coefficients) or on two grids (required for energy and virial
112      * calculations).
113      */
114     size_t spreadWorkGroupSize;
115
116     PmeKernelHandle splineKernelSingle;
117     PmeKernelHandle splineKernelThPerAtom4Single;
118     PmeKernelHandle spreadKernelSingle;
119     PmeKernelHandle spreadKernelThPerAtom4Single;
120     PmeKernelHandle splineAndSpreadKernelSingle;
121     PmeKernelHandle splineAndSpreadKernelThPerAtom4Single;
122     PmeKernelHandle splineAndSpreadKernelWriteSplinesSingle;
123     PmeKernelHandle splineAndSpreadKernelWriteSplinesThPerAtom4Single;
124     PmeKernelHandle splineKernelDual;
125     PmeKernelHandle splineKernelThPerAtom4Dual;
126     PmeKernelHandle spreadKernelDual;
127     PmeKernelHandle spreadKernelThPerAtom4Dual;
128     PmeKernelHandle splineAndSpreadKernelDual;
129     PmeKernelHandle splineAndSpreadKernelThPerAtom4Dual;
130     PmeKernelHandle splineAndSpreadKernelWriteSplinesDual;
131     PmeKernelHandle splineAndSpreadKernelWriteSplinesThPerAtom4Dual;
132     //@}
133
134     //@{
135     /** Same for gather: hardcoded X/Y unwrap parameters, order of 4, plus
136      * it can either reduce with previous forces in the host buffer, or ignore them.
137      * Also similarly to the gather we can use either order(4) or order*order (16) threads per atom
138      * and either recalculate the splines or read the ones written by the spread
139      * The kernels are templated separately for using one or two grids (required for
140      * calculating energies and virial).
141      */
142     size_t gatherWorkGroupSize;
143
144     PmeKernelHandle gatherKernelSingle;
145     PmeKernelHandle gatherKernelThPerAtom4Single;
146     PmeKernelHandle gatherKernelReadSplinesSingle;
147     PmeKernelHandle gatherKernelReadSplinesThPerAtom4Single;
148     PmeKernelHandle gatherKernelDual;
149     PmeKernelHandle gatherKernelThPerAtom4Dual;
150     PmeKernelHandle gatherKernelReadSplinesDual;
151     PmeKernelHandle gatherKernelReadSplinesThPerAtom4Dual;
152     //@}
153
154     //@{
155     /** Solve kernel doesn't care about the interpolation order, but can optionally
156      * compute energy and virial, and supports XYZ and YZX grid orderings.
157      * The kernels are templated separately for grids in state A and B.
158      */
159     size_t solveMaxWorkGroupSize;
160
161     PmeKernelHandle solveYZXKernelA;
162     PmeKernelHandle solveXYZKernelA;
163     PmeKernelHandle solveYZXEnergyKernelA;
164     PmeKernelHandle solveXYZEnergyKernelA;
165     PmeKernelHandle solveYZXKernelB;
166     PmeKernelHandle solveXYZKernelB;
167     PmeKernelHandle solveYZXEnergyKernelB;
168     PmeKernelHandle solveXYZEnergyKernelB;
169     //@}
170
171     PmeGpuProgramImpl() = delete;
172     //! Constructor for the given device
173     explicit PmeGpuProgramImpl(const DeviceContext& deviceContext);
174     // NOLINTNEXTLINE(performance-trivially-destructible)
175     ~PmeGpuProgramImpl();
176     GMX_DISALLOW_COPY_AND_ASSIGN(PmeGpuProgramImpl);
177
178     //! Return the warp size for which the kernels were compiled
179     int warpSize() const { return warpSize_; }
180
181 private:
182     // Compiles kernels, if supported. Called by the constructor.
183     void compileKernels(const DeviceInformation& deviceInfo);
184 };
185
186 #endif