1a130a6fb5ed3dc3b13d6021e7f328246e598bc1
[alexxy/gromacs.git] / src / gromacs / nbnxm / sycl / nbnxm_sycl.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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  */
35
36 /*! \internal \file
37  *  \brief
38  *  Data management and kernel launch functions for nbnxm sycl.
39  *
40  *  \ingroup module_nbnxm
41  */
42 #include "gmxpre.h"
43
44 #include "gromacs/nbnxm/gpu_common.h"
45 #include "gromacs/utility/exceptions.h"
46
47 #include "nbnxm_sycl_kernel.h"
48 #include "nbnxm_sycl_kernel_pruneonly.h"
49 #include "nbnxm_sycl_types.h"
50
51 namespace Nbnxm
52 {
53
54 /*! \brief
55  * Launch asynchronously the download of nonbonded forces from the GPU
56  * (and energies/shift forces if required).
57  */
58 void gpu_launch_cpyback(NbnxmGpu*                nb,
59                         struct nbnxn_atomdata_t* nbatom,
60                         const gmx::StepWorkload& stepWork,
61                         const AtomLocality       atomLocality)
62 {
63     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
64
65     const InteractionLocality iloc = gpuAtomToInteractionLocality(atomLocality);
66     GMX_ASSERT(iloc == InteractionLocality::Local
67                        || (iloc == InteractionLocality::NonLocal && nb->bNonLocalStreamDoneMarked == false),
68                "Non-local stream is indicating that the copy back event is enqueued at the "
69                "beginning of the copy back function.");
70
71     const DeviceStream& deviceStream = *nb->deviceStreams[iloc];
72     NBAtomData*         adat         = nb->atdat;
73
74     /* don't launch non-local copy-back if there was no non-local work to do */
75     if ((iloc == InteractionLocality::NonLocal) && !haveGpuShortRangeWork(*nb, iloc))
76     {
77         nb->bNonLocalStreamDoneMarked = false;
78         return;
79     }
80
81     /* local/nonlocal offset and length used for xq and f */
82     auto atomsRange = getGpuAtomRange(adat, atomLocality);
83
84     // With DD the local D2H transfer can only start after the non-local kernel has finished.
85     if (iloc == InteractionLocality::Local && nb->bNonLocalStreamDoneMarked)
86     {
87         nb->nonlocal_done.waitForEvent();
88         nb->bNonLocalStreamDoneMarked = false;
89     }
90
91     /* DtoH f
92      * Skip if buffer ops / reduction is offloaded to the GPU.
93      */
94     if (!stepWork.useGpuFBufferOps)
95     {
96         GMX_ASSERT(adat->f.elementSize() == sizeof(Float3),
97                    "The size of the force buffer element should be equal to the size of float3.");
98         copyFromDeviceBuffer(reinterpret_cast<Float3*>(nbatom->out[0].f.data()) + atomsRange.begin(),
99                              &adat->f,
100                              atomsRange.begin(),
101                              atomsRange.size(),
102                              deviceStream,
103                              GpuApiCallBehavior::Async,
104                              nullptr);
105     }
106
107     /* After the non-local D2H is launched the nonlocal_done event can be
108        recorded which signals that the local D2H can proceed. This event is not
109        placed after the non-local kernel because we want the non-local data
110        back first. */
111     if (iloc == InteractionLocality::NonLocal)
112     {
113         nb->nonlocal_done.markEvent(deviceStream);
114         nb->bNonLocalStreamDoneMarked = true;
115     }
116
117     /* only transfer energies in the local stream */
118     if (iloc == InteractionLocality::Local)
119     {
120         /* DtoH fshift when virial is needed */
121         if (stepWork.computeVirial)
122         {
123             GMX_ASSERT(sizeof(*nb->nbst.fShift) == adat->fShift.elementSize(),
124                        "Sizes of host- and device-side shift vector elements should be the same.");
125             copyFromDeviceBuffer(
126                     nb->nbst.fShift, &adat->fShift, 0, SHIFTS, deviceStream, GpuApiCallBehavior::Async, nullptr);
127         }
128
129         /* DtoH energies */
130         if (stepWork.computeEnergy)
131         {
132             GMX_ASSERT(sizeof(*nb->nbst.eLJ) == sizeof(float),
133                        "Sizes of host- and device-side LJ energy terms should be the same.");
134             copyFromDeviceBuffer(
135                     nb->nbst.eLJ, &adat->eLJ, 0, 1, deviceStream, GpuApiCallBehavior::Async, nullptr);
136             GMX_ASSERT(sizeof(*nb->nbst.eElec) == sizeof(float),
137                        "Sizes of host- and device-side electrostatic energy terms should be the "
138                        "same.");
139             copyFromDeviceBuffer(
140                     nb->nbst.eElec, &adat->eElec, 0, 1, deviceStream, GpuApiCallBehavior::Async, nullptr);
141         }
142     }
143 }
144
145 void gpu_launch_kernel_pruneonly(NbnxmGpu* nb, const InteractionLocality iloc, const int numParts)
146 {
147     gpu_plist* plist = nb->plist[iloc];
148
149     if (plist->haveFreshList)
150     {
151         GMX_ASSERT(numParts == 1, "With first pruning we expect 1 part");
152
153         /* Set rollingPruningNumParts to signal that it is not set */
154         plist->rollingPruningNumParts = 0;
155         plist->rollingPruningPart     = 0;
156     }
157     else
158     {
159         if (plist->rollingPruningNumParts == 0)
160         {
161             plist->rollingPruningNumParts = numParts;
162         }
163         else
164         {
165             GMX_ASSERT(numParts == plist->rollingPruningNumParts,
166                        "It is not allowed to change numParts in between list generation steps");
167         }
168     }
169
170     /* Use a local variable for part and update in plist, so we can return here
171      * without duplicating the part increment code.
172      */
173     const int part = plist->rollingPruningPart;
174
175     plist->rollingPruningPart++;
176     if (plist->rollingPruningPart >= plist->rollingPruningNumParts)
177     {
178         plist->rollingPruningPart = 0;
179     }
180
181     /* Compute the number of list entries to prune in this pass */
182     const int numSciInPart = (plist->nsci - part) / numParts;
183
184     /* Don't launch the kernel if there is no work to do */
185     if (numSciInPart <= 0)
186     {
187         plist->haveFreshList = false;
188         return;
189     }
190
191     launchNbnxmKernelPruneOnly(nb, iloc, numParts, part, numSciInPart);
192
193     if (plist->haveFreshList)
194     {
195         plist->haveFreshList = false;
196         nb->didPrune[iloc]   = true; // Mark that pruning has been done
197     }
198     else
199     {
200         nb->didRollingPrune[iloc] = true; // Mark that rolling pruning has been done
201     }
202 }
203
204 void gpu_launch_kernel(NbnxmGpu* nb, const gmx::StepWorkload& stepWork, const Nbnxm::InteractionLocality iloc)
205 {
206     const NBParamGpu* nbp   = nb->nbparam;
207     gpu_plist*        plist = nb->plist[iloc];
208
209     if (canSkipNonbondedWork(*nb, iloc))
210     {
211         plist->haveFreshList = false;
212         return;
213     }
214
215     if (nbp->useDynamicPruning && plist->haveFreshList)
216     {
217         // Prunes for rlistOuter and rlistInner, sets plist->haveFreshList=false
218         gpu_launch_kernel_pruneonly(nb, iloc, 1);
219     }
220
221     if (plist->nsci == 0)
222     {
223         /* Don't launch an empty local kernel */
224         return;
225     }
226
227     launchNbnxmKernel(nb, stepWork, iloc);
228 }
229
230 } // namespace Nbnxm