2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
38 * SYCL implementation of coordinate layout conversion kernel and its wrapper
40 * \ingroup module_nbnxm
44 #include "gromacs/gpu_utils/devicebuffer.h"
45 #include "gromacs/gpu_utils/gmxsycl.h"
46 #include "gromacs/nbnxm/grid.h"
47 #include "gromacs/nbnxm/nbnxm_gpu_buffer_ops_internal.h"
49 #include "nbnxm_sycl_types.h"
51 using cl::sycl::access::mode;
52 using cl::sycl::access::target;
57 /*! \brief SYCL kernel for transforming position coordinates from rvec to nbnxm layout.
59 * \param cgh SYCL's command group handler.
60 * \param[out] a_xq Coordinates buffer in nbnxm layout.
61 * \param[in] a_x Coordinates buffer.
62 * \param[in] a_atomIndex Atom index mapping.
63 * \param[in] a_numAtoms Array of number of atoms.
64 * \param[in] a_cellIndex Array of cell indices.
65 * \param[in] cellOffset First cell.
66 * \param[in] numAtomsPerCell Number of atoms per cell.
67 * \param[in] columnsOffset Index if the first column in the cell.
69 static auto nbnxmKernelTransformXToXq(cl::sycl::handler& cgh,
70 DeviceAccessor<Float4, mode::read_write> a_xq,
71 DeviceAccessor<Float3, mode::read> a_x,
72 DeviceAccessor<int, mode::read> a_atomIndex,
73 DeviceAccessor<int, mode::read> a_numAtoms,
74 DeviceAccessor<int, mode::read> a_cellIndex,
81 cgh.require(a_atomIndex);
82 cgh.require(a_numAtoms);
83 cgh.require(a_cellIndex);
85 return [=](cl::sycl::id<2> itemIdx) {
86 // Map cell-level parallelism to y component of block index.
87 const int cxy = itemIdx.get(1) + columnsOffset;
89 const int numAtoms = a_numAtoms[cxy];
90 const int offset = (cellOffset + a_cellIndex[cxy]) * numAtomsPerCell;
92 const int threadIndex = itemIdx.get(0);
94 // Perform layout conversion of each element.
95 if (threadIndex < numAtoms)
97 const float q = a_xq[threadIndex + offset][3];
98 const Float3 xNew = a_x[a_atomIndex[threadIndex + offset]];
99 a_xq[threadIndex + offset] = Float4(xNew[0], xNew[1], xNew[2], q);
104 // SYCL 1.2.1 requires providing a unique type for a kernel. Should not be needed for SYCL2020.
105 class NbnxmKernelTransformXToXqName;
107 void launchNbnxmKernelTransformXToXq(const Grid& grid,
109 DeviceBuffer<Float3> d_x,
110 const DeviceStream& deviceStream,
111 unsigned int numColumnsMax,
114 const unsigned int numColumns = grid.numColumns();
115 const unsigned int numAtomsMax = grid.numCellsColumnMax() * grid.numAtomsPerCell();
116 GMX_ASSERT(numColumns <= numColumnsMax, "Grid has more columns than allowed");
118 const cl::sycl::range<2> globalSize{ numAtomsMax, numColumns };
119 cl::sycl::queue q = deviceStream.stream();
121 q.submit([&](cl::sycl::handler& cgh) {
122 auto kernel = nbnxmKernelTransformXToXq(cgh,
129 grid.numAtomsPerCell(),
130 numColumnsMax * gridId);
131 cgh.parallel_for<NbnxmKernelTransformXToXqName>(globalSize, kernel);