Merge branch release-2019
[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, 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
49 #include <array>
50 #include <set>
51 #include <vector>
52
53 #if GMX_GPU == GMX_GPU_CUDA
54 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
55 #include "gromacs/gpu_utils/gpuregiontimer.cuh"
56 #elif GMX_GPU == GMX_GPU_OPENCL
57 #include "gromacs/gpu_utils/gpueventsynchronizer_ocl.h"
58 #include "gromacs/gpu_utils/gpuregiontimer_ocl.h"
59 #endif
60
61 #include "gromacs/timing/gpu_timing.h" // for gtPME_EVENT_COUNT
62
63 class GpuParallel3dFft;
64
65 /*! \internal \brief
66  * The main PME CUDA/OpenCL-specific host data structure, included in the PME GPU structure by the archSpecific pointer.
67  */
68 struct PmeGpuSpecific
69 {
70     /*! \brief The GPU stream where everything related to the PME happens. */
71     CommandStream pmeStream;
72
73     /*! \brief
74      * A handle to the GPU context.
75      * TODO: this is currently extracted from the implementation of pmeGpu->programHandle_,
76      * but should be a constructor parameter to PmeGpu, as well as PmeGpuProgram,
77      * managed by high-level code.
78      */
79     Context context;
80
81     /* Synchronization events */
82     /*! \brief Triggered after the grid has been copied to the host (after the spreading stage). */
83     GpuEventSynchronizer syncSpreadGridD2H;
84
85     /* Settings which are set at the start of the run */
86     /*! \brief A boolean which tells whether the complex and real grids for cu/clFFT are different or same. Currenty true. */
87     bool performOutOfPlaceFFT;
88     /*! \brief A boolean which tells if the GPU timing events are enabled.
89      *  False by default, can be enabled by setting the environment variable GMX_ENABLE_GPU_TIMING.
90      *  Note: will not be reliable when multiple GPU tasks are running concurrently on the same device context,
91      * as CUDA events on multiple streams are untrustworthy.
92      */
93     bool                                             useTiming;
94
95     //! Vector of FFT setups
96     std::vector<std::unique_ptr<GpuParallel3dFft > > fftSetup;
97
98     //! All the timers one might use
99     std::array<GpuRegionTimer, gtPME_EVENT_COUNT>    timingEvents;
100
101     //! Indices of timingEvents actually used
102     std::set<size_t>                                 activeTimers;
103
104     /* GPU arrays element counts (not the arrays sizes in bytes!).
105      * They might be larger than the actual meaningful data sizes.
106      * These are paired: the actual element count + the maximum element count that can fit in the current allocated memory.
107      * These integer pairs are mostly meaningful for the reallocateDeviceBuffer calls.
108      * As such, if DeviceBuffer is refactored into a class, they can be freely changed, too.
109      * The only exceptions are realGridSize and complexGridSize which are also used for grid clearing/copying.
110      * TODO: these should live in a clean buffered container type, and be refactored in the NB/cudautils as well.
111      */
112     /*! \brief The kernelParams.atoms.coordinates float element count (actual)*/
113     int coordinatesSize;
114     /*! \brief The kernelParams.atoms.coordinates float element count (reserved) */
115     int coordinatesSizeAlloc;
116     /*! \brief The kernelParams.atoms.forces float element count (actual) */
117     int forcesSize;
118     /*! \brief The kernelParams.atoms.forces float element count (reserved) */
119     int forcesSizeAlloc;
120     /*! \brief The kernelParams.atoms.gridlineIndices int element count (actual) */
121     int gridlineIndicesSize;
122     /*! \brief The kernelParams.atoms.gridlineIndices int element count (reserved) */
123     int gridlineIndicesSizeAlloc;
124     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (actual) */
125     int splineDataSize;
126     /*! \brief Both the kernelParams.atoms.theta and kernelParams.atoms.dtheta float element count (reserved) */
127     int splineDataSizeAlloc;
128     /*! \brief The kernelParams.atoms.coefficients float element count (actual) */
129     int coefficientsSize;
130     /*! \brief The kernelParams.atoms.coefficients float element count (reserved) */
131     int coefficientsSizeAlloc;
132     /*! \brief The kernelParams.grid.splineValuesArray float element count (actual) */
133     int splineValuesSize;
134     /*! \brief The kernelParams.grid.splineValuesArray float element count (reserved) */
135     int splineValuesSizeAlloc;
136     /*! \brief The kernelParams.grid.realGrid float element count (actual) */
137     int realGridSize;
138     /*! \brief The kernelParams.grid.realGrid float element count (reserved) */
139     int realGridSizeAlloc;
140     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (actual) */
141     int complexGridSize;
142     /*! \brief The kernelParams.grid.fourierGrid float (not float2!) element count (reserved) */
143     int complexGridSizeAlloc;
144 };
145
146 #endif