2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
36 * \brief Implements common util routines for different NBNXN GPU implementations
38 * \author Aleksei Iupinov <a.yupinov@gmail.com>
39 * \ingroup module_nbnxm
42 #ifndef GMX_NBNXM_GPU_COMMON_UTILS_H
43 #define GMX_NBNXM_GPU_COMMON_UTILS_H
47 #include "gromacs/listed_forces/gpubonded.h"
48 #include "gromacs/utility/exceptions.h"
49 #include "gromacs/utility/range.h"
50 #include "gromacs/nbnxm/nbnxm_gpu.h"
53 # include "cuda/nbnxm_cuda_types.h"
57 # include "opencl/nbnxm_ocl_types.h"
63 /*! \brief An early return condition for empty NB GPU workloads
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.
70 static inline bool canSkipNonbondedWork(const NbnxmGpu& nb, InteractionLocality iloc)
72 assert(nb.plist[iloc]);
73 return (iloc == InteractionLocality::NonLocal && nb.plist[iloc]->nsci == 0);
76 /*! \brief Convert atom locality to interaction locality.
78 * In the current implementation the this is straightforward conversion:
79 * local to local, non-local to non-local.
81 * \param[in] atomLocality Atom locality specifier
82 * \returns Interaction locality corresponding to the atom locality passed.
84 static inline InteractionLocality gpuAtomToInteractionLocality(const AtomLocality atomLocality)
87 /* determine interaction locality from atom locality */
88 if (atomLocality == AtomLocality::Local)
90 return InteractionLocality::Local;
92 else if (atomLocality == AtomLocality::NonLocal)
94 return InteractionLocality::NonLocal;
98 GMX_THROW(gmx::InconsistentInputError(
99 "Only Local and NonLocal atom locities can be converted to interaction locality."));
103 /*! \brief Returns true if there is GPU short-range work for the given interaction locality.
105 * Note that as, unlike nonbonded tasks, bonded tasks are not split into local/nonlocal,
106 * and therefore if there are GPU offloaded bonded interactions, this function will return
107 * true for all interaction localities.
109 * \param[inout] nb Pointer to the nonbonded GPU data structure
110 * \param[in] iLocality Interaction locality identifier
112 static inline bool haveGpuShortRangeWork(const NbnxmGpu& nb, const gmx::InteractionLocality iLocality)
114 return nb.haveWork[iLocality];
117 /*! \brief Calculate atom range and return start index and length.
119 * \param[in] atomData Atom descriptor data structure
120 * \param[in] atomLocality Atom locality specifier
121 * \returns Range of indexes for selected locality.
123 static inline gmx::Range<int> getGpuAtomRange(const NBAtomData* atomData, const AtomLocality atomLocality)
127 /* calculate the atom data index range based on locality */
128 if (atomLocality == AtomLocality::Local)
130 return gmx::Range<int>(0, atomData->numAtomsLocal);
132 else if (atomLocality == AtomLocality::NonLocal)
134 return gmx::Range<int>(atomData->numAtomsLocal, atomData->numAtoms);
138 GMX_THROW(gmx::InconsistentInputError(
139 "Only Local and NonLocal atom locities can be used to get atom ranges in NBNXM."));