Update instructions in containers.rst
[gromacs.git] / src / gromacs / nbnxm / cuda / nbnxm_buffer_ops_kernels.cuh
blob9a9ffc6c1ce8ca46bd4f726eb4fab006e1c2d2cd
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019,2020, 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  */
36 /*! \internal \file
37  *
38  * \brief
39  * CUDA kernels for GPU versions of copy_rvec_to_nbat_real and add_nbat_f_to_f.
40  *
41  *  \author Alan Gray <alang@nvidia.com>
42  *  \author Jon Vincent <jvincent@nvidia.com>
43  */
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.
49  *
50  * TODO:
51  *  - rename kernel so naming matches with the other NBNXM kernels;
52  *  - enable separate compilation unit
54  * \param[in]     numColumns          Extent of cell-level parallelism.
55  * \param[out]    gm_xq               Coordinates buffer in nbnxm layout.
56  * \tparam        setFillerCoords     Whether to set the coordinates of the filler particles.
57  * \param[in]     gm_x                Coordinates buffer.
58  * \param[in]     gm_atomIndex        Atom index mapping.
59  * \param[in]     gm_numAtoms         Array of number of atoms.
60  * \param[in]     gm_cellIndex        Array of cell indices.
61  * \param[in]     cellOffset          First cell.
62  * \param[in]     numAtomsPerCell     Number of atoms per cell.
63  */
64 template<bool setFillerCoords>
65 static __global__ void nbnxn_gpu_x_to_nbat_x_kernel(int numColumns,
66                                                     float4* __restrict__ gm_xq,
67                                                     const float3* __restrict__ gm_x,
68                                                     const int* __restrict__ gm_atomIndex,
69                                                     const int* __restrict__ gm_numAtoms,
70                                                     const int* __restrict__ gm_cellIndex,
71                                                     int cellOffset,
72                                                     int numAtomsPerCell)
76     const float farAway = -1000000.0f;
78     // Map cell-level parallelism to y component of CUDA block index.
79     int cxy = blockIdx.y;
81     if (cxy < numColumns)
82     {
84         const int numAtoms = gm_numAtoms[cxy];
85         const int offset   = (cellOffset + gm_cellIndex[cxy]) * numAtomsPerCell;
86         int       numAtomsRounded;
87         if (setFillerCoords)
88         {
89             // TODO: This can be done more efficiently
90             numAtomsRounded = (gm_cellIndex[cxy + 1] - gm_cellIndex[cxy]) * numAtomsPerCell;
91         }
92         else
93         {
94             // We fill only the real particle locations.
95             // We assume the filling entries at the end have been
96             // properly set before during pair-list generation.
97             numAtomsRounded = numAtoms;
98         }
100         const int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
102         // Destination address where x should be stored in nbnxm layout. We use this cast here to
103         // save only x, y and z components, not touching the w (q) component, which is pre-defined.
104         float3* gm_xqDest = (float3*)&gm_xq[threadIndex + offset];
106         // Perform layout conversion of each element.
107         if (threadIndex < numAtomsRounded)
108         {
109             if (threadIndex < numAtoms)
110             {
111                 *gm_xqDest = gm_x[gm_atomIndex[threadIndex + offset]];
112             }
113             else
114             {
115                 *gm_xqDest = make_float3(farAway);
116             }
117         }
118     }