2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
38 * Data management and kernel launch functions for nbnxm sycl.
40 * \ingroup module_nbnxm
44 #include "gromacs/nbnxm/gpu_common.h"
45 #include "gromacs/utility/exceptions.h"
47 #include "nbnxm_sycl_kernel.h"
48 #include "nbnxm_sycl_kernel_pruneonly.h"
49 #include "nbnxm_sycl_types.h"
55 void nbnxnInsertNonlocalGpuDependency(NbnxmGpu* nb, const InteractionLocality interactionLocality)
57 const DeviceStream& deviceStream = *nb->deviceStreams[interactionLocality];
58 if (nb->bUseTwoStreams)
60 if (interactionLocality == InteractionLocality::Local)
62 nb->misc_ops_and_local_H2D_done.markEvent(deviceStream);
66 nb->misc_ops_and_local_H2D_done.enqueueWaitEvent(deviceStream);
73 * Launch asynchronously the download of nonbonded forces from the GPU
74 * (and energies/shift forces if required).
76 void gpu_launch_cpyback(NbnxmGpu* nb,
77 struct nbnxn_atomdata_t* nbatom,
78 const gmx::StepWorkload& stepWork,
79 const AtomLocality atomLocality)
81 GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
83 const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
84 GMX_ASSERT(iloc == InteractionLocality::Local
85 || (iloc == InteractionLocality::NonLocal && nb->bNonLocalStreamDoneMarked == false),
86 "Non-local stream is indicating that the copy back event is enqueued at the "
87 "beginning of the copy back function.");
89 const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
90 sycl_atomdata_t* adat = nb->atdat;
92 /* don't launch non-local copy-back if there was no non-local work to do */
93 if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
95 nb->bNonLocalStreamDoneMarked = false;
99 int adatBegin, adatLen;
100 getGpuAtomRange(adat, atomLocality, &adatBegin, &adatLen);
102 // With DD the local D2H transfer can only start after the non-local kernel has finished.
103 if (iloc == InteractionLocality::Local && nb->bNonLocalStreamDoneMarked)
105 nb->nonlocal_done.waitForEvent();
106 nb->bNonLocalStreamDoneMarked = false;
110 * Skip if buffer ops / reduction is offloaded to the GPU.
112 if (!stepWork.useGpuFBufferOps)
114 GMX_ASSERT(adat->f.elementSize() == sizeof(Float3),
115 "The size of the force buffer element should be equal to the size of float3.");
116 copyFromDeviceBuffer(reinterpret_cast<Float3*>(nbatom->out[0].f.data()) + adatBegin,
121 GpuApiCallBehavior::Async,
125 /* After the non-local D2H is launched the nonlocal_done event can be
126 recorded which signals that the local D2H can proceed. This event is not
127 placed after the non-local kernel because we want the non-local data
129 if (iloc == InteractionLocality::NonLocal)
131 nb->nonlocal_done.markEvent(deviceStream);
132 nb->bNonLocalStreamDoneMarked = true;
135 /* only transfer energies in the local stream */
136 if (iloc == InteractionLocality::Local)
138 /* DtoH fshift when virial is needed */
139 if (stepWork.computeVirial)
141 GMX_ASSERT(sizeof(*nb->nbst.fshift) == adat->fShift.elementSize(),
142 "Sizes of host- and device-side shift vector elements should be the same.");
143 copyFromDeviceBuffer(
144 nb->nbst.fshift, &adat->fShift, 0, SHIFTS, deviceStream, GpuApiCallBehavior::Async, nullptr);
148 if (stepWork.computeEnergy)
150 GMX_ASSERT(sizeof(*nb->nbst.e_lj) == sizeof(float),
151 "Sizes of host- and device-side LJ energy terms should be the same.");
152 copyFromDeviceBuffer(
153 nb->nbst.e_lj, &adat->eLJ, 0, 1, deviceStream, GpuApiCallBehavior::Async, nullptr);
154 GMX_ASSERT(sizeof(*nb->nbst.e_el) == sizeof(float),
155 "Sizes of host- and device-side electrostatic energy terms should be the "
157 copyFromDeviceBuffer(
158 nb->nbst.e_el, &adat->eElec, 0, 1, deviceStream, GpuApiCallBehavior::Async, nullptr);
163 /*! \brief Launch asynchronously the xq buffer host to device copy. */
164 void gpu_copy_xq_to_gpu(NbnxmGpu* nb, const nbnxn_atomdata_t* nbatom, const AtomLocality atomLocality)
166 GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
167 validateGpuAtomLocality(atomLocality);
169 const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
171 sycl_atomdata_t* adat = nb->atdat;
172 gpu_plist* plist = nb->plist[iloc];
173 const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
175 /* Don't launch the non-local H2D copy if there is no dependent
176 work to do: neither non-local nor other (e.g. bonded) work
177 to do that has as input the nbnxn coordinates.
178 Doing the same for the local kernel is more complicated, since the
179 local part of the force array also depends on the non-local kernel.
180 So to avoid complicating the code and to reduce the risk of bugs,
181 we always call the local local x+q copy (and the rest of the local
182 work in nbnxn_gpu_launch_kernel().
184 if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
186 plist->haveFreshList = false;
188 // The event is marked for Local interactions unconditionally,
189 // so it has to be released here because of the early return
190 // for NonLocal interactions.
191 nb->misc_ops_and_local_H2D_done.reset();
196 int adatBegin, adatLen;
197 getGpuAtomRange(adat, atomLocality, &adatBegin, &adatLen);
200 GMX_ASSERT(adat->xq.elementSize() == sizeof(Float4),
201 "The size of the xyzq buffer element should be equal to the size of float4.");
202 copyToDeviceBuffer(&adat->xq,
203 reinterpret_cast<const Float4*>(nbatom->x().data()) + adatBegin,
207 GpuApiCallBehavior::Async,
210 /* No need to enforce stream synchronization with events like we do in CUDA/OpenCL.
211 * Runtime should do the scheduling correctly based on data dependencies.
212 * But for consistency's sake, we do it anyway. */
213 /* When we get here all misc operations issued in the local stream as well as
214 * the local xq H2D are done, so we record that in the local stream and wait for it in the
215 * nonlocal one. This wait needs to precede any PP tasks, bonded or nonbonded, that may
216 * compute on interactions between local and nonlocal atoms. */
217 nbnxnInsertNonlocalGpuDependency(nb, iloc);
220 void gpu_launch_kernel_pruneonly(NbnxmGpu* nb, const InteractionLocality iloc, const int numParts)
222 gpu_plist* plist = nb->plist[iloc];
224 if (plist->haveFreshList)
226 GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
228 /* Set rollingPruningNumParts to signal that it is not set */
229 plist->rollingPruningNumParts = 0;
230 plist->rollingPruningPart = 0;
234 if (plist->rollingPruningNumParts == 0)
236 plist->rollingPruningNumParts = numParts;
240 GMX_ASSERT(numParts == plist->rollingPruningNumParts,
241 "It is not allowed to change numParts in between list generation steps");
245 /* Use a local variable for part and update in plist, so we can return here
246 * without duplicating the part increment code.
248 const int part = plist->rollingPruningPart;
250 plist->rollingPruningPart++;
251 if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
253 plist->rollingPruningPart = 0;
256 /* Compute the number of list entries to prune in this pass */
257 const int numSciInPart = (plist->nsci - part) / numParts;
259 /* Don't launch the kernel if there is no work to do */
260 if (numSciInPart <= 0)
262 plist->haveFreshList = false;
266 launchNbnxmKernelPruneOnly(nb, iloc, numParts, part, numSciInPart);
268 if (plist->haveFreshList)
270 plist->haveFreshList = false;
271 nb->didPrune[iloc] = true; // Mark that pruning has been done
275 nb->didRollingPrune[iloc] = true; // Mark that rolling pruning has been done
279 void gpu_launch_kernel(NbnxmGpu* nb, const gmx::StepWorkload& stepWork, const Nbnxm::InteractionLocality iloc)
281 const NBParamGpu* nbp = nb->nbparam;
282 gpu_plist* plist = nb->plist[iloc];
284 if (canSkipNonbondedWork(*nb, iloc))
286 plist->haveFreshList = false;
290 if (nbp->useDynamicPruning && plist->haveFreshList)
292 // Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
293 gpu_launch_kernel_pruneonly(nb, iloc, 1);
296 if (plist->nsci == 0)
298 /* Don't launch an empty local kernel */
302 launchNbnxmKernel(nb, stepWork, iloc);