185fea7de3420607514eb54726087591e39f63a2
[alexxy/gromacs.git] / src / gromacs / ewald / pme_gpu_internal.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020 by the GROMACS development team.
5  * Copyright (c) 2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  *
39  * \brief This file contains internal function implementations
40  * for performing the PME calculations on GPU.
41  *
42  * Note that this file is compiled as regular C++ source in OpenCL builds, but
43  * it is treated as CUDA source in CUDA-enabled GPU builds.
44  *
45  * \author Aleksei Iupinov <a.yupinov@gmail.com>
46  * \ingroup module_ewald
47  */
48
49 #include "gmxpre.h"
50
51 #include "pme_gpu_internal.h"
52
53 #include "config.h"
54
55 #include <list>
56 #include <memory>
57 #include <string>
58
59 #include "gromacs/ewald/ewald_utils.h"
60 #include "gromacs/fft/gpu_3dfft.h"
61 #include "gromacs/gpu_utils/device_context.h"
62 #include "gromacs/gpu_utils/device_stream.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/gpu_utils/pmalloc.h"
65 #if GMX_GPU_SYCL
66 #    include "gromacs/gpu_utils/syclutils.h"
67 #endif
68 #include "gromacs/hardware/device_information.h"
69 #include "gromacs/math/invertmatrix.h"
70 #include "gromacs/math/units.h"
71 #include "gromacs/timing/gpu_timing.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/logger.h"
76 #include "gromacs/utility/stringutil.h"
77 #include "gromacs/ewald/pme.h"
78
79 #if GMX_GPU_CUDA
80 #    include "pme.cuh"
81 #endif
82
83 #include "pme_gpu_calculate_splines.h"
84 #include "pme_gpu_constants.h"
85 #include "pme_gpu_program_impl.h"
86 #include "pme_gpu_timings.h"
87 #include "pme_gpu_types.h"
88 #include "pme_gpu_types_host.h"
89 #include "pme_gpu_types_host_impl.h"
90 #include "pme_grid.h"
91 #include "pme_internal.h"
92 #include "pme_solve.h"
93
94 /*! \brief
95  * CUDA only
96  * Atom limit above which it is advantageous to turn on the
97  * recalculating of the splines in the gather and using less threads per atom in the spline and spread
98  */
99 constexpr int c_pmeGpuPerformanceAtomLimit = 23000;
100
101 /*! \internal \brief
102  * Wrapper for getting a pointer to the plain C++ part of the GPU kernel parameters structure.
103  *
104  * \param[in] pmeGpu  The PME GPU structure.
105  * \returns The pointer to the kernel parameters.
106  */
107 static PmeGpuKernelParamsBase* pme_gpu_get_kernel_params_base_ptr(const PmeGpu* pmeGpu)
108 {
109     // reinterpret_cast is needed because the derived CUDA structure is not known in this file
110     auto* kernelParamsPtr = reinterpret_cast<PmeGpuKernelParamsBase*>(pmeGpu->kernelParams.get());
111     return kernelParamsPtr;
112 }
113
114 /*! \brief
115  * Atom data block size (in terms of number of atoms).
116  * This is the least common multiple of number of atoms processed by
117  * a single block/workgroup of the spread and gather kernels.
118  * The GPU atom data buffers must be padded, which means that
119  * the numbers of atoms used for determining the size of the memory
120  * allocation must be divisible by this.
121  */
122 constexpr int c_pmeAtomDataBlockSize = 64;
123
124 int pme_gpu_get_atom_data_block_size()
125 {
126     return c_pmeAtomDataBlockSize;
127 }
128
129 void pme_gpu_synchronize(const PmeGpu* pmeGpu)
130 {
131     pmeGpu->archSpecific->pmeStream_.synchronize();
132 }
133
134 void pme_gpu_alloc_energy_virial(PmeGpu* pmeGpu)
135 {
136     const size_t energyAndVirialSize = c_virialAndEnergyCount * sizeof(float);
137
138     GMX_ASSERT(
139             pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
140             "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
141
142     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
143     {
144         allocateDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy[gridIndex],
145                              c_virialAndEnergyCount,
146                              pmeGpu->archSpecific->deviceContext_);
147         pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_virialAndEnergy[gridIndex]), energyAndVirialSize);
148     }
149 }
150
151 void pme_gpu_free_energy_virial(PmeGpu* pmeGpu)
152 {
153     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
154     {
155         freeDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy[gridIndex]);
156         pfree(pmeGpu->staging.h_virialAndEnergy[gridIndex]);
157         pmeGpu->staging.h_virialAndEnergy[gridIndex] = nullptr;
158     }
159 }
160
161 void pme_gpu_clear_energy_virial(const PmeGpu* pmeGpu)
162 {
163     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
164     {
165         clearDeviceBufferAsync(&pmeGpu->kernelParams->constants.d_virialAndEnergy[gridIndex],
166                                0,
167                                c_virialAndEnergyCount,
168                                pmeGpu->archSpecific->pmeStream_);
169     }
170 }
171
172 void pme_gpu_realloc_and_copy_bspline_values(PmeGpu* pmeGpu, const int gridIndex)
173 {
174     GMX_ASSERT(
175             pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
176             "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
177     GMX_ASSERT(gridIndex < pmeGpu->common->ngrids,
178                "Invalid combination of gridIndex and number of grids");
179
180     const int splineValuesOffset[DIM] = { 0,
181                                           pmeGpu->kernelParams->grid.realGridSize[XX],
182                                           pmeGpu->kernelParams->grid.realGridSize[XX]
183                                                   + pmeGpu->kernelParams->grid.realGridSize[YY] };
184     memcpy(&pmeGpu->kernelParams->grid.splineValuesOffset, &splineValuesOffset, sizeof(splineValuesOffset));
185
186     const int newSplineValuesSize = pmeGpu->kernelParams->grid.realGridSize[XX]
187                                     + pmeGpu->kernelParams->grid.realGridSize[YY]
188                                     + pmeGpu->kernelParams->grid.realGridSize[ZZ];
189     const bool shouldRealloc = (newSplineValuesSize > pmeGpu->archSpecific->splineValuesSize[gridIndex]);
190     reallocateDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli[gridIndex],
191                            newSplineValuesSize,
192                            &pmeGpu->archSpecific->splineValuesSize[gridIndex],
193                            &pmeGpu->archSpecific->splineValuesCapacity[gridIndex],
194                            pmeGpu->archSpecific->deviceContext_);
195     if (shouldRealloc)
196     {
197         /* Reallocate the host buffer */
198         pfree(pmeGpu->staging.h_splineModuli[gridIndex]);
199         pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_splineModuli[gridIndex]),
200                 newSplineValuesSize * sizeof(float));
201     }
202     for (int i = 0; i < DIM; i++)
203     {
204         memcpy(pmeGpu->staging.h_splineModuli[gridIndex] + splineValuesOffset[i],
205                pmeGpu->common->bsp_mod[i].data(),
206                pmeGpu->common->bsp_mod[i].size() * sizeof(float));
207     }
208     /* TODO: pin original buffer instead! */
209     copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli[gridIndex],
210                        pmeGpu->staging.h_splineModuli[gridIndex],
211                        0,
212                        newSplineValuesSize,
213                        pmeGpu->archSpecific->pmeStream_,
214                        pmeGpu->settings.transferKind,
215                        nullptr);
216 }
217
218 void pme_gpu_free_bspline_values(const PmeGpu* pmeGpu)
219 {
220     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
221     {
222         pfree(pmeGpu->staging.h_splineModuli[gridIndex]);
223         freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli[gridIndex]);
224     }
225 }
226
227 void pme_gpu_realloc_forces(PmeGpu* pmeGpu)
228 {
229     const size_t newForcesSize = pmeGpu->nAtomsAlloc;
230     GMX_ASSERT(newForcesSize > 0, "Bad number of atoms in PME GPU");
231     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces,
232                            newForcesSize,
233                            &pmeGpu->archSpecific->forcesSize,
234                            &pmeGpu->archSpecific->forcesSizeAlloc,
235                            pmeGpu->archSpecific->deviceContext_);
236     pmeGpu->staging.h_forces.reserveWithPadding(pmeGpu->nAtomsAlloc);
237     pmeGpu->staging.h_forces.resizeWithPadding(pmeGpu->kernelParams->atoms.nAtoms);
238 }
239
240 void pme_gpu_free_forces(const PmeGpu* pmeGpu)
241 {
242     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces);
243 }
244
245 void pme_gpu_copy_input_forces(PmeGpu* pmeGpu)
246 {
247     GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
248     copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces,
249                        pmeGpu->staging.h_forces.data(),
250                        0,
251                        pmeGpu->kernelParams->atoms.nAtoms,
252                        pmeGpu->archSpecific->pmeStream_,
253                        pmeGpu->settings.transferKind,
254                        nullptr);
255 }
256
257 void pme_gpu_copy_output_forces(PmeGpu* pmeGpu)
258 {
259     GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
260     copyFromDeviceBuffer(pmeGpu->staging.h_forces.data(),
261                          &pmeGpu->kernelParams->atoms.d_forces,
262                          0,
263                          pmeGpu->kernelParams->atoms.nAtoms,
264                          pmeGpu->archSpecific->pmeStream_,
265                          pmeGpu->settings.transferKind,
266                          nullptr);
267 }
268
269 void pme_gpu_realloc_and_copy_input_coefficients(const PmeGpu* pmeGpu,
270                                                  const float*  h_coefficients,
271                                                  const int     gridIndex)
272 {
273     GMX_ASSERT(h_coefficients, "Bad host-side charge buffer in PME GPU");
274     const size_t newCoefficientsSize = pmeGpu->nAtomsAlloc;
275     GMX_ASSERT(newCoefficientsSize > 0, "Bad number of atoms in PME GPU");
276     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex],
277                            newCoefficientsSize,
278                            &pmeGpu->archSpecific->coefficientsSize[gridIndex],
279                            &pmeGpu->archSpecific->coefficientsCapacity[gridIndex],
280                            pmeGpu->archSpecific->deviceContext_);
281     copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex],
282                        const_cast<float*>(h_coefficients),
283                        0,
284                        pmeGpu->kernelParams->atoms.nAtoms,
285                        pmeGpu->archSpecific->pmeStream_,
286                        pmeGpu->settings.transferKind,
287                        nullptr);
288
289     const size_t paddingIndex = pmeGpu->kernelParams->atoms.nAtoms;
290     const size_t paddingCount = pmeGpu->nAtomsAlloc - paddingIndex;
291     if (paddingCount > 0)
292     {
293         clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex],
294                                paddingIndex,
295                                paddingCount,
296                                pmeGpu->archSpecific->pmeStream_);
297     }
298 }
299
300 void pme_gpu_free_coefficients(const PmeGpu* pmeGpu)
301 {
302     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
303     {
304         freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex]);
305     }
306 }
307
308 void pme_gpu_realloc_spline_data(PmeGpu* pmeGpu)
309 {
310     const int order             = pmeGpu->common->pme_order;
311     const int newSplineDataSize = DIM * order * pmeGpu->nAtomsAlloc;
312     GMX_ASSERT(newSplineDataSize > 0, "Bad number of atoms in PME GPU");
313     /* Two arrays of the same size */
314     const bool shouldRealloc        = (newSplineDataSize > pmeGpu->archSpecific->splineDataSize);
315     int        currentSizeTemp      = pmeGpu->archSpecific->splineDataSize;
316     int        currentSizeTempAlloc = pmeGpu->archSpecific->splineDataSizeAlloc;
317     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta,
318                            newSplineDataSize,
319                            &currentSizeTemp,
320                            &currentSizeTempAlloc,
321                            pmeGpu->archSpecific->deviceContext_);
322     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta,
323                            newSplineDataSize,
324                            &pmeGpu->archSpecific->splineDataSize,
325                            &pmeGpu->archSpecific->splineDataSizeAlloc,
326                            pmeGpu->archSpecific->deviceContext_);
327     // the host side reallocation
328     if (shouldRealloc)
329     {
330         pfree(pmeGpu->staging.h_theta);
331         pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_theta), newSplineDataSize * sizeof(float));
332         pfree(pmeGpu->staging.h_dtheta);
333         pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_dtheta), newSplineDataSize * sizeof(float));
334     }
335 }
336
337 void pme_gpu_free_spline_data(const PmeGpu* pmeGpu)
338 {
339     /* Two arrays of the same size */
340     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta);
341     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta);
342     pfree(pmeGpu->staging.h_theta);
343     pfree(pmeGpu->staging.h_dtheta);
344 }
345
346 void pme_gpu_realloc_grid_indices(PmeGpu* pmeGpu)
347 {
348     const size_t newIndicesSize = DIM * pmeGpu->nAtomsAlloc;
349     GMX_ASSERT(newIndicesSize > 0, "Bad number of atoms in PME GPU");
350     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices,
351                            newIndicesSize,
352                            &pmeGpu->archSpecific->gridlineIndicesSize,
353                            &pmeGpu->archSpecific->gridlineIndicesSizeAlloc,
354                            pmeGpu->archSpecific->deviceContext_);
355     pfree(pmeGpu->staging.h_gridlineIndices);
356     pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_gridlineIndices), newIndicesSize * sizeof(int));
357 }
358
359 void pme_gpu_free_grid_indices(const PmeGpu* pmeGpu)
360 {
361     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices);
362     pfree(pmeGpu->staging.h_gridlineIndices);
363 }
364
365 void pme_gpu_realloc_grids(PmeGpu* pmeGpu)
366 {
367     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
368
369     const int newRealGridSize = kernelParamsPtr->grid.realGridSizePadded[XX]
370                                 * kernelParamsPtr->grid.realGridSizePadded[YY]
371                                 * kernelParamsPtr->grid.realGridSizePadded[ZZ];
372     const int newComplexGridSize = kernelParamsPtr->grid.complexGridSizePadded[XX]
373                                    * kernelParamsPtr->grid.complexGridSizePadded[YY]
374                                    * kernelParamsPtr->grid.complexGridSizePadded[ZZ] * 2;
375     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
376     {
377         // Multiplied by 2 because we count complex grid size for complex numbers, but all allocations/pointers are float
378         if (pmeGpu->archSpecific->performOutOfPlaceFFT)
379         {
380             /* 2 separate grids */
381             reallocateDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid[gridIndex],
382                                    newComplexGridSize,
383                                    &pmeGpu->archSpecific->complexGridSize[gridIndex],
384                                    &pmeGpu->archSpecific->complexGridCapacity[gridIndex],
385                                    pmeGpu->archSpecific->deviceContext_);
386             reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid[gridIndex],
387                                    newRealGridSize,
388                                    &pmeGpu->archSpecific->realGridSize[gridIndex],
389                                    &pmeGpu->archSpecific->realGridCapacity[gridIndex],
390                                    pmeGpu->archSpecific->deviceContext_);
391         }
392         else
393         {
394             /* A single buffer so that any grid will fit */
395             const int newGridsSize = std::max(newRealGridSize, newComplexGridSize);
396             reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid[gridIndex],
397                                    newGridsSize,
398                                    &pmeGpu->archSpecific->realGridSize[gridIndex],
399                                    &pmeGpu->archSpecific->realGridCapacity[gridIndex],
400                                    pmeGpu->archSpecific->deviceContext_);
401             kernelParamsPtr->grid.d_fourierGrid[gridIndex] = kernelParamsPtr->grid.d_realGrid[gridIndex];
402             pmeGpu->archSpecific->complexGridSize[gridIndex] =
403                     pmeGpu->archSpecific->realGridSize[gridIndex];
404             // the size might get used later for copying the grid
405         }
406     }
407 }
408
409 void pme_gpu_free_grids(const PmeGpu* pmeGpu)
410 {
411     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
412     {
413         if (pmeGpu->archSpecific->performOutOfPlaceFFT)
414         {
415             freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_fourierGrid[gridIndex]);
416         }
417         freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid[gridIndex]);
418     }
419 }
420
421 void pme_gpu_clear_grids(const PmeGpu* pmeGpu)
422 {
423     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
424     {
425         clearDeviceBufferAsync(&pmeGpu->kernelParams->grid.d_realGrid[gridIndex],
426                                0,
427                                pmeGpu->archSpecific->realGridSize[gridIndex],
428                                pmeGpu->archSpecific->pmeStream_);
429     }
430 }
431
432 void pme_gpu_realloc_and_copy_fract_shifts(PmeGpu* pmeGpu)
433 {
434     pme_gpu_free_fract_shifts(pmeGpu);
435
436     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
437
438     const int nx                  = kernelParamsPtr->grid.realGridSize[XX];
439     const int ny                  = kernelParamsPtr->grid.realGridSize[YY];
440     const int nz                  = kernelParamsPtr->grid.realGridSize[ZZ];
441     const int cellCount           = c_pmeNeighborUnitcellCount;
442     const int gridDataOffset[DIM] = { 0, cellCount * nx, cellCount * (nx + ny) };
443
444     memcpy(kernelParamsPtr->grid.tablesOffsets, &gridDataOffset, sizeof(gridDataOffset));
445
446     const int newFractShiftsSize = cellCount * (nx + ny + nz);
447
448     initParamLookupTable(&kernelParamsPtr->grid.d_fractShiftsTable,
449                          &kernelParamsPtr->fractShiftsTableTexture,
450                          pmeGpu->common->fsh.data(),
451                          newFractShiftsSize,
452                          pmeGpu->archSpecific->deviceContext_);
453
454     initParamLookupTable(&kernelParamsPtr->grid.d_gridlineIndicesTable,
455                          &kernelParamsPtr->gridlineIndicesTableTexture,
456                          pmeGpu->common->nn.data(),
457                          newFractShiftsSize,
458                          pmeGpu->archSpecific->deviceContext_);
459 }
460
461 void pme_gpu_free_fract_shifts(const PmeGpu* pmeGpu)
462 {
463     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
464 #if GMX_GPU_CUDA
465     destroyParamLookupTable(&kernelParamsPtr->grid.d_fractShiftsTable,
466                             &kernelParamsPtr->fractShiftsTableTexture);
467     destroyParamLookupTable(&kernelParamsPtr->grid.d_gridlineIndicesTable,
468                             &kernelParamsPtr->gridlineIndicesTableTexture);
469 #elif GMX_GPU_OPENCL || GMX_GPU_SYCL
470     freeDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable);
471     freeDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable);
472 #endif
473 }
474
475 bool pme_gpu_stream_query(const PmeGpu* pmeGpu)
476 {
477     return haveStreamTasksCompleted(pmeGpu->archSpecific->pmeStream_);
478 }
479
480 void pme_gpu_copy_input_gather_grid(const PmeGpu* pmeGpu, const float* h_grid, const int gridIndex)
481 {
482     copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid[gridIndex],
483                        h_grid,
484                        0,
485                        pmeGpu->archSpecific->realGridSize[gridIndex],
486                        pmeGpu->archSpecific->pmeStream_,
487                        pmeGpu->settings.transferKind,
488                        nullptr);
489 }
490
491 void pme_gpu_copy_output_spread_grid(const PmeGpu* pmeGpu, float* h_grid, const int gridIndex)
492 {
493     copyFromDeviceBuffer(h_grid,
494                          &pmeGpu->kernelParams->grid.d_realGrid[gridIndex],
495                          0,
496                          pmeGpu->archSpecific->realGridSize[gridIndex],
497                          pmeGpu->archSpecific->pmeStream_,
498                          pmeGpu->settings.transferKind,
499                          nullptr);
500     pmeGpu->archSpecific->syncSpreadGridD2H.markEvent(pmeGpu->archSpecific->pmeStream_);
501 }
502
503 void pme_gpu_copy_output_spread_atom_data(const PmeGpu* pmeGpu)
504 {
505     const size_t splinesCount    = DIM * pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order;
506     auto*        kernelParamsPtr = pmeGpu->kernelParams.get();
507     copyFromDeviceBuffer(pmeGpu->staging.h_dtheta,
508                          &kernelParamsPtr->atoms.d_dtheta,
509                          0,
510                          splinesCount,
511                          pmeGpu->archSpecific->pmeStream_,
512                          pmeGpu->settings.transferKind,
513                          nullptr);
514     copyFromDeviceBuffer(pmeGpu->staging.h_theta,
515                          &kernelParamsPtr->atoms.d_theta,
516                          0,
517                          splinesCount,
518                          pmeGpu->archSpecific->pmeStream_,
519                          pmeGpu->settings.transferKind,
520                          nullptr);
521     copyFromDeviceBuffer(pmeGpu->staging.h_gridlineIndices,
522                          &kernelParamsPtr->atoms.d_gridlineIndices,
523                          0,
524                          kernelParamsPtr->atoms.nAtoms * DIM,
525                          pmeGpu->archSpecific->pmeStream_,
526                          pmeGpu->settings.transferKind,
527                          nullptr);
528 }
529
530 void pme_gpu_copy_input_gather_atom_data(const PmeGpu* pmeGpu)
531 {
532     const size_t splinesCount    = DIM * pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order;
533     auto*        kernelParamsPtr = pmeGpu->kernelParams.get();
534
535     // TODO: could clear only the padding and not the whole thing, but this is a test-exclusive code anyway
536     clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_gridlineIndices,
537                            0,
538                            pmeGpu->nAtomsAlloc * DIM,
539                            pmeGpu->archSpecific->pmeStream_);
540     clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_dtheta,
541                            0,
542                            pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM,
543                            pmeGpu->archSpecific->pmeStream_);
544     clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_theta,
545                            0,
546                            pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM,
547                            pmeGpu->archSpecific->pmeStream_);
548
549     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_dtheta,
550                        pmeGpu->staging.h_dtheta,
551                        0,
552                        splinesCount,
553                        pmeGpu->archSpecific->pmeStream_,
554                        pmeGpu->settings.transferKind,
555                        nullptr);
556     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_theta,
557                        pmeGpu->staging.h_theta,
558                        0,
559                        splinesCount,
560                        pmeGpu->archSpecific->pmeStream_,
561                        pmeGpu->settings.transferKind,
562                        nullptr);
563     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_gridlineIndices,
564                        pmeGpu->staging.h_gridlineIndices,
565                        0,
566                        kernelParamsPtr->atoms.nAtoms * DIM,
567                        pmeGpu->archSpecific->pmeStream_,
568                        pmeGpu->settings.transferKind,
569                        nullptr);
570 }
571
572 void pme_gpu_sync_spread_grid(const PmeGpu* pmeGpu)
573 {
574     pmeGpu->archSpecific->syncSpreadGridD2H.waitForEvent();
575 }
576
577 /*! \brief Internal GPU initialization for PME.
578  *
579  * \param[in]  pmeGpu         GPU PME data.
580  * \param[in]  deviceContext  GPU context.
581  * \param[in]  deviceStream   GPU stream.
582  */
583 static void pme_gpu_init_internal(PmeGpu* pmeGpu, const DeviceContext& deviceContext, const DeviceStream& deviceStream)
584 {
585     /* Allocate the target-specific structures */
586     pmeGpu->archSpecific.reset(new PmeGpuSpecific(deviceContext, deviceStream));
587     pmeGpu->kernelParams.reset(new PmeGpuKernelParams());
588
589     pmeGpu->archSpecific->performOutOfPlaceFFT = true;
590     /* This should give better performance, according to the cuFFT documentation.
591      * The performance seems to be the same though.
592      * TODO: PME could also try to pick up nice grid sizes (with factors of 2, 3, 5, 7).
593      */
594
595 #if GMX_GPU_CUDA
596     pmeGpu->maxGridWidthX = deviceContext.deviceInfo().prop.maxGridSize[0];
597 #else
598     // Use this path for any non-CUDA GPU acceleration
599     // TODO: is there no really global work size limit in OpenCL?
600     pmeGpu->maxGridWidthX = INT32_MAX / 2;
601 #endif
602 }
603
604 void pme_gpu_reinit_3dfft(const PmeGpu* pmeGpu)
605 {
606     if (pme_gpu_settings(pmeGpu).performGPUFFT)
607     {
608         pmeGpu->archSpecific->fftSetup.resize(0);
609         const bool         performOutOfPlaceFFT      = pmeGpu->archSpecific->performOutOfPlaceFFT;
610         const bool         allocateGrid              = false;
611         MPI_Comm           comm                      = MPI_COMM_NULL;
612         std::array<int, 1> gridOffsetsInXForEachRank = { 0 };
613         std::array<int, 1> gridOffsetsInYForEachRank = { 0 };
614 #if GMX_GPU_CUDA
615         const gmx::FftBackend backend = gmx::FftBackend::Cufft;
616 #elif GMX_GPU_OPENCL
617         const gmx::FftBackend backend = gmx::FftBackend::Ocl;
618 #elif GMX_GPU_SYCL
619         const gmx::FftBackend backend = gmx::FftBackend::Sycl;
620 #else
621         GMX_RELEASE_ASSERT(false, "Unknown GPU backend");
622         const gmx::FftBackend backend = gmx::FftBackend::Count;
623 #endif
624
625         PmeGpuGridParams& grid = pme_gpu_get_kernel_params_base_ptr(pmeGpu)->grid;
626         for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
627         {
628             pmeGpu->archSpecific->fftSetup.push_back(
629                     std::make_unique<gmx::Gpu3dFft>(backend,
630                                                     allocateGrid,
631                                                     comm,
632                                                     gridOffsetsInXForEachRank,
633                                                     gridOffsetsInYForEachRank,
634                                                     grid.realGridSize[ZZ],
635                                                     performOutOfPlaceFFT,
636                                                     pmeGpu->archSpecific->deviceContext_,
637                                                     pmeGpu->archSpecific->pmeStream_,
638                                                     grid.realGridSize,
639                                                     grid.realGridSizePadded,
640                                                     grid.complexGridSizePadded,
641                                                     &(grid.d_realGrid[gridIndex]),
642                                                     &(grid.d_fourierGrid[gridIndex])));
643         }
644     }
645 }
646
647 void pme_gpu_destroy_3dfft(const PmeGpu* pmeGpu)
648 {
649     pmeGpu->archSpecific->fftSetup.resize(0);
650 }
651
652 void pme_gpu_getEnergyAndVirial(const gmx_pme_t& pme, const float lambda, PmeOutput* output)
653 {
654     const PmeGpu* pmeGpu = pme.gpu;
655
656     GMX_ASSERT(lambda == 1.0 || pmeGpu->common->ngrids == 2,
657                "Invalid combination of lambda and number of grids");
658
659     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
660     {
661         for (int j = 0; j < c_virialAndEnergyCount; j++)
662         {
663             GMX_ASSERT(std::isfinite(pmeGpu->staging.h_virialAndEnergy[gridIndex][j]),
664                        "PME GPU produces incorrect energy/virial.");
665         }
666     }
667     for (int dim1 = 0; dim1 < DIM; dim1++)
668     {
669         for (int dim2 = 0; dim2 < DIM; dim2++)
670         {
671             output->coulombVirial_[dim1][dim2] = 0;
672         }
673     }
674     output->coulombEnergy_ = 0;
675     float scale            = 1.0;
676     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
677     {
678         if (pmeGpu->common->ngrids == 2)
679         {
680             scale = gridIndex == 0 ? (1.0 - lambda) : lambda;
681         }
682         output->coulombVirial_[XX][XX] +=
683                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][0];
684         output->coulombVirial_[YY][YY] +=
685                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][1];
686         output->coulombVirial_[ZZ][ZZ] +=
687                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][2];
688         output->coulombVirial_[XX][YY] +=
689                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][3];
690         output->coulombVirial_[YY][XX] +=
691                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][3];
692         output->coulombVirial_[XX][ZZ] +=
693                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][4];
694         output->coulombVirial_[ZZ][XX] +=
695                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][4];
696         output->coulombVirial_[YY][ZZ] +=
697                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][5];
698         output->coulombVirial_[ZZ][YY] +=
699                 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][5];
700         output->coulombEnergy_ += scale * 0.5F * pmeGpu->staging.h_virialAndEnergy[gridIndex][6];
701     }
702     if (pmeGpu->common->ngrids > 1)
703     {
704         output->coulombDvdl_ = 0.5F
705                                * (pmeGpu->staging.h_virialAndEnergy[FEP_STATE_B][6]
706                                   - pmeGpu->staging.h_virialAndEnergy[FEP_STATE_A][6]);
707     }
708 }
709
710 /*! \brief Sets the force-related members in \p output
711  *
712  * \param[in]   pmeGpu      PME GPU data structure
713  * \param[out]  output      Pointer to PME output data structure
714  */
715 static void pme_gpu_getForceOutput(PmeGpu* pmeGpu, PmeOutput* output)
716 {
717     output->haveForceOutput_ = !pmeGpu->settings.useGpuForceReduction;
718     if (output->haveForceOutput_)
719     {
720         output->forces_ = pmeGpu->staging.h_forces;
721     }
722 }
723
724 PmeOutput pme_gpu_getOutput(const gmx_pme_t& pme, const bool computeEnergyAndVirial, const real lambdaQ)
725 {
726     PmeGpu* pmeGpu = pme.gpu;
727
728     PmeOutput output;
729
730     pme_gpu_getForceOutput(pmeGpu, &output);
731
732     if (computeEnergyAndVirial)
733     {
734         if (pme_gpu_settings(pmeGpu).performGPUSolve)
735         {
736             pme_gpu_getEnergyAndVirial(pme, lambdaQ, &output);
737         }
738         else
739         {
740             get_pme_ener_vir_q(pme.solve_work, pme.nthread, &output);
741         }
742     }
743     return output;
744 }
745
746 void pme_gpu_update_input_box(PmeGpu gmx_unused* pmeGpu, const matrix gmx_unused box)
747 {
748 #if GMX_DOUBLE
749     GMX_THROW(gmx::NotImplementedError("PME is implemented for single-precision only on GPU"));
750 #else
751     matrix scaledBox;
752     pmeGpu->common->boxScaler->scaleBox(box, scaledBox);
753     auto* kernelParamsPtr              = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
754     kernelParamsPtr->current.boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
755     GMX_ASSERT(kernelParamsPtr->current.boxVolume != 0.0F, "Zero volume of the unit cell");
756     matrix recipBox;
757     gmx::invertBoxMatrix(scaledBox, recipBox);
758
759     /* The GPU recipBox is transposed as compared to the CPU recipBox.
760      * Spread uses matrix columns (while solve and gather use rows).
761      * There is no particular reason for this; it might be further rethought/optimized for better access patterns.
762      */
763     const real newRecipBox[DIM][DIM] = { { recipBox[XX][XX], recipBox[YY][XX], recipBox[ZZ][XX] },
764                                          { 0.0, recipBox[YY][YY], recipBox[ZZ][YY] },
765                                          { 0.0, 0.0, recipBox[ZZ][ZZ] } };
766     memcpy(kernelParamsPtr->current.recipBox, newRecipBox, sizeof(matrix));
767 #endif
768 }
769
770 /*! \brief \libinternal
771  * (Re-)initializes all the PME GPU data related to the grid size and cut-off.
772  *
773  * \param[in] pmeGpu            The PME GPU structure.
774  */
775 static void pme_gpu_reinit_grids(PmeGpu* pmeGpu)
776 {
777     auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
778
779     GMX_ASSERT(
780             pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
781             "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
782
783     kernelParamsPtr->grid.ewaldFactor =
784             (M_PI * M_PI) / (pmeGpu->common->ewaldcoeff_q * pmeGpu->common->ewaldcoeff_q);
785     /* The grid size variants */
786     for (int i = 0; i < DIM; i++)
787     {
788         kernelParamsPtr->grid.realGridSize[i] = pmeGpu->common->nk[i];
789         kernelParamsPtr->grid.realGridSizeFP[i] =
790                 static_cast<float>(kernelParamsPtr->grid.realGridSize[i]);
791         kernelParamsPtr->grid.realGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
792
793         // The complex grid currently uses no padding;
794         // if it starts to do so, then another test should be added for that
795         kernelParamsPtr->grid.complexGridSize[i]       = kernelParamsPtr->grid.realGridSize[i];
796         kernelParamsPtr->grid.complexGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
797     }
798     /* FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
799     if (!pme_gpu_settings(pmeGpu).performGPUFFT)
800     {
801         // This allows for GPU spreading grid and CPU fftgrid to have the same layout, so that we can copy the data directly
802         kernelParamsPtr->grid.realGridSizePadded[ZZ] =
803                 (kernelParamsPtr->grid.realGridSize[ZZ] / 2 + 1) * 2;
804     }
805     /* GPU FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
806     kernelParamsPtr->grid.complexGridSize[ZZ] /= 2;
807     kernelParamsPtr->grid.complexGridSize[ZZ]++;
808     kernelParamsPtr->grid.complexGridSizePadded[ZZ] = kernelParamsPtr->grid.complexGridSize[ZZ];
809
810     pme_gpu_realloc_and_copy_fract_shifts(pmeGpu);
811     for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
812     {
813         pme_gpu_realloc_and_copy_bspline_values(pmeGpu, gridIndex);
814     }
815
816     pme_gpu_realloc_grids(pmeGpu);
817     pme_gpu_reinit_3dfft(pmeGpu);
818 }
819
820 /* Several GPU functions that refer to the CPU PME data live here.
821  * We would like to keep these away from the GPU-framework specific code for clarity,
822  * as well as compilation issues with MPI.
823  */
824
825 /*! \brief \libinternal
826  * Copies everything useful from the PME CPU to the PME GPU structure.
827  * The goal is to minimize interaction with the PME CPU structure in the GPU code.
828  *
829  * \param[in] pme         The PME structure.
830  */
831 static void pme_gpu_copy_common_data_from(const gmx_pme_t* pme)
832 {
833     /* TODO: Consider refactoring the CPU PME code to use the same structure,
834      * so that this function becomes 2 lines */
835     PmeGpu* pmeGpu               = pme->gpu;
836     pmeGpu->common->ngrids       = pme->bFEP_q ? 2 : 1;
837     pmeGpu->common->epsilon_r    = pme->epsilon_r;
838     pmeGpu->common->ewaldcoeff_q = pme->ewaldcoeff_q;
839     pmeGpu->common->nk[XX]       = pme->nkx;
840     pmeGpu->common->nk[YY]       = pme->nky;
841     pmeGpu->common->nk[ZZ]       = pme->nkz;
842     pmeGpu->common->pme_order    = pme->pme_order;
843     if (pmeGpu->common->pme_order != c_pmeGpuOrder)
844     {
845         GMX_THROW(gmx::NotImplementedError("pme_order != 4 is not implemented!"));
846     }
847     for (int i = 0; i < DIM; i++)
848     {
849         pmeGpu->common->bsp_mod[i].assign(pme->bsp_mod[i], pme->bsp_mod[i] + pmeGpu->common->nk[i]);
850     }
851     const int cellCount = c_pmeNeighborUnitcellCount;
852     pmeGpu->common->fsh.resize(0);
853     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshx, pme->fshx + cellCount * pme->nkx);
854     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshy, pme->fshy + cellCount * pme->nky);
855     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshz, pme->fshz + cellCount * pme->nkz);
856     pmeGpu->common->nn.resize(0);
857     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnx, pme->nnx + cellCount * pme->nkx);
858     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nny, pme->nny + cellCount * pme->nky);
859     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnz, pme->nnz + cellCount * pme->nkz);
860     pmeGpu->common->runMode       = pme->runMode;
861     pmeGpu->common->isRankPmeOnly = !pme->bPPnode;
862     pmeGpu->common->boxScaler     = pme->boxScaler.get();
863 }
864
865 /*! \libinternal \brief
866  * uses heuristics to select the best performing PME gather and scatter kernels
867  *
868  * \param[in,out] pmeGpu         The PME GPU structure.
869  */
870 static void pme_gpu_select_best_performing_pme_spreadgather_kernels(PmeGpu* pmeGpu)
871 {
872     if (GMX_GPU_CUDA && pmeGpu->kernelParams->atoms.nAtoms > c_pmeGpuPerformanceAtomLimit)
873     {
874         pmeGpu->settings.threadsPerAtom     = ThreadsPerAtom::Order;
875         pmeGpu->settings.recalculateSplines = true;
876     }
877     else
878     {
879         pmeGpu->settings.threadsPerAtom     = ThreadsPerAtom::OrderSquared;
880         pmeGpu->settings.recalculateSplines = false;
881     }
882 }
883
884
885 /*! \libinternal \brief
886  * Initializes the PME GPU data at the beginning of the run.
887  * TODO: this should become PmeGpu::PmeGpu()
888  *
889  * \param[in,out] pme            The PME structure.
890  * \param[in]     deviceContext  The GPU context.
891  * \param[in]     deviceStream   The GPU stream.
892  * \param[in,out] pmeGpuProgram  The handle to the program/kernel data created outside (e.g. in unit tests/runner)
893  */
894 static void pme_gpu_init(gmx_pme_t*           pme,
895                          const DeviceContext& deviceContext,
896                          const DeviceStream&  deviceStream,
897                          const PmeGpuProgram* pmeGpuProgram)
898 {
899     pme->gpu       = new PmeGpu();
900     PmeGpu* pmeGpu = pme->gpu;
901     changePinningPolicy(&pmeGpu->staging.h_forces, pme_get_pinning_policy());
902     pmeGpu->common = std::make_shared<PmeShared>();
903
904     /* These settings are set here for the whole run; dynamic ones are set in pme_gpu_reinit() */
905     /* A convenience variable. */
906     pmeGpu->settings.useDecomposition = (pme->nnodes != 1);
907     /* TODO: CPU gather with GPU spread is broken due to different theta/dtheta layout. */
908     pmeGpu->settings.performGPUGather = true;
909     // By default GPU-side reduction is off (explicitly set here for tests, otherwise reset per-step)
910     pmeGpu->settings.useGpuForceReduction = false;
911
912     pme_gpu_set_testing(pmeGpu, false);
913
914     GMX_ASSERT(pmeGpuProgram != nullptr, "GPU kernels must be already compiled");
915     pmeGpu->programHandle_ = pmeGpuProgram;
916
917     pmeGpu->initializedClfftLibrary_ = std::make_unique<gmx::ClfftInitializer>();
918
919     pme_gpu_init_internal(pmeGpu, deviceContext, deviceStream);
920
921     pme_gpu_copy_common_data_from(pme);
922     pme_gpu_alloc_energy_virial(pmeGpu);
923
924     GMX_ASSERT(pmeGpu->common->epsilon_r != 0.0F, "PME GPU: bad electrostatic coefficient");
925
926     auto* kernelParamsPtr               = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
927     kernelParamsPtr->constants.elFactor = gmx::c_one4PiEps0 / pmeGpu->common->epsilon_r;
928 }
929
930 void pme_gpu_get_real_grid_sizes(const PmeGpu* pmeGpu, gmx::IVec* gridSize, gmx::IVec* paddedGridSize)
931 {
932     GMX_ASSERT(gridSize != nullptr, "");
933     GMX_ASSERT(paddedGridSize != nullptr, "");
934     GMX_ASSERT(pmeGpu != nullptr, "");
935     auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
936     for (int i = 0; i < DIM; i++)
937     {
938         (*gridSize)[i]       = kernelParamsPtr->grid.realGridSize[i];
939         (*paddedGridSize)[i] = kernelParamsPtr->grid.realGridSizePadded[i];
940     }
941 }
942
943 void pme_gpu_reinit(gmx_pme_t*           pme,
944                     const DeviceContext* deviceContext,
945                     const DeviceStream*  deviceStream,
946                     const PmeGpuProgram* pmeGpuProgram)
947 {
948     GMX_ASSERT(pme != nullptr, "Need valid PME object");
949
950     if (!pme->gpu)
951     {
952         GMX_RELEASE_ASSERT(deviceContext != nullptr,
953                            "Device context can not be nullptr when setting up PME on GPU.");
954         GMX_RELEASE_ASSERT(deviceStream != nullptr,
955                            "Device stream can not be nullptr when setting up PME on GPU.");
956         /* First-time initialization */
957         pme_gpu_init(pme, *deviceContext, *deviceStream, pmeGpuProgram);
958     }
959     else
960     {
961         /* After this call nothing in the GPU code should refer to the gmx_pme_t *pme itself - until the next pme_gpu_reinit */
962         pme_gpu_copy_common_data_from(pme);
963     }
964     /* GPU FFT will only get used for a single rank.*/
965     pme->gpu->settings.performGPUFFT =
966             (pme->gpu->common->runMode == PmeRunMode::GPU) && !pme->gpu->settings.useDecomposition;
967     pme->gpu->settings.performGPUSolve = (pme->gpu->common->runMode == PmeRunMode::GPU);
968
969     /* Reinit active timers */
970     pme_gpu_reinit_timings(pme->gpu);
971
972     pme_gpu_reinit_grids(pme->gpu);
973     // Note: if timing the reinit launch overhead becomes more relevant
974     // (e.g. with regulat PP-PME re-balancing), we should pass wcycle here.
975     pme_gpu_reinit_computation(pme, nullptr);
976     /* Clear the previous box - doesn't hurt, and forces the PME CPU recipbox
977      * update for mixed mode on grid switch. TODO: use shared recipbox field.
978      */
979     std::memset(pme->gpu->common->previousBox, 0, sizeof(pme->gpu->common->previousBox));
980 }
981
982 void pme_gpu_destroy(PmeGpu* pmeGpu)
983 {
984     /* Free lots of data */
985     pme_gpu_free_energy_virial(pmeGpu);
986     pme_gpu_free_bspline_values(pmeGpu);
987     pme_gpu_free_forces(pmeGpu);
988     pme_gpu_free_coefficients(pmeGpu);
989     pme_gpu_free_spline_data(pmeGpu);
990     pme_gpu_free_grid_indices(pmeGpu);
991     pme_gpu_free_fract_shifts(pmeGpu);
992     pme_gpu_free_grids(pmeGpu);
993
994     pme_gpu_destroy_3dfft(pmeGpu);
995
996     delete pmeGpu;
997 }
998
999 void pme_gpu_reinit_atoms(PmeGpu* pmeGpu, const int nAtoms, const real* chargesA, const real* chargesB)
1000 {
1001     auto* kernelParamsPtr         = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
1002     kernelParamsPtr->atoms.nAtoms = nAtoms;
1003     const int  block_size         = pme_gpu_get_atom_data_block_size();
1004     const int  nAtomsNewPadded    = ((nAtoms + block_size - 1) / block_size) * block_size;
1005     const bool haveToRealloc      = (pmeGpu->nAtomsAlloc < nAtomsNewPadded);
1006     pmeGpu->nAtomsAlloc           = nAtomsNewPadded;
1007
1008 #if GMX_DOUBLE
1009     GMX_RELEASE_ASSERT(false, "Only single precision supported");
1010     GMX_UNUSED_VALUE(charges);
1011 #else
1012     int gridIndex = 0;
1013     /* Could also be checked for haveToRealloc, but the copy always needs to be performed */
1014     pme_gpu_realloc_and_copy_input_coefficients(pmeGpu, reinterpret_cast<const float*>(chargesA), gridIndex);
1015     gridIndex++;
1016     if (chargesB != nullptr)
1017     {
1018         pme_gpu_realloc_and_copy_input_coefficients(
1019                 pmeGpu, reinterpret_cast<const float*>(chargesB), gridIndex);
1020     }
1021     else
1022     {
1023         /* Fill the second set of coefficients with chargesA as well to be able to avoid
1024          * conditionals in the GPU kernels */
1025         /* FIXME: This should be avoided by making a separate templated version of the
1026          * relevant kernel(s) (probably only pme_gather_kernel). That would require a
1027          * reduction of the current number of templated parameters of that kernel. */
1028         pme_gpu_realloc_and_copy_input_coefficients(
1029                 pmeGpu, reinterpret_cast<const float*>(chargesA), gridIndex);
1030     }
1031 #endif
1032
1033     if (haveToRealloc)
1034     {
1035         pme_gpu_realloc_forces(pmeGpu);
1036         pme_gpu_realloc_spline_data(pmeGpu);
1037         pme_gpu_realloc_grid_indices(pmeGpu);
1038     }
1039     pme_gpu_select_best_performing_pme_spreadgather_kernels(pmeGpu);
1040 }
1041
1042 /*! \internal \brief
1043  * Returns raw timing event from the corresponding GpuRegionTimer (if timings are enabled).
1044  * In CUDA result can be nullptr stub, per GpuRegionTimer implementation.
1045  *
1046  * \param[in] pmeGpu         The PME GPU data structure.
1047  * \param[in] pmeStageId     The PME GPU stage gtPME_ index from the enum in src/gromacs/timing/gpu_timing.h
1048  */
1049 static CommandEvent* pme_gpu_fetch_timing_event(const PmeGpu* pmeGpu, PmeStage pmeStageId)
1050 {
1051     CommandEvent* timingEvent = nullptr;
1052     if (pme_gpu_timings_enabled(pmeGpu))
1053     {
1054         GMX_ASSERT(pmeStageId < PmeStage::Count, "Wrong PME GPU timing event index");
1055         timingEvent = pmeGpu->archSpecific->timingEvents[pmeStageId].fetchNextEvent();
1056     }
1057     return timingEvent;
1058 }
1059
1060 void pme_gpu_3dfft(const PmeGpu* pmeGpu, gmx_fft_direction dir, const int grid_index)
1061 {
1062     PmeStage timerId = (dir == GMX_FFT_REAL_TO_COMPLEX) ? PmeStage::FftTransformR2C
1063                                                         : PmeStage::FftTransformC2R;
1064
1065     pme_gpu_start_timing(pmeGpu, timerId);
1066     pmeGpu->archSpecific->fftSetup[grid_index]->perform3dFft(
1067             dir, pme_gpu_fetch_timing_event(pmeGpu, timerId));
1068     pme_gpu_stop_timing(pmeGpu, timerId);
1069 }
1070
1071 /*! \brief
1072  * Given possibly large \p blockCount, returns a compact 1D or 2D grid for kernel scheduling,
1073  * to minimize number of unused blocks.
1074  */
1075 std::pair<int, int> inline pmeGpuCreateGrid(const PmeGpu* pmeGpu, int blockCount)
1076 {
1077     // How many maximum widths in X do we need (hopefully just one)
1078     const int minRowCount = (blockCount + pmeGpu->maxGridWidthX - 1) / pmeGpu->maxGridWidthX;
1079     // Trying to make things even
1080     const int colCount = (blockCount + minRowCount - 1) / minRowCount;
1081     GMX_ASSERT((colCount * minRowCount - blockCount) >= 0, "pmeGpuCreateGrid: totally wrong");
1082     GMX_ASSERT((colCount * minRowCount - blockCount) < minRowCount,
1083                "pmeGpuCreateGrid: excessive blocks");
1084     return std::pair<int, int>(colCount, minRowCount);
1085 }
1086
1087 /*! \brief
1088  * Returns a pointer to appropriate spline and spread kernel based on the input bool values
1089  *
1090  * \param[in]  pmeGpu                   The PME GPU structure.
1091  * \param[in]  threadsPerAtom           Controls whether we should use order or order*order threads per atom
1092  * \param[in]  writeSplinesToGlobal     bool controlling if we should write spline data to global memory
1093  * \param[in]  numGrids                 Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1094  *
1095  * \return Pointer to CUDA kernel
1096  */
1097 static auto selectSplineAndSpreadKernelPtr(const PmeGpu*  pmeGpu,
1098                                            ThreadsPerAtom threadsPerAtom,
1099                                            bool           writeSplinesToGlobal,
1100                                            const int      numGrids)
1101 {
1102     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1103     if (writeSplinesToGlobal)
1104     {
1105         if (threadsPerAtom == ThreadsPerAtom::Order)
1106         {
1107             if (numGrids == 2)
1108             {
1109                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesThPerAtom4Dual;
1110             }
1111             else
1112             {
1113                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesThPerAtom4Single;
1114             }
1115         }
1116         else
1117         {
1118             if (numGrids == 2)
1119             {
1120                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesDual;
1121             }
1122             else
1123             {
1124                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesSingle;
1125             }
1126         }
1127     }
1128     else
1129     {
1130         if (threadsPerAtom == ThreadsPerAtom::Order)
1131         {
1132             if (numGrids == 2)
1133             {
1134                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Dual;
1135             }
1136             else
1137             {
1138                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Single;
1139             }
1140         }
1141         else
1142         {
1143             if (numGrids == 2)
1144             {
1145                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelDual;
1146             }
1147             else
1148             {
1149                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelSingle;
1150             }
1151         }
1152     }
1153
1154     return kernelPtr;
1155 }
1156
1157 /*! \brief
1158  * Returns a pointer to appropriate spline kernel based on the input bool values
1159  *
1160  * \param[in]  pmeGpu                   The PME GPU structure.
1161  * \param[in]  threadsPerAtom           Controls whether we should use order or order*order threads per atom
1162  * \param[in]  writeSplinesToGlobal     bool controlling if we should write spline data to global memory
1163  * \param[in]  numGrids                 Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1164  *
1165  * \return Pointer to CUDA kernel
1166  */
1167 static auto selectSplineKernelPtr(const PmeGpu*   pmeGpu,
1168                                   ThreadsPerAtom  threadsPerAtom,
1169                                   bool gmx_unused writeSplinesToGlobal,
1170                                   const int       numGrids)
1171 {
1172     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1173     GMX_ASSERT(
1174             writeSplinesToGlobal,
1175             "Spline data should always be written to global memory when just calculating splines");
1176
1177     if (threadsPerAtom == ThreadsPerAtom::Order)
1178     {
1179         if (numGrids == 2)
1180         {
1181             kernelPtr = pmeGpu->programHandle_->impl_->splineKernelThPerAtom4Dual;
1182         }
1183         else
1184         {
1185             kernelPtr = pmeGpu->programHandle_->impl_->splineKernelThPerAtom4Single;
1186         }
1187     }
1188     else
1189     {
1190         if (numGrids == 2)
1191         {
1192             kernelPtr = pmeGpu->programHandle_->impl_->splineKernelDual;
1193         }
1194         else
1195         {
1196             kernelPtr = pmeGpu->programHandle_->impl_->splineKernelSingle;
1197         }
1198     }
1199     return kernelPtr;
1200 }
1201
1202 /*! \brief
1203  * Returns a pointer to appropriate spread kernel based on the input bool values
1204  *
1205  * \param[in]  pmeGpu                   The PME GPU structure.
1206  * \param[in]  threadsPerAtom           Controls whether we should use order or order*order threads per atom
1207  * \param[in]  writeSplinesToGlobal     bool controlling if we should write spline data to global memory
1208  * \param[in]  numGrids                 Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1209  *
1210  * \return Pointer to CUDA kernel
1211  */
1212 static auto selectSpreadKernelPtr(const PmeGpu*  pmeGpu,
1213                                   ThreadsPerAtom threadsPerAtom,
1214                                   bool           writeSplinesToGlobal,
1215                                   const int      numGrids)
1216 {
1217     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1218     if (writeSplinesToGlobal)
1219     {
1220         if (threadsPerAtom == ThreadsPerAtom::Order)
1221         {
1222             if (numGrids == 2)
1223             {
1224                 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelThPerAtom4Dual;
1225             }
1226             else
1227             {
1228                 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelThPerAtom4Single;
1229             }
1230         }
1231         else
1232         {
1233             if (numGrids == 2)
1234             {
1235                 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelDual;
1236             }
1237             else
1238             {
1239                 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelSingle;
1240             }
1241         }
1242     }
1243     else
1244     {
1245         /* if we are not saving the spline data we need to recalculate it
1246            using the spline and spread Kernel */
1247         if (threadsPerAtom == ThreadsPerAtom::Order)
1248         {
1249             if (numGrids == 2)
1250             {
1251                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Dual;
1252             }
1253             else
1254             {
1255                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Single;
1256             }
1257         }
1258         else
1259         {
1260             if (numGrids == 2)
1261             {
1262                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelDual;
1263             }
1264             else
1265             {
1266                 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelSingle;
1267             }
1268         }
1269     }
1270     return kernelPtr;
1271 }
1272
1273 void pme_gpu_spread(const PmeGpu*         pmeGpu,
1274                     GpuEventSynchronizer* xReadyOnDevice,
1275                     real**                h_grids,
1276                     bool                  computeSplines,
1277                     bool                  spreadCharges,
1278                     const real            lambda)
1279 {
1280     GMX_ASSERT(
1281             pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
1282             "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
1283
1284     GMX_ASSERT(computeSplines || spreadCharges,
1285                "PME spline/spread kernel has invalid input (nothing to do)");
1286     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1287     GMX_ASSERT(kernelParamsPtr->atoms.nAtoms > 0, "No atom data in PME GPU spread");
1288
1289     const size_t blockSize = pmeGpu->programHandle_->impl_->spreadWorkGroupSize;
1290
1291     const int order = pmeGpu->common->pme_order;
1292     GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1293     const bool writeGlobal = pmeGpu->settings.copyAllOutputs;
1294     const int  threadsPerAtom =
1295             (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
1296     const bool recalculateSplines = pmeGpu->settings.recalculateSplines;
1297
1298     GMX_ASSERT(!GMX_GPU_OPENCL || pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::OrderSquared,
1299                "Only 16 threads per atom supported in OpenCL");
1300     GMX_ASSERT(!GMX_GPU_OPENCL || !recalculateSplines,
1301                "Recalculating splines not supported in OpenCL");
1302
1303     const int atomsPerBlock = blockSize / threadsPerAtom;
1304
1305     // TODO: pick smaller block size in runtime if needed
1306     // (e.g. on 660 Ti where 50% occupancy is ~25% faster than 100% occupancy with RNAse (~17.8k atoms))
1307     // If doing so, change atomsPerBlock in the kernels as well.
1308     // TODO: test varying block sizes on modern arch-s as well
1309     // TODO: also consider using cudaFuncSetCacheConfig() for preferring shared memory on older architectures
1310     //(for spline data mostly)
1311     GMX_ASSERT(!(c_pmeAtomDataBlockSize % atomsPerBlock),
1312                "inconsistent atom data padding vs. spreading block size");
1313
1314     // Ensure that coordinates are ready on the device before launching spread;
1315     // only needed with CUDA on PP+PME ranks, not on separate PME ranks, in unit tests
1316     // nor in OpenCL as these cases use a single stream (hence xReadyOnDevice == nullptr).
1317     GMX_ASSERT(!GMX_GPU_CUDA || xReadyOnDevice != nullptr || pmeGpu->common->isRankPmeOnly
1318                        || pme_gpu_settings(pmeGpu).copyAllOutputs,
1319                "Need a valid coordinate synchronizer on PP+PME ranks with CUDA.");
1320
1321     if (xReadyOnDevice)
1322     {
1323         xReadyOnDevice->enqueueWaitEvent(pmeGpu->archSpecific->pmeStream_);
1324     }
1325
1326     const int blockCount = pmeGpu->nAtomsAlloc / atomsPerBlock;
1327     auto      dimGrid    = pmeGpuCreateGrid(pmeGpu, blockCount);
1328
1329     if (pmeGpu->common->ngrids == 1)
1330     {
1331         kernelParamsPtr->current.scale = 1.0;
1332     }
1333     else
1334     {
1335         kernelParamsPtr->current.scale = 1.0 - lambda;
1336     }
1337
1338     KernelLaunchConfig config;
1339     config.blockSize[0] = order;
1340     config.blockSize[1] = (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? 1 : order);
1341     config.blockSize[2] = atomsPerBlock;
1342     config.gridSize[0]  = dimGrid.first;
1343     config.gridSize[1]  = dimGrid.second;
1344
1345     PmeStage                           timingId;
1346     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1347     if (computeSplines)
1348     {
1349         if (spreadCharges)
1350         {
1351             timingId  = PmeStage::SplineAndSpread;
1352             kernelPtr = selectSplineAndSpreadKernelPtr(pmeGpu,
1353                                                        pmeGpu->settings.threadsPerAtom,
1354                                                        writeGlobal || (!recalculateSplines),
1355                                                        pmeGpu->common->ngrids);
1356         }
1357         else
1358         {
1359             timingId  = PmeStage::Spline;
1360             kernelPtr = selectSplineKernelPtr(pmeGpu,
1361                                               pmeGpu->settings.threadsPerAtom,
1362                                               writeGlobal || (!recalculateSplines),
1363                                               pmeGpu->common->ngrids);
1364         }
1365     }
1366     else
1367     {
1368         timingId  = PmeStage::Spread;
1369         kernelPtr = selectSpreadKernelPtr(pmeGpu,
1370                                           pmeGpu->settings.threadsPerAtom,
1371                                           writeGlobal || (!recalculateSplines),
1372                                           pmeGpu->common->ngrids);
1373     }
1374
1375
1376     pme_gpu_start_timing(pmeGpu, timingId);
1377     auto* timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1378 #if c_canEmbedBuffers
1379     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1380 #else
1381     const auto kernelArgs =
1382             prepareGpuKernelArguments(kernelPtr,
1383                                       config,
1384                                       kernelParamsPtr,
1385                                       &kernelParamsPtr->atoms.d_theta,
1386                                       &kernelParamsPtr->atoms.d_dtheta,
1387                                       &kernelParamsPtr->atoms.d_gridlineIndices,
1388                                       &kernelParamsPtr->grid.d_realGrid[FEP_STATE_A],
1389                                       &kernelParamsPtr->grid.d_realGrid[FEP_STATE_B],
1390                                       &kernelParamsPtr->grid.d_fractShiftsTable,
1391                                       &kernelParamsPtr->grid.d_gridlineIndicesTable,
1392                                       &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_A],
1393                                       &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_B],
1394                                       &kernelParamsPtr->atoms.d_coordinates);
1395 #endif
1396
1397     launchGpuKernel(
1398             kernelPtr, config, pmeGpu->archSpecific->pmeStream_, timingEvent, "PME spline/spread", kernelArgs);
1399     pme_gpu_stop_timing(pmeGpu, timingId);
1400
1401     const auto& settings    = pmeGpu->settings;
1402     const bool copyBackGrid = spreadCharges && (!settings.performGPUFFT || settings.copyAllOutputs);
1403     if (copyBackGrid)
1404     {
1405         for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
1406         {
1407             float* h_grid = h_grids[gridIndex];
1408             pme_gpu_copy_output_spread_grid(pmeGpu, h_grid, gridIndex);
1409         }
1410     }
1411     const bool copyBackAtomData =
1412             computeSplines && (!settings.performGPUGather || settings.copyAllOutputs);
1413     if (copyBackAtomData)
1414     {
1415         pme_gpu_copy_output_spread_atom_data(pmeGpu);
1416     }
1417 }
1418
1419 void pme_gpu_solve(const PmeGpu* pmeGpu,
1420                    const int     gridIndex,
1421                    t_complex*    h_grid,
1422                    GridOrdering  gridOrdering,
1423                    bool          computeEnergyAndVirial)
1424 {
1425     GMX_ASSERT(
1426             pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
1427             "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
1428     GMX_ASSERT(gridIndex < pmeGpu->common->ngrids,
1429                "Invalid combination of gridIndex and number of grids");
1430
1431     const auto& settings               = pmeGpu->settings;
1432     const bool  copyInputAndOutputGrid = !settings.performGPUFFT || settings.copyAllOutputs;
1433
1434     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1435
1436     float* h_gridFloat = reinterpret_cast<float*>(h_grid);
1437     if (copyInputAndOutputGrid)
1438     {
1439         copyToDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid[gridIndex],
1440                            h_gridFloat,
1441                            0,
1442                            pmeGpu->archSpecific->complexGridSize[gridIndex],
1443                            pmeGpu->archSpecific->pmeStream_,
1444                            pmeGpu->settings.transferKind,
1445                            nullptr);
1446     }
1447
1448     int majorDim = -1, middleDim = -1, minorDim = -1;
1449     switch (gridOrdering)
1450     {
1451         case GridOrdering::YZX:
1452             majorDim  = YY;
1453             middleDim = ZZ;
1454             minorDim  = XX;
1455             break;
1456
1457         case GridOrdering::XYZ:
1458             majorDim  = XX;
1459             middleDim = YY;
1460             minorDim  = ZZ;
1461             break;
1462
1463         default: GMX_ASSERT(false, "Implement grid ordering here and below for the kernel launch");
1464     }
1465
1466     const int maxBlockSize = pmeGpu->programHandle_->impl_->solveMaxWorkGroupSize;
1467
1468     const int gridLineSize      = pmeGpu->kernelParams->grid.complexGridSize[minorDim];
1469     const int gridLinesPerBlock = std::max(maxBlockSize / gridLineSize, 1);
1470     const int blocksPerGridLine = (gridLineSize + maxBlockSize - 1) / maxBlockSize;
1471     int       cellsPerBlock;
1472     if (blocksPerGridLine == 1)
1473     {
1474         cellsPerBlock = gridLineSize * gridLinesPerBlock;
1475     }
1476     else
1477     {
1478         cellsPerBlock = (gridLineSize + blocksPerGridLine - 1) / blocksPerGridLine;
1479     }
1480     const int warpSize  = pmeGpu->programHandle_->warpSize();
1481     const int blockSize = (cellsPerBlock + warpSize - 1) / warpSize * warpSize;
1482
1483     static_assert(!GMX_GPU_CUDA || c_solveMaxWarpsPerBlock / 2 >= 4,
1484                   "The CUDA solve energy kernels needs at least 4 warps. "
1485                   "Here we launch at least half of the max warps.");
1486
1487     KernelLaunchConfig config;
1488     config.blockSize[0] = blockSize;
1489     config.gridSize[0]  = blocksPerGridLine;
1490     // rounding up to full warps so that shuffle operations produce defined results
1491     config.gridSize[1] = (pmeGpu->kernelParams->grid.complexGridSize[middleDim] + gridLinesPerBlock - 1)
1492                          / gridLinesPerBlock;
1493     config.gridSize[2] = pmeGpu->kernelParams->grid.complexGridSize[majorDim];
1494
1495     PmeStage                           timingId  = PmeStage::Solve;
1496     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1497     if (gridOrdering == GridOrdering::YZX)
1498     {
1499         if (gridIndex == 0)
1500         {
1501             kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveYZXEnergyKernelA
1502                                                : pmeGpu->programHandle_->impl_->solveYZXKernelA;
1503         }
1504         else
1505         {
1506             kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveYZXEnergyKernelB
1507                                                : pmeGpu->programHandle_->impl_->solveYZXKernelB;
1508         }
1509     }
1510     else if (gridOrdering == GridOrdering::XYZ)
1511     {
1512         if (gridIndex == 0)
1513         {
1514             kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveXYZEnergyKernelA
1515                                                : pmeGpu->programHandle_->impl_->solveXYZKernelA;
1516         }
1517         else
1518         {
1519             kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveXYZEnergyKernelB
1520                                                : pmeGpu->programHandle_->impl_->solveXYZKernelB;
1521         }
1522     }
1523
1524     pme_gpu_start_timing(pmeGpu, timingId);
1525     auto* timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1526 #if c_canEmbedBuffers
1527     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1528 #else
1529     const auto kernelArgs =
1530             prepareGpuKernelArguments(kernelPtr,
1531                                       config,
1532                                       kernelParamsPtr,
1533                                       &kernelParamsPtr->grid.d_splineModuli[gridIndex],
1534                                       &kernelParamsPtr->constants.d_virialAndEnergy[gridIndex],
1535                                       &kernelParamsPtr->grid.d_fourierGrid[gridIndex]);
1536 #endif
1537     launchGpuKernel(kernelPtr, config, pmeGpu->archSpecific->pmeStream_, timingEvent, "PME solve", kernelArgs);
1538     pme_gpu_stop_timing(pmeGpu, timingId);
1539
1540     if (computeEnergyAndVirial)
1541     {
1542         copyFromDeviceBuffer(pmeGpu->staging.h_virialAndEnergy[gridIndex],
1543                              &kernelParamsPtr->constants.d_virialAndEnergy[gridIndex],
1544                              0,
1545                              c_virialAndEnergyCount,
1546                              pmeGpu->archSpecific->pmeStream_,
1547                              pmeGpu->settings.transferKind,
1548                              nullptr);
1549     }
1550
1551     if (copyInputAndOutputGrid)
1552     {
1553         copyFromDeviceBuffer(h_gridFloat,
1554                              &kernelParamsPtr->grid.d_fourierGrid[gridIndex],
1555                              0,
1556                              pmeGpu->archSpecific->complexGridSize[gridIndex],
1557                              pmeGpu->archSpecific->pmeStream_,
1558                              pmeGpu->settings.transferKind,
1559                              nullptr);
1560     }
1561 }
1562
1563 /*! \brief
1564  * Returns a pointer to appropriate gather kernel based on the inputvalues
1565  *
1566  * \param[in]  pmeGpu                   The PME GPU structure.
1567  * \param[in]  threadsPerAtom           Controls whether we should use order or order*order threads per atom
1568  * \param[in]  readSplinesFromGlobal    bool controlling if we should write spline data to global memory
1569  * \param[in]  numGrids                 Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1570  *
1571  * \return Pointer to CUDA kernel
1572  */
1573 inline auto selectGatherKernelPtr(const PmeGpu*  pmeGpu,
1574                                   ThreadsPerAtom threadsPerAtom,
1575                                   bool           readSplinesFromGlobal,
1576                                   const int      numGrids)
1577
1578 {
1579     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1580
1581     if (readSplinesFromGlobal)
1582     {
1583         if (threadsPerAtom == ThreadsPerAtom::Order)
1584         {
1585             if (numGrids == 2)
1586             {
1587                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesThPerAtom4Dual;
1588             }
1589             else
1590             {
1591                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesThPerAtom4Single;
1592             }
1593         }
1594         else
1595         {
1596             if (numGrids == 2)
1597             {
1598                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesDual;
1599             }
1600             else
1601             {
1602                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesSingle;
1603             }
1604         }
1605     }
1606     else
1607     {
1608         if (threadsPerAtom == ThreadsPerAtom::Order)
1609         {
1610             if (numGrids == 2)
1611             {
1612                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelThPerAtom4Dual;
1613             }
1614             else
1615             {
1616                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelThPerAtom4Single;
1617             }
1618         }
1619         else
1620         {
1621             if (numGrids == 2)
1622             {
1623                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelDual;
1624             }
1625             else
1626             {
1627                 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelSingle;
1628             }
1629         }
1630     }
1631     return kernelPtr;
1632 }
1633
1634 void pme_gpu_gather(PmeGpu* pmeGpu, real** h_grids, const float lambda)
1635 {
1636     GMX_ASSERT(
1637             pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
1638             "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
1639
1640     const auto& settings = pmeGpu->settings;
1641
1642     if (!settings.performGPUFFT || settings.copyAllOutputs)
1643     {
1644         for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
1645         {
1646             float* h_grid = const_cast<float*>(h_grids[gridIndex]);
1647             pme_gpu_copy_input_gather_grid(pmeGpu, h_grid, gridIndex);
1648         }
1649     }
1650
1651     if (settings.copyAllOutputs)
1652     {
1653         pme_gpu_copy_input_gather_atom_data(pmeGpu);
1654     }
1655
1656     /* Set if we have unit tests */
1657     const bool   readGlobal = pmeGpu->settings.copyAllOutputs;
1658     const size_t blockSize  = pmeGpu->programHandle_->impl_->gatherWorkGroupSize;
1659     const int    order      = pmeGpu->common->pme_order;
1660     GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1661     const int threadsPerAtom =
1662             (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
1663     const bool recalculateSplines = pmeGpu->settings.recalculateSplines;
1664
1665     GMX_ASSERT(!GMX_GPU_OPENCL || pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::OrderSquared,
1666                "Only 16 threads per atom supported in OpenCL");
1667     GMX_ASSERT(!GMX_GPU_OPENCL || !recalculateSplines,
1668                "Recalculating splines not supported in OpenCL");
1669
1670     const int atomsPerBlock = blockSize / threadsPerAtom;
1671
1672     GMX_ASSERT(!(c_pmeAtomDataBlockSize % atomsPerBlock),
1673                "inconsistent atom data padding vs. gathering block size");
1674
1675     const int blockCount = pmeGpu->nAtomsAlloc / atomsPerBlock;
1676     auto      dimGrid    = pmeGpuCreateGrid(pmeGpu, blockCount);
1677
1678     KernelLaunchConfig config;
1679     config.blockSize[0] = order;
1680     config.blockSize[1] = (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? 1 : order);
1681     config.blockSize[2] = atomsPerBlock;
1682     config.gridSize[0]  = dimGrid.first;
1683     config.gridSize[1]  = dimGrid.second;
1684
1685     // TODO test different cache configs
1686
1687     PmeStage                           timingId = PmeStage::Gather;
1688     PmeGpuProgramImpl::PmeKernelHandle kernelPtr =
1689             selectGatherKernelPtr(pmeGpu,
1690                                   pmeGpu->settings.threadsPerAtom,
1691                                   readGlobal || (!recalculateSplines),
1692                                   pmeGpu->common->ngrids);
1693     // TODO design kernel selection getters and make PmeGpu a friend of PmeGpuProgramImpl
1694
1695     pme_gpu_start_timing(pmeGpu, timingId);
1696     auto* timingEvent     = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1697     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1698     if (pmeGpu->common->ngrids == 1)
1699     {
1700         kernelParamsPtr->current.scale = 1.0;
1701     }
1702     else
1703     {
1704         kernelParamsPtr->current.scale = 1.0 - lambda;
1705     }
1706
1707 #if c_canEmbedBuffers
1708     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1709 #else
1710     const auto kernelArgs =
1711             prepareGpuKernelArguments(kernelPtr,
1712                                       config,
1713                                       kernelParamsPtr,
1714                                       &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_A],
1715                                       &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_B],
1716                                       &kernelParamsPtr->grid.d_realGrid[FEP_STATE_A],
1717                                       &kernelParamsPtr->grid.d_realGrid[FEP_STATE_B],
1718                                       &kernelParamsPtr->atoms.d_theta,
1719                                       &kernelParamsPtr->atoms.d_dtheta,
1720                                       &kernelParamsPtr->atoms.d_gridlineIndices,
1721                                       &kernelParamsPtr->atoms.d_forces);
1722 #endif
1723     launchGpuKernel(kernelPtr, config, pmeGpu->archSpecific->pmeStream_, timingEvent, "PME gather", kernelArgs);
1724     pme_gpu_stop_timing(pmeGpu, timingId);
1725
1726     if (pmeGpu->settings.useGpuForceReduction)
1727     {
1728         pmeGpu->archSpecific->pmeForcesReady.markEvent(pmeGpu->archSpecific->pmeStream_);
1729     }
1730     else
1731     {
1732         pme_gpu_copy_output_forces(pmeGpu);
1733     }
1734 }
1735
1736 DeviceBuffer<gmx::RVec> pme_gpu_get_kernelparam_forces(const PmeGpu* pmeGpu)
1737 {
1738     if (pmeGpu && pmeGpu->kernelParams)
1739     {
1740         return pmeGpu->kernelParams->atoms.d_forces;
1741     }
1742     else
1743     {
1744         return DeviceBuffer<gmx::RVec>{};
1745     }
1746 }
1747
1748 void pme_gpu_set_kernelparam_coordinates(const PmeGpu* pmeGpu, DeviceBuffer<gmx::RVec> d_x)
1749 {
1750     GMX_ASSERT(pmeGpu && pmeGpu->kernelParams,
1751                "PME GPU device buffer can not be set in non-GPU builds or before the GPU PME was "
1752                "initialized.");
1753
1754     GMX_ASSERT(checkDeviceBuffer(d_x, pmeGpu->kernelParams->atoms.nAtoms),
1755                "The device-side buffer can not be set.");
1756
1757     pmeGpu->kernelParams->atoms.d_coordinates = d_x;
1758 }
1759
1760 GpuEventSynchronizer* pme_gpu_get_forces_ready_synchronizer(const PmeGpu* pmeGpu)
1761 {
1762     if (pmeGpu && pmeGpu->kernelParams)
1763     {
1764         return &pmeGpu->archSpecific->pmeForcesReady;
1765     }
1766     else
1767     {
1768         return nullptr;
1769     }
1770 }