Avoid allocating SYCL buffer on each call to PME solve
[alexxy/gromacs.git] / src / gromacs / ewald / pme_solve_sycl.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 Implements PME GPU Fourier grid solving in SYCL.
38  *
39  *  \author Mark Abraham <mark.j.abraham@gmail.com>
40  */
41
42 #include "gmxpre.h"
43
44 #include "pme_solve_sycl.h"
45
46 #include <cassert>
47
48 #include "gromacs/gpu_utils/gmxsycl.h"
49 #include "gromacs/gpu_utils/sycl_kernel_utils.h"
50 #include "gromacs/math/units.h"
51
52 #include "pme_gpu_constants.h"
53
54 using cl::sycl::access::mode;
55
56 /*! \brief
57  * PME complex grid solver kernel function.
58  *
59  * \tparam     gridOrdering             Specifies the dimension ordering of the complex grid.
60  * \tparam     computeEnergyAndVirial   Tells if the reciprocal energy and virial should be
61  *                                        computed.
62  * \tparam     subGroupSize             Describes the width of a SYCL subgroup
63  */
64 template<GridOrdering gridOrdering, bool computeEnergyAndVirial, int subGroupSize>
65 auto makeSolveKernel(cl::sycl::handler&                cgh,
66                      DeviceAccessor<float, mode::read> a_splineModuli,
67                      SolveKernelParams                 solveKernelParams,
68                      OptionalAccessor<float, mode::read_write, computeEnergyAndVirial> a_virialAndEnergy,
69                      DeviceAccessor<float, mode::read_write> a_fourierGrid)
70 {
71     a_splineModuli.bind(cgh);
72     if constexpr (computeEnergyAndVirial)
73     {
74         a_virialAndEnergy.bind(cgh);
75     }
76     a_fourierGrid.bind(cgh);
77
78     /* Reduce 7 outputs per warp in the shared memory */
79     const int stride =
80             8; // this is c_virialAndEnergyCount==7 rounded up to power of 2 for convenience, hence the assert
81     static_assert(c_virialAndEnergyCount == 7);
82     const int reductionBufferSize = c_solveMaxWarpsPerBlock * stride;
83     cl::sycl::accessor<float, 1, mode::read_write, cl::sycl::target::local> sm_virialAndEnergy(
84             cl::sycl::range<1>(reductionBufferSize), cgh);
85
86     /* Each thread works on one cell of the Fourier space complex 3D grid (gm_grid).
87      * Each block handles up to c_solveMaxWarpsPerBlock * subGroupSize cells -
88      * depending on the grid contiguous dimension size,
89      * that can range from a part of a single gridline to several complete gridlines.
90      */
91     return [=](cl::sycl::nd_item<3> itemIdx) [[intel::reqd_sub_group_size(subGroupSize)]]
92     {
93         /* This kernel supports 2 different grid dimension orderings: YZX and XYZ */
94         int majorDim, middleDim, minorDim;
95         switch (gridOrdering)
96         {
97             case GridOrdering::YZX:
98                 majorDim  = YY;
99                 middleDim = ZZ;
100                 minorDim  = XX;
101                 break;
102
103             case GridOrdering::XYZ:
104                 majorDim  = XX;
105                 middleDim = YY;
106                 minorDim  = ZZ;
107                 break;
108
109             default: assert(false);
110         }
111
112         /* Global memory pointers */
113         const float* __restrict__ gm_splineValueMajor =
114                 a_splineModuli.get_pointer() + solveKernelParams.splineValuesOffset[majorDim];
115         const float* __restrict__ gm_splineValueMiddle =
116                 a_splineModuli.get_pointer() + solveKernelParams.splineValuesOffset[middleDim];
117         const float* __restrict__ gm_splineValueMinor =
118                 a_splineModuli.get_pointer() + solveKernelParams.splineValuesOffset[minorDim];
119         // The Fourier grid is allocated as float values, even though
120         // it logically contains complex values. (It also can be
121         // the same memory as the real grid for in-place transforms.)
122         // The buffer underlying the accessor may have a size that is
123         // larger than the active grid, because it is allocated with
124         // reallocateDeviceBuffer. The size of that larger-than-needed
125         // grid can be an odd number of floats, even though actual
126         // grid code only accesses up to an even number of floats. If
127         // we would use the reinterpet method of the accessor to
128         // convert from float to float2, runtime boundary checks can
129         // fail because of this mismatch. So, we extract the
130         // underlying global_ptr and use that to construct
131         // cl::sycl::float2 values when needed.
132         cl::sycl::global_ptr<float> gm_fourierGrid = a_fourierGrid.get_pointer();
133
134         /* Various grid sizes and indices */
135         const int localOffsetMinor = 0, localOffsetMajor = 0, localOffsetMiddle = 0;
136         const int localSizeMinor   = solveKernelParams.complexGridSizePadded[minorDim];
137         const int localSizeMiddle  = solveKernelParams.complexGridSizePadded[middleDim];
138         const int localCountMiddle = solveKernelParams.complexGridSize[middleDim];
139         const int localCountMinor  = solveKernelParams.complexGridSize[minorDim];
140         const int nMajor           = solveKernelParams.realGridSize[majorDim];
141         const int nMiddle          = solveKernelParams.realGridSize[middleDim];
142         const int nMinor           = solveKernelParams.realGridSize[minorDim];
143         const int maxkMajor        = (nMajor + 1) / 2;  // X or Y
144         const int maxkMiddle       = (nMiddle + 1) / 2; // Y OR Z => only check for !YZX
145         const int maxkMinor        = (nMinor + 1) / 2;  // Z or X => only check for YZX
146
147         const int threadLocalId     = itemIdx.get_local_linear_id();
148         const int gridLineSize      = localCountMinor;
149         const int gridLineIndex     = threadLocalId / gridLineSize;
150         const int gridLineCellIndex = threadLocalId - gridLineSize * gridLineIndex;
151         const int gridLinesPerBlock =
152                 cl::sycl::max(itemIdx.get_local_range(2) / size_t(gridLineSize), size_t(1));
153         const int activeWarps = (itemIdx.get_local_range(2) / subGroupSize);
154         const int indexMinor = itemIdx.get_group(2) * itemIdx.get_local_range(2) + gridLineCellIndex;
155         const int indexMiddle = itemIdx.get_group(1) * gridLinesPerBlock + gridLineIndex;
156         const int indexMajor  = itemIdx.get_group(0);
157
158         /* Optional outputs */
159         float energy = 0.0F;
160         float virxx  = 0.0F;
161         float virxy  = 0.0F;
162         float virxz  = 0.0F;
163         float viryy  = 0.0F;
164         float viryz  = 0.0F;
165         float virzz  = 0.0F;
166
167         assert(indexMajor < solveKernelParams.complexGridSize[majorDim]);
168         if ((indexMiddle < localCountMiddle) & (indexMinor < localCountMinor)
169             & (gridLineIndex < gridLinesPerBlock))
170         {
171             /* The offset should be equal to the global thread index for coalesced access */
172             const int gridThreadIndex =
173                     (indexMajor * localSizeMiddle + indexMiddle) * localSizeMinor + indexMinor;
174
175             const int kMajor = indexMajor + localOffsetMajor;
176             /* Checking either X in XYZ, or Y in YZX cases */
177             const float mMajor = (kMajor < maxkMajor) ? kMajor : (kMajor - nMajor);
178
179             const int kMiddle = indexMiddle + localOffsetMiddle;
180             float     mMiddle = kMiddle;
181             /* Checking Y in XYZ case */
182             if (gridOrdering == GridOrdering::XYZ)
183             {
184                 mMiddle = (kMiddle < maxkMiddle) ? kMiddle : (kMiddle - nMiddle);
185             }
186             const int kMinor = localOffsetMinor + indexMinor;
187             float     mMinor = kMinor;
188             /* Checking X in YZX case */
189             if (gridOrdering == GridOrdering::YZX)
190             {
191                 mMinor = (kMinor < maxkMinor) ? kMinor : (kMinor - nMinor);
192             }
193             /* We should skip the k-space point (0,0,0) */
194             const bool notZeroPoint = (kMinor > 0) | (kMajor > 0) | (kMiddle > 0);
195
196             float mX, mY, mZ;
197             switch (gridOrdering)
198             {
199                 case GridOrdering::YZX:
200                     mX = mMinor;
201                     mY = mMajor;
202                     mZ = mMiddle;
203                     break;
204
205                 case GridOrdering::XYZ:
206                     mX = mMajor;
207                     mY = mMiddle;
208                     mZ = mMinor;
209                     break;
210
211                 default: assert(false);
212             }
213
214             /* 0.5 correction factor for the first and last components of a Z dimension */
215             float corner_fac = 1.0F;
216             switch (gridOrdering)
217             {
218                 case GridOrdering::YZX:
219                     if ((kMiddle == 0) | (kMiddle == maxkMiddle))
220                     {
221                         corner_fac = 0.5F;
222                     }
223                     break;
224
225                 case GridOrdering::XYZ:
226                     if ((kMinor == 0) | (kMinor == maxkMinor))
227                     {
228                         corner_fac = 0.5F;
229                     }
230                     break;
231
232                 default: assert(false);
233             }
234
235             if (notZeroPoint)
236             {
237                 const float mhxk = mX * solveKernelParams.recipBox[XX][XX];
238                 const float mhyk = mX * solveKernelParams.recipBox[XX][YY]
239                                    + mY * solveKernelParams.recipBox[YY][YY];
240                 const float mhzk = mX * solveKernelParams.recipBox[XX][ZZ]
241                                    + mY * solveKernelParams.recipBox[YY][ZZ]
242                                    + mZ * solveKernelParams.recipBox[ZZ][ZZ];
243
244                 const float m2k = mhxk * mhxk + mhyk * mhyk + mhzk * mhzk;
245                 assert(m2k != 0.0F);
246                 float denom = m2k * float(M_PI) * solveKernelParams.boxVolume * gm_splineValueMajor[kMajor]
247                               * gm_splineValueMiddle[kMiddle] * gm_splineValueMinor[kMinor];
248                 assert(sycl_2020::isfinite(denom));
249                 assert(denom != 0.0F);
250
251                 const float tmp1   = cl::sycl::exp(-solveKernelParams.ewaldFactor * m2k);
252                 const float etermk = solveKernelParams.elFactor * tmp1 / denom;
253
254                 // sycl::float2::load and store are buggy in hipSYCL,
255                 // but can probably be used after resolution of
256                 // https://github.com/illuhad/hipSYCL/issues/647
257                 cl::sycl::float2 gridValue;
258                 sycl_2020::loadToVec(
259                         gridThreadIndex, cl::sycl::global_ptr<const float>(gm_fourierGrid), &gridValue);
260                 const cl::sycl::float2 oldGridValue = gridValue;
261                 gridValue *= etermk;
262                 sycl_2020::storeFromVec(gridValue, gridThreadIndex, gm_fourierGrid);
263
264                 if (computeEnergyAndVirial)
265                 {
266                     const float tmp1k = 2.0F * cl::sycl::dot(gridValue, oldGridValue);
267
268                     float vfactor = (solveKernelParams.ewaldFactor + 1.0F / m2k) * 2.0F;
269                     float ets2    = corner_fac * tmp1k;
270                     energy        = ets2;
271
272                     float ets2vf = ets2 * vfactor;
273
274                     virxx = ets2vf * mhxk * mhxk - ets2;
275                     virxy = ets2vf * mhxk * mhyk;
276                     virxz = ets2vf * mhxk * mhzk;
277                     viryy = ets2vf * mhyk * mhyk - ets2;
278                     viryz = ets2vf * mhyk * mhzk;
279                     virzz = ets2vf * mhzk * mhzk - ets2;
280                 }
281             }
282         }
283
284         /* Optional energy/virial reduction */
285         if constexpr (computeEnergyAndVirial)
286         {
287             /* A tricky shuffle reduction inspired by reduce_force_j_warp_shfl.
288              * The idea is to reduce 7 energy/virial components into a single variable (aligned by
289              * 8). We will reduce everything into virxx.
290              */
291
292             /* We can only reduce warp-wise */
293             const int width = subGroupSize;
294             static_assert(subGroupSize >= 8);
295
296             sycl_2020::sub_group sg = itemIdx.get_sub_group();
297
298             /* Making pair sums */
299             virxx += sycl_2020::shift_left(sg, virxx, 1);
300             viryy += sycl_2020::shift_right(sg, viryy, 1);
301             virzz += sycl_2020::shift_left(sg, virzz, 1);
302             virxy += sycl_2020::shift_right(sg, virxy, 1);
303             virxz += sycl_2020::shift_left(sg, virxz, 1);
304             viryz += sycl_2020::shift_right(sg, viryz, 1);
305             energy += sycl_2020::shift_left(sg, energy, 1);
306             if (threadLocalId & 1)
307             {
308                 virxx = viryy; // virxx now holds virxx and viryy pair sums
309                 virzz = virxy; // virzz now holds virzz and virxy pair sums
310                 virxz = viryz; // virxz now holds virxz and viryz pair sums
311             }
312
313             /* Making quad sums */
314             virxx += sycl_2020::shift_left(sg, virxx, 2);
315             virzz += sycl_2020::shift_right(sg, virzz, 2);
316             virxz += sycl_2020::shift_left(sg, virxz, 2);
317             energy += sycl_2020::shift_right(sg, energy, 2);
318             if (threadLocalId & 2)
319             {
320                 virxx = virzz; // virxx now holds quad sums of virxx, virxy, virzz and virxy
321                 virxz = energy; // virxz now holds quad sums of virxz, viryz, energy and unused paddings
322             }
323
324             /* Making octet sums */
325             virxx += sycl_2020::shift_left(sg, virxx, 4);
326             virxz += sycl_2020::shift_right(sg, virxz, 4);
327             if (threadLocalId & 4)
328             {
329                 virxx = virxz; // virxx now holds all 7 components' octet sums + unused paddings
330             }
331
332             /* We only need to reduce virxx now */
333 #pragma unroll
334             for (int delta = 8; delta < width; delta <<= 1)
335             {
336                 virxx += sycl_2020::shift_left(sg, virxx, delta);
337             }
338             /* Now first 7 threads of each warp have the full output contributions in virxx */
339
340             const int  componentIndex      = threadLocalId & (subGroupSize - 1);
341             const bool validComponentIndex = (componentIndex < c_virialAndEnergyCount);
342
343             if (validComponentIndex)
344             {
345                 const int warpIndex = threadLocalId / subGroupSize;
346                 sm_virialAndEnergy[warpIndex * stride + componentIndex] = virxx;
347             }
348             itemIdx.barrier(cl::sycl::access::fence_space::local_space);
349
350             /* Reduce to the single warp size */
351             const int targetIndex = threadLocalId;
352 #pragma unroll
353             for (int reductionStride = reductionBufferSize >> 1; reductionStride >= subGroupSize;
354                  reductionStride >>= 1)
355             {
356                 const int sourceIndex = targetIndex + reductionStride;
357                 if ((targetIndex < reductionStride) & (sourceIndex < activeWarps * stride))
358                 {
359                     sm_virialAndEnergy[targetIndex] += sm_virialAndEnergy[sourceIndex];
360                 }
361                 itemIdx.barrier(cl::sycl::access::fence_space::local_space);
362             }
363
364             /* Now use shuffle again */
365             /* NOTE: This reduction assumes there are at least 4 warps (asserted).
366              *       To use fewer warps, add to the conditional:
367              *       && threadLocalId < activeWarps * stride
368              */
369             assert(activeWarps * stride >= subGroupSize);
370             if (threadLocalId < subGroupSize)
371             {
372                 float output = sm_virialAndEnergy[threadLocalId];
373 #pragma unroll
374                 for (int delta = stride; delta < subGroupSize; delta <<= 1)
375                 {
376                     output += sycl_2020::shift_left(sg, output, delta);
377                 }
378                 /* Final output */
379                 if (validComponentIndex)
380                 {
381                     assert(sycl_2020::isfinite(output));
382                     atomicFetchAdd(a_virialAndEnergy[componentIndex], output);
383                 }
384             }
385         }
386     };
387 }
388
389 template<GridOrdering gridOrdering, bool computeEnergyAndVirial, int gridIndex, int subGroupSize>
390 PmeSolveKernel<gridOrdering, computeEnergyAndVirial, gridIndex, subGroupSize>::PmeSolveKernel()
391 {
392     reset();
393 }
394
395 template<GridOrdering gridOrdering, bool computeEnergyAndVirial, int gridIndex, int subGroupSize>
396 void PmeSolveKernel<gridOrdering, computeEnergyAndVirial, gridIndex, subGroupSize>::setArg(size_t argIndex,
397                                                                                            void* arg)
398 {
399     if (argIndex == 0)
400     {
401         auto* params = reinterpret_cast<PmeGpuKernelParams*>(arg);
402
403         constParams_                             = &params->constants;
404         gridParams_                              = &params->grid;
405         solveKernelParams_.ewaldFactor           = params->grid.ewaldFactor;
406         solveKernelParams_.realGridSize          = params->grid.realGridSize;
407         solveKernelParams_.complexGridSize       = params->grid.complexGridSize;
408         solveKernelParams_.complexGridSizePadded = params->grid.complexGridSizePadded;
409         solveKernelParams_.splineValuesOffset    = params->grid.splineValuesOffset;
410         solveKernelParams_.recipBox[XX]          = params->current.recipBox[XX];
411         solveKernelParams_.recipBox[YY]          = params->current.recipBox[YY];
412         solveKernelParams_.recipBox[ZZ]          = params->current.recipBox[ZZ];
413         solveKernelParams_.boxVolume             = params->current.boxVolume;
414         solveKernelParams_.elFactor              = params->constants.elFactor;
415     }
416     else
417     {
418         GMX_RELEASE_ASSERT(argIndex == 0, "Trying to pass too many args to the solve kernel");
419     }
420 }
421
422 template<GridOrdering gridOrdering, bool computeEnergyAndVirial, int gridIndex, int subGroupSize>
423 cl::sycl::event PmeSolveKernel<gridOrdering, computeEnergyAndVirial, gridIndex, subGroupSize>::launch(
424         const KernelLaunchConfig& config,
425         const DeviceStream&       deviceStream)
426 {
427     GMX_RELEASE_ASSERT(gridParams_, "Can not launch the kernel before setting its args");
428     GMX_RELEASE_ASSERT(constParams_, "Can not launch the kernel before setting its args");
429
430     using KernelNameType = PmeSolveKernel<gridOrdering, computeEnergyAndVirial, gridIndex, subGroupSize>;
431
432     // SYCL has different multidimensional layout than OpenCL/CUDA.
433     const cl::sycl::range<3> localSize{ config.blockSize[2], config.blockSize[1], config.blockSize[0] };
434     const cl::sycl::range<3> groupRange{ config.gridSize[2], config.gridSize[1], config.gridSize[0] };
435     const cl::sycl::nd_range<3> range{ groupRange * localSize, localSize };
436
437     cl::sycl::queue q = deviceStream.stream();
438
439     cl::sycl::event e = q.submit([&](cl::sycl::handler& cgh) {
440         auto kernel = makeSolveKernel<gridOrdering, computeEnergyAndVirial, subGroupSize>(
441                 cgh,
442                 gridParams_->d_splineModuli[gridIndex],
443                 solveKernelParams_,
444                 constParams_->d_virialAndEnergy[gridIndex],
445                 gridParams_->d_fourierGrid[gridIndex]);
446         cgh.parallel_for<KernelNameType>(range, kernel);
447     });
448
449     // Delete set args, so we don't forget to set them before the next launch.
450     reset();
451
452     return e;
453 }
454
455 template<GridOrdering gridOrdering, bool computeEnergyAndVirial, int gridIndex, int subGroupSize>
456 void PmeSolveKernel<gridOrdering, computeEnergyAndVirial, gridIndex, subGroupSize>::reset()
457 {
458     gridParams_  = nullptr;
459     constParams_ = nullptr;
460 }
461
462 //! Kernel class instantiations
463 /* Disable the "explicit template instantiation 'PmeSplineAndSpreadKernel<...>' will emit a vtable in every
464  * translation unit [-Wweak-template-vtables]" warning.
465  * It is only explicitly instantiated in this translation unit, so we should be safe.
466  */
467 #ifdef __clang__
468 #    pragma clang diagnostic push
469 #    pragma clang diagnostic ignored "-Wweak-template-vtables"
470 #endif
471
472 #define INSTANTIATE(subGroupSize)                                             \
473     template class PmeSolveKernel<GridOrdering::XYZ, false, 0, subGroupSize>; \
474     template class PmeSolveKernel<GridOrdering::XYZ, true, 0, subGroupSize>;  \
475     template class PmeSolveKernel<GridOrdering::YZX, false, 0, subGroupSize>; \
476     template class PmeSolveKernel<GridOrdering::YZX, true, 0, subGroupSize>;  \
477     template class PmeSolveKernel<GridOrdering::XYZ, false, 1, subGroupSize>; \
478     template class PmeSolveKernel<GridOrdering::XYZ, true, 1, subGroupSize>;  \
479     template class PmeSolveKernel<GridOrdering::YZX, false, 1, subGroupSize>; \
480     template class PmeSolveKernel<GridOrdering::YZX, true, 1, subGroupSize>;
481
482 #if GMX_SYCL_DPCPP
483 INSTANTIATE(16);
484 #elif GMX_SYCL_HIPSYCL
485 INSTANTIATE(32);
486 INSTANTIATE(64);
487 #endif
488
489 #ifdef __clang__
490 #    pragma clang diagnostic pop
491 #endif