f1b1a6db8178a2a20aeb71ae4b251d815edf08a8
[alexxy/gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_cuda_types.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2012, The GROMACS development team.
6  * Copyright (c) 2013-2019,2020,2021, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37
38 /*! \internal \file
39  *  \brief
40  *  Data types used internally in the nbnxn_cuda module.
41  *
42  *  \author Szilárd Páll <pall.szilard@gmail.com>
43  *  \ingroup module_nbnxm
44  */
45
46 #ifndef NBNXM_CUDA_TYPES_H
47 #define NBNXM_CUDA_TYPES_H
48
49 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/gpu_utils/devicebuffer.h"
52 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
53 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
54 #include "gromacs/gpu_utils/gputraits.cuh"
55 #include "gromacs/mdtypes/interaction_const.h"
56 #include "gromacs/nbnxm/gpu_types_common.h"
57 #include "gromacs/nbnxm/nbnxm.h"
58 #include "gromacs/nbnxm/pairlist.h"
59 #include "gromacs/timing/gpu_timing.h"
60 #include "gromacs/utility/enumerationhelpers.h"
61
62 /* TODO: consider moving this to kernel_utils */
63 /* Convenience defines */
64 /*! \brief cluster size = number of atoms per cluster. */
65 static constexpr int c_clSize = c_nbnxnGpuClusterSize;
66
67 /* All structs prefixed with "cu_" hold data used in GPU calculations and
68  * are passed to the kernels, except cu_timers_t. */
69 /*! \cond */
70 typedef struct cu_atomdata cu_atomdata_t;
71 /*! \endcond */
72
73
74 /** \internal
75  * \brief Staging area for temporary data downloaded from the GPU.
76  *
77  *  The energies/shift forces get downloaded here first, before getting added
78  *  to the CPU-side aggregate values.
79  */
80 struct nb_staging_t
81 {
82     //! LJ energy
83     float* e_lj = nullptr;
84     //! electrostatic energy
85     float* e_el = nullptr;
86     //! shift forces
87     Float3* fshift = nullptr;
88 };
89
90 /** \internal
91  * \brief Nonbonded atom data - both inputs and outputs.
92  */
93 struct cu_atomdata
94 {
95     //! number of atoms
96     int natoms;
97     //! number of local atoms
98     int natoms_local;
99     //! allocation size for the atom data (xq, f)
100     int nalloc;
101
102     //! atom coordinates + charges, size natoms
103     DeviceBuffer<Float4> xq;
104     //! force output array, size natoms
105     DeviceBuffer<Float3> f;
106
107     //! LJ energy output, size 1
108     DeviceBuffer<float> e_lj;
109     //! Electrostatics energy input, size 1
110     DeviceBuffer<float> e_el;
111
112     //! shift forces
113     DeviceBuffer<Float3> fshift;
114
115     //! number of atom types
116     int ntypes;
117     //! atom type indices, size natoms
118     DeviceBuffer<int> atom_types;
119     //! sqrt(c6),sqrt(c12) size natoms
120     DeviceBuffer<Float2> lj_comb;
121
122     //! shifts
123     DeviceBuffer<Float3> shift_vec;
124     //! true if the shift vector has been uploaded
125     bool bShiftVecUploaded;
126 };
127
128 /** \internal
129  * \brief Typedef of actual timer type.
130  */
131 typedef struct Nbnxm::gpu_timers_t cu_timers_t;
132
133 /*! \internal
134  * \brief Main data structure for CUDA nonbonded force calculations.
135  */
136 struct NbnxmGpu
137 {
138     /*! \brief GPU device context.
139      *
140      * \todo Make it constant reference, once NbnxmGpu is a proper class.
141      */
142     const DeviceContext* deviceContext_;
143     /*! \brief true if doing both local/non-local NB work on GPU */
144     bool bUseTwoStreams = false;
145     //! true indicates that the nonlocal_done event was marked
146     bool bNonLocalStreamDoneMarked = false;
147
148     /*! \brief atom data */
149     cu_atomdata_t* atdat = nullptr;
150     /*! \brief array of atom indices */
151     int* atomIndices = nullptr;
152     /*! \brief size of atom indices */
153     int atomIndicesSize = 0;
154     /*! \brief size of atom indices allocated in device buffer */
155     int atomIndicesSize_alloc = 0;
156     /*! \brief x buf ops num of atoms */
157     int* cxy_na = nullptr;
158     /*! \brief number of elements in cxy_na */
159     int ncxy_na = 0;
160     /*! \brief number of elements allocated allocated in device buffer */
161     int ncxy_na_alloc = 0;
162     /*! \brief x buf ops cell index mapping */
163     int* cxy_ind = nullptr;
164     /*! \brief number of elements in cxy_ind */
165     int ncxy_ind = 0;
166     /*! \brief number of elements allocated allocated in device buffer */
167     int ncxy_ind_alloc = 0;
168     /*! \brief parameters required for the non-bonded calc. */
169     NBParamGpu* nbparam = nullptr;
170     /*! \brief pair-list data structures (local and non-local) */
171     gmx::EnumerationArray<Nbnxm::InteractionLocality, Nbnxm::gpu_plist*> plist = { { nullptr } };
172     /*! \brief staging area where fshift/energies get downloaded */
173     nb_staging_t nbst;
174     /*! \brief local and non-local GPU streams */
175     gmx::EnumerationArray<Nbnxm::InteractionLocality, const DeviceStream*> deviceStreams;
176
177     /*! \brief Event triggered when the non-local non-bonded
178      * kernel is done (and the local transfer can proceed) */
179     GpuEventSynchronizer nonlocal_done;
180     /*! \brief Event triggered when the tasks issued in the local
181      * stream that need to precede the non-local force or buffer
182      * operation calculations are done (e.g. f buffer 0-ing, local
183      * x/q H2D, buffer op initialization in local stream that is
184      * required also by nonlocal stream ) */
185     GpuEventSynchronizer misc_ops_and_local_H2D_done;
186
187     /*! \brief True if there is work for the current domain in the
188      * respective locality.
189      *
190      * This includes local/nonlocal GPU work, either bonded or
191      * nonbonded, scheduled to be executed in the current
192      * domain. As long as bonded work is not split up into
193      * local/nonlocal, if there is bonded GPU work, both flags
194      * will be true. */
195     gmx::EnumerationArray<Nbnxm::InteractionLocality, bool> haveWork = { { false } };
196
197     /* NOTE: With current CUDA versions (<=5.0) timing doesn't work with multiple
198      * concurrent streams, so we won't time if both l/nl work is done on GPUs.
199      * Timer init/uninit is still done even with timing off so only the condition
200      * setting bDoTime needs to be change if this CUDA "feature" gets fixed. */
201     /*! \brief True if event-based timing is enabled. */
202     bool bDoTime = false;
203     /*! \brief CUDA event-based timers. */
204     cu_timers_t* timers = nullptr;
205     /*! \brief Timing data. TODO: deprecate this and query timers for accumulated data instead */
206     gmx_wallclock_gpu_nbnxn_t* timings = nullptr;
207 };
208
209 #endif /* NBNXN_CUDA_TYPES_H */