Unify gpu_copy_xq_to_gpu(...) function
[alexxy/gromacs.git] / src / gromacs / nbnxm / gpu_common_utils.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2017,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 Implements common util routines for different NBNXN GPU implementations
37  *
38  * \author Aleksei Iupinov <a.yupinov@gmail.com>
39  * \ingroup module_nbnxm
40  */
41
42 #ifndef GMX_NBNXM_GPU_COMMON_UTILS_H
43 #define GMX_NBNXM_GPU_COMMON_UTILS_H
44
45 #include "config.h"
46
47 #include "gromacs/listed_forces/gpubonded.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/range.h"
50 #include "gromacs/nbnxm/nbnxm_gpu.h"
51
52 #if GMX_GPU_CUDA
53 #    include "cuda/nbnxm_cuda_types.h"
54 #endif
55
56 #if GMX_GPU_OPENCL
57 #    include "opencl/nbnxm_ocl_types.h"
58 #endif
59
60 namespace Nbnxm
61 {
62
63 /*! \brief An early return condition for empty NB GPU workloads
64  *
65  * This is currently used for non-local kernels/transfers only.
66  * Skipping the local kernel is more complicated, since the
67  * local part of the force array also depends on the non-local kernel.
68  * The skip of the local kernel is taken care of separately.
69  */
70 static inline bool canSkipNonbondedWork(const NbnxmGpu& nb, InteractionLocality iloc)
71 {
72     assert(nb.plist[iloc]);
73     return (iloc == InteractionLocality::NonLocal && nb.plist[iloc]->nsci == 0);
74 }
75
76 /*! \brief Check that atom locality values are valid for the GPU module.
77  *
78  *  In the GPU module atom locality "all" is not supported, the local and
79  *  non-local ranges are treated separately.
80  *
81  *  \param[in] atomLocality atom locality specifier
82  */
83 static inline void validateGpuAtomLocality(const AtomLocality atomLocality)
84 {
85     std::string str = gmx::formatString(
86             "Invalid atom locality passed (%d); valid here is only "
87             "local (%d) or nonlocal (%d)",
88             static_cast<int>(atomLocality),
89             static_cast<int>(AtomLocality::Local),
90             static_cast<int>(AtomLocality::NonLocal));
91
92     GMX_ASSERT(atomLocality == AtomLocality::Local || atomLocality == AtomLocality::NonLocal, str.c_str());
93 }
94
95 /*! \brief Convert atom locality to interaction locality.
96  *
97  *  In the current implementation the this is straightforward conversion:
98  *  local to local, non-local to non-local.
99  *
100  *  \param[in] atomLocality Atom locality specifier
101  *  \returns                Interaction locality corresponding to the atom locality passed.
102  */
103 static inline InteractionLocality gpuAtomToInteractionLocality(const AtomLocality atomLocality)
104 {
105     validateGpuAtomLocality(atomLocality);
106
107     /* determine interaction locality from atom locality */
108     if (atomLocality == AtomLocality::Local)
109     {
110         return InteractionLocality::Local;
111     }
112     else if (atomLocality == AtomLocality::NonLocal)
113     {
114         return InteractionLocality::NonLocal;
115     }
116     else
117     {
118         gmx_incons("Wrong locality");
119     }
120 }
121
122 /*! \brief Returns true if there is GPU short-range work for the given interaction locality.
123  *
124  * Note that as, unlike nonbonded tasks, bonded tasks are not split into local/nonlocal,
125  * and therefore if there are GPU offloaded bonded interactions, this function will return
126  * true for all interaction localities.
127  *
128  * \param[inout]  nb        Pointer to the nonbonded GPU data structure
129  * \param[in]     iLocality Interaction locality identifier
130  */
131 static inline bool haveGpuShortRangeWork(const NbnxmGpu& nb, const gmx::InteractionLocality iLocality)
132 {
133     return nb.haveWork[iLocality];
134 }
135
136 /*! \brief Calculate atom range and return start index and length.
137  *
138  * \param[in] atomData Atom descriptor data structure
139  * \param[in] atomLocality Atom locality specifier
140  * \returns Range of indexes for selected locality.
141  */
142 static inline gmx::Range<int> getGpuAtomRange(const NBAtomData* atomData, const AtomLocality atomLocality)
143 {
144     assert(atomData);
145     validateGpuAtomLocality(atomLocality);
146
147     /* calculate the atom data index range based on locality */
148     if (atomLocality == AtomLocality::Local)
149     {
150         return gmx::Range<int>(0, atomData->numAtomsLocal);
151     }
152     else
153     {
154         return gmx::Range<int>(atomData->numAtomsLocal, atomData->numAtoms);
155     }
156 }
157
158 } // namespace Nbnxm
159
160 #endif