Move GPU 3D FFT code to fft module
[alexxy/gromacs.git] / src / gromacs / ewald / pme_gpu_types_host_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
36 /*! \internal \file
37  * \brief Defines the host-side PME GPU data structure, which is dependent on the GPU types.
38  * It's included by pointer in the general PmeGpu host structure in pme_gpu_types_host.h.
39  *
40  * \author Aleksei Iupinov <a.yupinov@gmail.com>
41  * \ingroup module_ewald
42  */
43
44 #ifndef PMEGPUTYPESHOSTIMPL_H
45 #define PMEGPUTYPESHOSTIMPL_H
46
47 #include "config.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49
50 #include <array>
51 #include <set>
52 #include <vector>
53
54 #if GMX_GPU_CUDA
55 #    include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
56 #    include "gromacs/gpu_utils/gpuregiontimer.cuh"
57 #elif GMX_GPU_OPENCL
58 #    include "gromacs/gpu_utils/gpueventsynchronizer_ocl.h"
59 #    include "gromacs/gpu_utils/gpuregiontimer_ocl.h"
60 #elif GMX_GPU_SYCL
61 #    include "gromacs/gpu_utils/gpueventsynchronizer_sycl.h"
62 #    include "gromacs/gpu_utils/gpuregiontimer_sycl.h"
63 #endif
64
65 #include "gromacs/fft/gpu_3dfft.h"
66 #include "gromacs/timing/gpu_timing.h" // for gtPME_EVENT_COUNT
67
68 #ifndef NUMFEPSTATES
69 //! Number of FEP states.
70 #    define NUMFEPSTATES 2
71 #endif
72
73 namespace gmx
74 {
75 class Gpu3dFft;
76 } // namespace gmx
77
78 /*! \internal \brief
79  * The main PME CUDA/OpenCL-specific host data structure, included in the PME GPU structure by the archSpecific pointer.
80  */
81 struct PmeGpuSpecific
82 {
83     /*! \brief Constructor
84      *
85      * \param[in] deviceContext  GPU device context
86      * \param[in] pmeStream      GPU pme stream.
87      */
88     PmeGpuSpecific(const DeviceContext& deviceContext, const DeviceStream& pmeStream) :
89         deviceContext_(deviceContext), pmeStream_(pmeStream)
90     {
91     }
92
93     /*! \brief
94      * A handle to the GPU context.
95      * TODO: this is currently extracted from the implementation of pmeGpu->programHandle_,
96      * but should be a constructor parameter to PmeGpu, as well as PmeGpuProgram,
97      * managed by high-level code.
98      */
99     const DeviceContext& deviceContext_;
100
101     /*! \brief The GPU stream where everything related to the PME happens. */
102     const DeviceStream& pmeStream_;
103
104     /* Synchronization events */
105     /*! \brief Triggered after the PME Force Calculations have been completed */
106     GpuEventSynchronizer pmeForcesReady;
107     /*! \brief Triggered after the grid has been copied to the host (after the spreading stage). */
108     GpuEventSynchronizer syncSpreadGridD2H;
109
110     /* Settings which are set at the start of the run */
111     /*! \brief A boolean which tells whether the complex and real grids for cu/clFFT are different or same. Currenty true. */
112     bool performOutOfPlaceFFT = false;
113     /*! \brief A boolean which tells if the GPU timing events are enabled.
114      *  False by default, can be enabled by setting the environment variable GMX_ENABLE_GPU_TIMING.
115      *  Note: will not be reliable when multiple GPU tasks are running concurrently on the same
116      * device context, as CUDA events on multiple streams are untrustworthy.
117      */
118     bool useTiming = false;
119
120     //! Vector of FFT setups
121     std::vector<std::unique_ptr<gmx::Gpu3dFft>> fftSetup;
122
123     //! All the timers one might use
124     gmx::EnumerationArray<PmeStage, GpuRegionTimer> timingEvents;
125
126     //! Indices of timingEvents actually used
127     std::set<PmeStage> activeTimers;
128
129     /* GPU arrays element counts (not the arrays sizes in bytes!).
130      * They might be larger than the actual meaningful data sizes.
131      * These are paired: the actual element count + the maximum element count that can fit in the current allocated memory.
132      * These integer pairs are mostly meaningful for the reallocateDeviceBuffer calls.
133      * As such, if DeviceBuffer is refactored into a class, they can be freely changed, too.
134      * The only exceptions are realGridSize and complexGridSize which are also used for grid clearing/copying.
135      * TODO: these should live in a clean buffered container type, and be refactored in the NB/cudautils as well.
136      */
137     /*! \brief The kernelParams.atoms.coordinates float element count (actual)*/
138     int coordinatesSize = 0;
139     /*! \brief The kernelParams.atoms.coordinates float element count (reserved) */
140     int coordinatesSizeAlloc = 0;
141     /*! \brief The kernelParams.atoms.forces float element count (actual) */
142     int forcesSize = 0;
143     /*! \brief The kernelParams.atoms.forces float element count (reserved) */
144     int forcesSizeAlloc = 0;
145     /*! \brief The kernelParams.atoms.gridlineIndices int element count (actual) */
146     int gridlineIndicesSize = 0;
147     /*! \brief The kernelParams.atoms.gridlineIndices int element count (reserved) */
148     int gridlineIndicesSizeAlloc = 0;
149     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (actual) */
150     int splineDataSize = 0;
151     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (reserved) */
152     int splineDataSizeAlloc = 0;
153     /*! \brief The kernelParams.atoms.coefficients float element count (actual) */
154     int coefficientsSize[NUMFEPSTATES] = { 0, 0 };
155     /*! \brief The kernelParams.atoms.coefficients float element count (reserved) */
156     int coefficientsCapacity[NUMFEPSTATES] = { 0, 0 };
157     /*! \brief The kernelParams.grid.splineValuesArray float element count (actual) */
158     int splineValuesSize[NUMFEPSTATES] = { 0, 0 };
159     /*! \brief The kernelParams.grid.splineValuesArray float element count (reserved) */
160     int splineValuesCapacity[NUMFEPSTATES] = { 0, 0 };
161     /*! \brief The kernelParams.grid.realGrid float element count (actual) */
162     int realGridSize[NUMFEPSTATES] = { 0, 0 };
163     /*! \brief The kernelParams.grid.realGrid float element count (reserved) */
164     int realGridCapacity[NUMFEPSTATES] = { 0, 0 };
165     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (actual) */
166     int complexGridSize[NUMFEPSTATES] = { 0, 0 };
167     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (reserved) */
168     int complexGridCapacity[NUMFEPSTATES] = { 0, 0 };
169 };
170
171 #endif