2 * This file is part of the GROMACS molecular simulation package.
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.
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.
39 * CUDA kernels for GPU versions of copy_rvec_to_nbat_real and add_nbat_f_to_f.
41 * \author Alan Gray <alang@nvidia.com>
42 * \author Jon Vincent <jvincent@nvidia.com>
45 #include "gromacs/gpu_utils/vectype_ops.cuh"
46 #include "gromacs/nbnxm/nbnxm.h"
48 /*! \brief CUDA kernel for transforming position coordinates from rvec to nbnxm layout.
51 * - improve/simplify/document use of cxy_na and na_round
52 * - rename kernel so naming matches with the other NBNXM kernels;
53 * - enable separate compilation unit
55 * \param[in] numColumns extent of cell-level parallelism
56 * \param[out] gm_coordinatesNbnxm coordinates buffer in nbnxm layout
57 * \param[in] setFillerCoords tells whether to set the coordinates of the filler particles
58 * \param[in] gm_coordinatesRvec coordinates buffer in rvec format
59 * \param[in] gm_atomIndex atom index mapping
60 * \param[in] gm_numAtoms array of number of atoms
61 * \param[in] gm_cellIndex array of cell indices
62 * \param[in] cellOffset first cell
63 * \param[in] numAtomsPerCell number of atoms per cell
65 __global__ void nbnxn_gpu_x_to_nbat_x_kernel(int numColumns,
66 float* __restrict__ gm_coordinatesNbnxm,
68 const rvec* __restrict__ gm_coordinatesRvec,
69 const int* __restrict__ gm_atomIndex,
70 const int* __restrict__ gm_numAtoms,
71 const int* __restrict__ gm_cellIndex,
76 __global__ void nbnxn_gpu_x_to_nbat_x_kernel(int numColumns,
77 float* __restrict__ gm_coordinatesNbnxm,
79 const rvec* __restrict__ gm_coordinatesRvec,
80 const int* __restrict__ gm_atomIndex,
81 const int* __restrict__ gm_numAtoms,
82 const int* __restrict__ gm_cellIndex,
88 const float farAway = -1000000.0f;
90 /* map cell-level parallelism to y component of CUDA block index */
96 int na = gm_numAtoms[cxy];
97 int a0 = (cellOffset + gm_cellIndex[cxy]) * numAtomsPerCell;
101 // TODO: This can be done more efficiently
102 na_round = (gm_cellIndex[cxy + 1] - gm_cellIndex[cxy]) * numAtomsPerCell;
106 /* We fill only the real particle locations.
107 * We assume the filling entries at the end have been
108 * properly set before during pair-list generation.
113 /* map parallelism within a cell to x component of CUDA block index linearized
114 * with threads within a block */
116 i = blockIdx.x * blockDim.x + threadIdx.x;
118 j0 = a0 * STRIDE_XYZQ;
120 // destination address where x shoud be stored in nbnxm layout
121 float3* gm_coordinatesDest = (float3*)&gm_coordinatesNbnxm[j0 + 4 * i];
123 /* perform conversion of each element */
128 *gm_coordinatesDest = *((float3*)gm_coordinatesRvec[gm_atomIndex[a0 + i]]);
132 *gm_coordinatesDest = make_float3(farAway);
138 /*! \brief CUDA kernel to sum up the force components
140 * \tparam accumulateForce If the initial forces in \p gm_fTotal should be saved.
141 * \tparam addPmeForce Whether the PME force should be added to the total.
143 * \param[in] gm_forcesNbnxm Non-bonded forces in nbnxm format.
144 * \param[in] gm_forcesPme PME forces.
145 * \param[in,out] gm_forcesTotal Force buffer to be reduced into.
146 * \param[in] cell Cell index mapping.
147 * \param[in] atomStart Start atom index.
148 * \param[in] numAtoms Number of atoms.
150 template<bool accumulateForce, bool addPmeForce>
151 __global__ void nbnxn_gpu_add_nbat_f_to_f_kernel(const float3* __restrict__ gb_forcesNbnxm,
152 const float3* __restrict__ gm_forcesPme,
153 float3* gm_forcesTotal,
154 const int* __restrict__ gm_cell,
157 template<bool accumulateForce, bool addPmeForce>
158 __global__ void nbnxn_gpu_add_nbat_f_to_f_kernel(const float3* __restrict__ gb_forcesNbnxm,
159 const float3* __restrict__ gm_forcesPme,
160 float3* gm_forcesTotal,
161 const int* __restrict__ gm_cell,
166 /* map particle-level parallelism to 1D CUDA thread and block index */
167 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
169 /* perform addition for each particle*/
170 if (threadIndex < numAtoms)
173 int i = gm_cell[atomStart + threadIndex];
174 float3* gm_forcesDest = (float3*)&gm_forcesTotal[atomStart + threadIndex];
179 temp = *gm_forcesDest;
180 temp += gb_forcesNbnxm[i];
184 temp = gb_forcesNbnxm[i];
188 temp += gm_forcesPme[atomStart + threadIndex];
190 *gm_forcesDest = temp;