a7d6e1e96309ff60e7ab1ce1c4381fa78660214f
[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, 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  *
38  * \brief This file contains internal function implementations
39  * for performing the PME calculations on GPU.
40  *
41  * Note that this file is compiled as regular C++ source in OpenCL builds, but
42  * it is treated as CUDA source in CUDA-enabled GPU builds.
43  *
44  * \author Aleksei Iupinov <a.yupinov@gmail.com>
45  * \ingroup module_ewald
46  */
47
48 #include "gmxpre.h"
49
50 #include "pme_gpu_internal.h"
51
52 #include "config.h"
53
54 #include <list>
55 #include <memory>
56 #include <string>
57
58 #include "gromacs/ewald/ewald_utils.h"
59 #include "gromacs/gpu_utils/gpu_utils.h"
60 #include "gromacs/math/invertmatrix.h"
61 #include "gromacs/math/units.h"
62 #include "gromacs/timing/gpu_timing.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/fatalerror.h"
65 #include "gromacs/utility/gmxassert.h"
66 #include "gromacs/utility/logger.h"
67 #include "gromacs/utility/stringutil.h"
68
69 #if GMX_GPU == GMX_GPU_CUDA
70 #    include "gromacs/gpu_utils/pmalloc_cuda.h"
71
72 #    include "pme.cuh"
73 #elif GMX_GPU == GMX_GPU_OPENCL
74 #    include "gromacs/gpu_utils/gmxopencl.h"
75 #endif
76
77 #include "gromacs/ewald/pme.h"
78
79 #include "pme_gpu_3dfft.h"
80 #include "pme_gpu_constants.h"
81 #include "pme_gpu_program_impl.h"
82 #include "pme_gpu_timings.h"
83 #include "pme_gpu_types.h"
84 #include "pme_gpu_types_host.h"
85 #include "pme_gpu_types_host_impl.h"
86 #include "pme_gpu_utils.h"
87 #include "pme_grid.h"
88 #include "pme_internal.h"
89 #include "pme_solve.h"
90
91 /*! \brief
92  * CUDA only
93  * Atom limit above which it is advantageous to turn on the
94  * recalcuating of the splines in the gather and using less threads per atom in the spline and spread
95  */
96 constexpr int c_pmeGpuPerformanceAtomLimit = 23000;
97
98 /*! \internal \brief
99  * Wrapper for getting a pointer to the plain C++ part of the GPU kernel parameters structure.
100  *
101  * \param[in] pmeGpu  The PME GPU structure.
102  * \returns The pointer to the kernel parameters.
103  */
104 static PmeGpuKernelParamsBase* pme_gpu_get_kernel_params_base_ptr(const PmeGpu* pmeGpu)
105 {
106     // reinterpret_cast is needed because the derived CUDA structure is not known in this file
107     auto* kernelParamsPtr = reinterpret_cast<PmeGpuKernelParamsBase*>(pmeGpu->kernelParams.get());
108     return kernelParamsPtr;
109 }
110
111 int pme_gpu_get_atom_data_alignment(const PmeGpu* /*unused*/)
112 {
113     // TODO: this can be simplified, as c_pmeAtomDataAlignment is now constant
114     if (c_usePadding)
115     {
116         return c_pmeAtomDataAlignment;
117     }
118     else
119     {
120         return 0;
121     }
122 }
123
124 int pme_gpu_get_atoms_per_warp(const PmeGpu* pmeGpu)
125 {
126     if (pmeGpu->settings.useOrderThreadsPerAtom)
127     {
128         return pmeGpu->programHandle_->impl_->warpSize / c_pmeSpreadGatherThreadsPerAtom4ThPerAtom;
129     }
130     else
131     {
132         return pmeGpu->programHandle_->impl_->warpSize / c_pmeSpreadGatherThreadsPerAtom;
133     }
134 }
135
136 void pme_gpu_synchronize(const PmeGpu* pmeGpu)
137 {
138     gpuStreamSynchronize(pmeGpu->archSpecific->pmeStream);
139 }
140
141 void pme_gpu_alloc_energy_virial(PmeGpu* pmeGpu)
142 {
143     const size_t energyAndVirialSize = c_virialAndEnergyCount * sizeof(float);
144     allocateDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy, c_virialAndEnergyCount,
145                          pmeGpu->archSpecific->context);
146     pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_virialAndEnergy), energyAndVirialSize);
147 }
148
149 void pme_gpu_free_energy_virial(PmeGpu* pmeGpu)
150 {
151     freeDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy);
152     pfree(pmeGpu->staging.h_virialAndEnergy);
153     pmeGpu->staging.h_virialAndEnergy = nullptr;
154 }
155
156 void pme_gpu_clear_energy_virial(const PmeGpu* pmeGpu)
157 {
158     clearDeviceBufferAsync(&pmeGpu->kernelParams->constants.d_virialAndEnergy, 0,
159                            c_virialAndEnergyCount, pmeGpu->archSpecific->pmeStream);
160 }
161
162 void pme_gpu_realloc_and_copy_bspline_values(PmeGpu* pmeGpu)
163 {
164     const int splineValuesOffset[DIM] = { 0, pmeGpu->kernelParams->grid.realGridSize[XX],
165                                           pmeGpu->kernelParams->grid.realGridSize[XX]
166                                                   + pmeGpu->kernelParams->grid.realGridSize[YY] };
167     memcpy(&pmeGpu->kernelParams->grid.splineValuesOffset, &splineValuesOffset, sizeof(splineValuesOffset));
168
169     const int newSplineValuesSize = pmeGpu->kernelParams->grid.realGridSize[XX]
170                                     + pmeGpu->kernelParams->grid.realGridSize[YY]
171                                     + pmeGpu->kernelParams->grid.realGridSize[ZZ];
172     const bool shouldRealloc = (newSplineValuesSize > pmeGpu->archSpecific->splineValuesSize);
173     reallocateDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli, newSplineValuesSize,
174                            &pmeGpu->archSpecific->splineValuesSize,
175                            &pmeGpu->archSpecific->splineValuesSizeAlloc, pmeGpu->archSpecific->context);
176     if (shouldRealloc)
177     {
178         /* Reallocate the host buffer */
179         pfree(pmeGpu->staging.h_splineModuli);
180         pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_splineModuli),
181                 newSplineValuesSize * sizeof(float));
182     }
183     for (int i = 0; i < DIM; i++)
184     {
185         memcpy(pmeGpu->staging.h_splineModuli + splineValuesOffset[i],
186                pmeGpu->common->bsp_mod[i].data(), pmeGpu->common->bsp_mod[i].size() * sizeof(float));
187     }
188     /* TODO: pin original buffer instead! */
189     copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli, pmeGpu->staging.h_splineModuli,
190                        0, newSplineValuesSize, pmeGpu->archSpecific->pmeStream,
191                        pmeGpu->settings.transferKind, nullptr);
192 }
193
194 void pme_gpu_free_bspline_values(const PmeGpu* pmeGpu)
195 {
196     pfree(pmeGpu->staging.h_splineModuli);
197     freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli);
198 }
199
200 void pme_gpu_realloc_forces(PmeGpu* pmeGpu)
201 {
202     const size_t newForcesSize = pmeGpu->nAtomsAlloc * DIM;
203     GMX_ASSERT(newForcesSize > 0, "Bad number of atoms in PME GPU");
204     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces, newForcesSize,
205                            &pmeGpu->archSpecific->forcesSize,
206                            &pmeGpu->archSpecific->forcesSizeAlloc, pmeGpu->archSpecific->context);
207     pmeGpu->staging.h_forces.reserveWithPadding(pmeGpu->nAtomsAlloc);
208     pmeGpu->staging.h_forces.resizeWithPadding(pmeGpu->kernelParams->atoms.nAtoms);
209 }
210
211 void pme_gpu_free_forces(const PmeGpu* pmeGpu)
212 {
213     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces);
214 }
215
216 void pme_gpu_copy_input_forces(PmeGpu* pmeGpu)
217 {
218     GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
219     float* h_forcesFloat = reinterpret_cast<float*>(pmeGpu->staging.h_forces.data());
220     copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces, h_forcesFloat, 0,
221                        DIM * pmeGpu->kernelParams->atoms.nAtoms, pmeGpu->archSpecific->pmeStream,
222                        pmeGpu->settings.transferKind, nullptr);
223 }
224
225 void pme_gpu_copy_output_forces(PmeGpu* pmeGpu)
226 {
227     GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
228     float* h_forcesFloat = reinterpret_cast<float*>(pmeGpu->staging.h_forces.data());
229     copyFromDeviceBuffer(h_forcesFloat, &pmeGpu->kernelParams->atoms.d_forces, 0,
230                          DIM * pmeGpu->kernelParams->atoms.nAtoms, pmeGpu->archSpecific->pmeStream,
231                          pmeGpu->settings.transferKind, nullptr);
232 }
233
234 void pme_gpu_realloc_coordinates(const PmeGpu* pmeGpu)
235 {
236     const size_t newCoordinatesSize = pmeGpu->nAtomsAlloc * DIM;
237     GMX_ASSERT(newCoordinatesSize > 0, "Bad number of atoms in PME GPU");
238     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates, newCoordinatesSize,
239                            &pmeGpu->archSpecific->coordinatesSize,
240                            &pmeGpu->archSpecific->coordinatesSizeAlloc, pmeGpu->archSpecific->context);
241     if (c_usePadding)
242     {
243         const size_t paddingIndex = DIM * pmeGpu->kernelParams->atoms.nAtoms;
244         const size_t paddingCount = DIM * pmeGpu->nAtomsAlloc - paddingIndex;
245         if (paddingCount > 0)
246         {
247             clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coordinates, paddingIndex,
248                                    paddingCount, pmeGpu->archSpecific->pmeStream);
249         }
250     }
251 }
252
253 void pme_gpu_free_coordinates(const PmeGpu* pmeGpu)
254 {
255     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates);
256 }
257
258 void pme_gpu_realloc_and_copy_input_coefficients(const PmeGpu* pmeGpu, const float* h_coefficients)
259 {
260     GMX_ASSERT(h_coefficients, "Bad host-side charge buffer in PME GPU");
261     const size_t newCoefficientsSize = pmeGpu->nAtomsAlloc;
262     GMX_ASSERT(newCoefficientsSize > 0, "Bad number of atoms in PME GPU");
263     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients, newCoefficientsSize,
264                            &pmeGpu->archSpecific->coefficientsSize,
265                            &pmeGpu->archSpecific->coefficientsSizeAlloc, pmeGpu->archSpecific->context);
266     copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients,
267                        const_cast<float*>(h_coefficients), 0, pmeGpu->kernelParams->atoms.nAtoms,
268                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
269     if (c_usePadding)
270     {
271         const size_t paddingIndex = pmeGpu->kernelParams->atoms.nAtoms;
272         const size_t paddingCount = pmeGpu->nAtomsAlloc - paddingIndex;
273         if (paddingCount > 0)
274         {
275             clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coefficients, paddingIndex,
276                                    paddingCount, pmeGpu->archSpecific->pmeStream);
277         }
278     }
279 }
280
281 void pme_gpu_free_coefficients(const PmeGpu* pmeGpu)
282 {
283     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients);
284 }
285
286 void pme_gpu_realloc_spline_data(PmeGpu* pmeGpu)
287 {
288     const int    order        = pmeGpu->common->pme_order;
289     const int    alignment    = pme_gpu_get_atoms_per_warp(pmeGpu);
290     const size_t nAtomsPadded = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
291     const int    newSplineDataSize = DIM * order * nAtomsPadded;
292     GMX_ASSERT(newSplineDataSize > 0, "Bad number of atoms in PME GPU");
293     /* Two arrays of the same size */
294     const bool shouldRealloc        = (newSplineDataSize > pmeGpu->archSpecific->splineDataSize);
295     int        currentSizeTemp      = pmeGpu->archSpecific->splineDataSize;
296     int        currentSizeTempAlloc = pmeGpu->archSpecific->splineDataSizeAlloc;
297     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta, newSplineDataSize,
298                            &currentSizeTemp, &currentSizeTempAlloc, pmeGpu->archSpecific->context);
299     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta, newSplineDataSize,
300                            &pmeGpu->archSpecific->splineDataSize,
301                            &pmeGpu->archSpecific->splineDataSizeAlloc, pmeGpu->archSpecific->context);
302     // the host side reallocation
303     if (shouldRealloc)
304     {
305         pfree(pmeGpu->staging.h_theta);
306         pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_theta), newSplineDataSize * sizeof(float));
307         pfree(pmeGpu->staging.h_dtheta);
308         pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_dtheta), newSplineDataSize * sizeof(float));
309     }
310 }
311
312 void pme_gpu_free_spline_data(const PmeGpu* pmeGpu)
313 {
314     /* Two arrays of the same size */
315     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta);
316     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta);
317     pfree(pmeGpu->staging.h_theta);
318     pfree(pmeGpu->staging.h_dtheta);
319 }
320
321 void pme_gpu_realloc_grid_indices(PmeGpu* pmeGpu)
322 {
323     const size_t newIndicesSize = DIM * pmeGpu->nAtomsAlloc;
324     GMX_ASSERT(newIndicesSize > 0, "Bad number of atoms in PME GPU");
325     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices, newIndicesSize,
326                            &pmeGpu->archSpecific->gridlineIndicesSize,
327                            &pmeGpu->archSpecific->gridlineIndicesSizeAlloc, pmeGpu->archSpecific->context);
328     pfree(pmeGpu->staging.h_gridlineIndices);
329     pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_gridlineIndices), newIndicesSize * sizeof(int));
330 }
331
332 void pme_gpu_free_grid_indices(const PmeGpu* pmeGpu)
333 {
334     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices);
335     pfree(pmeGpu->staging.h_gridlineIndices);
336 }
337
338 void pme_gpu_realloc_grids(PmeGpu* pmeGpu)
339 {
340     auto*     kernelParamsPtr = pmeGpu->kernelParams.get();
341     const int newRealGridSize = kernelParamsPtr->grid.realGridSizePadded[XX]
342                                 * kernelParamsPtr->grid.realGridSizePadded[YY]
343                                 * kernelParamsPtr->grid.realGridSizePadded[ZZ];
344     const int newComplexGridSize = kernelParamsPtr->grid.complexGridSizePadded[XX]
345                                    * kernelParamsPtr->grid.complexGridSizePadded[YY]
346                                    * kernelParamsPtr->grid.complexGridSizePadded[ZZ] * 2;
347     // Multiplied by 2 because we count complex grid size for complex numbers, but all allocations/pointers are float
348     if (pmeGpu->archSpecific->performOutOfPlaceFFT)
349     {
350         /* 2 separate grids */
351         reallocateDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid, newComplexGridSize,
352                                &pmeGpu->archSpecific->complexGridSize,
353                                &pmeGpu->archSpecific->complexGridSizeAlloc, pmeGpu->archSpecific->context);
354         reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid, newRealGridSize,
355                                &pmeGpu->archSpecific->realGridSize,
356                                &pmeGpu->archSpecific->realGridSizeAlloc, pmeGpu->archSpecific->context);
357     }
358     else
359     {
360         /* A single buffer so that any grid will fit */
361         const int newGridsSize = std::max(newRealGridSize, newComplexGridSize);
362         reallocateDeviceBuffer(
363                 &kernelParamsPtr->grid.d_realGrid, newGridsSize, &pmeGpu->archSpecific->realGridSize,
364                 &pmeGpu->archSpecific->realGridSizeAlloc, pmeGpu->archSpecific->context);
365         kernelParamsPtr->grid.d_fourierGrid   = kernelParamsPtr->grid.d_realGrid;
366         pmeGpu->archSpecific->complexGridSize = pmeGpu->archSpecific->realGridSize;
367         // the size might get used later for copying the grid
368     }
369 }
370
371 void pme_gpu_free_grids(const PmeGpu* pmeGpu)
372 {
373     if (pmeGpu->archSpecific->performOutOfPlaceFFT)
374     {
375         freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_fourierGrid);
376     }
377     freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid);
378 }
379
380 void pme_gpu_clear_grids(const PmeGpu* pmeGpu)
381 {
382     clearDeviceBufferAsync(&pmeGpu->kernelParams->grid.d_realGrid, 0,
383                            pmeGpu->archSpecific->realGridSize, pmeGpu->archSpecific->pmeStream);
384 }
385
386 void pme_gpu_realloc_and_copy_fract_shifts(PmeGpu* pmeGpu)
387 {
388     pme_gpu_free_fract_shifts(pmeGpu);
389
390     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
391
392     const int nx                  = kernelParamsPtr->grid.realGridSize[XX];
393     const int ny                  = kernelParamsPtr->grid.realGridSize[YY];
394     const int nz                  = kernelParamsPtr->grid.realGridSize[ZZ];
395     const int cellCount           = c_pmeNeighborUnitcellCount;
396     const int gridDataOffset[DIM] = { 0, cellCount * nx, cellCount * (nx + ny) };
397
398     memcpy(kernelParamsPtr->grid.tablesOffsets, &gridDataOffset, sizeof(gridDataOffset));
399
400     const int newFractShiftsSize = cellCount * (nx + ny + nz);
401
402 #if GMX_GPU == GMX_GPU_CUDA
403     initParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable, kernelParamsPtr->fractShiftsTableTexture,
404                          pmeGpu->common->fsh.data(), newFractShiftsSize);
405
406     initParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
407                          kernelParamsPtr->gridlineIndicesTableTexture, pmeGpu->common->nn.data(),
408                          newFractShiftsSize);
409 #elif GMX_GPU == GMX_GPU_OPENCL
410     // No dedicated texture routines....
411     allocateDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable, newFractShiftsSize,
412                          pmeGpu->archSpecific->context);
413     allocateDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable, newFractShiftsSize,
414                          pmeGpu->archSpecific->context);
415     copyToDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable, pmeGpu->common->fsh.data(), 0,
416                        newFractShiftsSize, pmeGpu->archSpecific->pmeStream,
417                        GpuApiCallBehavior::Async, nullptr);
418     copyToDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable, pmeGpu->common->nn.data(), 0,
419                        newFractShiftsSize, pmeGpu->archSpecific->pmeStream,
420                        GpuApiCallBehavior::Async, nullptr);
421 #endif
422 }
423
424 void pme_gpu_free_fract_shifts(const PmeGpu* pmeGpu)
425 {
426     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
427 #if GMX_GPU == GMX_GPU_CUDA
428     destroyParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable,
429                             kernelParamsPtr->fractShiftsTableTexture);
430     destroyParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
431                             kernelParamsPtr->gridlineIndicesTableTexture);
432 #elif GMX_GPU == GMX_GPU_OPENCL
433     freeDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable);
434     freeDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable);
435 #endif
436 }
437
438 bool pme_gpu_stream_query(const PmeGpu* pmeGpu)
439 {
440     return haveStreamTasksCompleted(pmeGpu->archSpecific->pmeStream);
441 }
442
443 void pme_gpu_copy_input_gather_grid(const PmeGpu* pmeGpu, float* h_grid)
444 {
445     copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid, h_grid, 0, pmeGpu->archSpecific->realGridSize,
446                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
447 }
448
449 void pme_gpu_copy_output_spread_grid(const PmeGpu* pmeGpu, float* h_grid)
450 {
451     copyFromDeviceBuffer(h_grid, &pmeGpu->kernelParams->grid.d_realGrid, 0,
452                          pmeGpu->archSpecific->realGridSize, pmeGpu->archSpecific->pmeStream,
453                          pmeGpu->settings.transferKind, nullptr);
454     pmeGpu->archSpecific->syncSpreadGridD2H.markEvent(pmeGpu->archSpecific->pmeStream);
455 }
456
457 void pme_gpu_copy_output_spread_atom_data(const PmeGpu* pmeGpu)
458 {
459     const int    alignment       = pme_gpu_get_atoms_per_warp(pmeGpu);
460     const size_t nAtomsPadded    = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
461     const size_t splinesCount    = DIM * nAtomsPadded * pmeGpu->common->pme_order;
462     auto*        kernelParamsPtr = pmeGpu->kernelParams.get();
463     copyFromDeviceBuffer(pmeGpu->staging.h_dtheta, &kernelParamsPtr->atoms.d_dtheta, 0, splinesCount,
464                          pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
465     copyFromDeviceBuffer(pmeGpu->staging.h_theta, &kernelParamsPtr->atoms.d_theta, 0, splinesCount,
466                          pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
467     copyFromDeviceBuffer(pmeGpu->staging.h_gridlineIndices, &kernelParamsPtr->atoms.d_gridlineIndices,
468                          0, kernelParamsPtr->atoms.nAtoms * DIM, pmeGpu->archSpecific->pmeStream,
469                          pmeGpu->settings.transferKind, nullptr);
470 }
471
472 void pme_gpu_copy_input_gather_atom_data(const PmeGpu* pmeGpu)
473 {
474     const int    alignment       = pme_gpu_get_atoms_per_warp(pmeGpu);
475     const size_t nAtomsPadded    = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
476     const size_t splinesCount    = DIM * nAtomsPadded * pmeGpu->common->pme_order;
477     auto*        kernelParamsPtr = pmeGpu->kernelParams.get();
478     if (c_usePadding)
479     {
480         // TODO: could clear only the padding and not the whole thing, but this is a test-exclusive code anyway
481         clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_gridlineIndices, 0,
482                                pmeGpu->nAtomsAlloc * DIM, pmeGpu->archSpecific->pmeStream);
483         clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_dtheta, 0,
484                                pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM,
485                                pmeGpu->archSpecific->pmeStream);
486         clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_theta, 0,
487                                pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM,
488                                pmeGpu->archSpecific->pmeStream);
489     }
490     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_dtheta, pmeGpu->staging.h_dtheta, 0, splinesCount,
491                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
492     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_theta, pmeGpu->staging.h_theta, 0, splinesCount,
493                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
494     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_gridlineIndices, pmeGpu->staging.h_gridlineIndices,
495                        0, kernelParamsPtr->atoms.nAtoms * DIM, pmeGpu->archSpecific->pmeStream,
496                        pmeGpu->settings.transferKind, nullptr);
497 }
498
499 void pme_gpu_sync_spread_grid(const PmeGpu* pmeGpu)
500 {
501     pmeGpu->archSpecific->syncSpreadGridD2H.waitForEvent();
502 }
503
504 void pme_gpu_init_internal(PmeGpu* pmeGpu)
505 {
506 #if GMX_GPU == GMX_GPU_CUDA
507     // Prepare to use the device that this PME task was assigned earlier.
508     // Other entities, such as CUDA timing events, are known to implicitly use the device context.
509     CU_RET_ERR(cudaSetDevice(pmeGpu->deviceInfo->id), "Switching to PME CUDA device");
510 #endif
511
512     /* Allocate the target-specific structures */
513     pmeGpu->archSpecific.reset(new PmeGpuSpecific());
514     pmeGpu->kernelParams.reset(new PmeGpuKernelParams());
515
516     pmeGpu->archSpecific->performOutOfPlaceFFT = true;
517     /* This should give better performance, according to the cuFFT documentation.
518      * The performance seems to be the same though.
519      * TODO: PME could also try to pick up nice grid sizes (with factors of 2, 3, 5, 7).
520      */
521
522     // TODO: this is just a convenient reuse because programHandle_ currently is in charge of creating context
523     pmeGpu->archSpecific->context = pmeGpu->programHandle_->impl_->context;
524
525     // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?) and reuse in NB
526     if (GMX_GPU == GMX_GPU_CUDA)
527     {
528         /* WARNING: CUDA timings are incorrect with multiple streams.
529          *          This is the main reason why they are disabled by default.
530          */
531         // TODO: Consider turning on by default when we can detect nr of streams.
532         pmeGpu->archSpecific->useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
533     }
534     else if (GMX_GPU == GMX_GPU_OPENCL)
535     {
536         pmeGpu->archSpecific->useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
537     }
538
539 #if GMX_GPU == GMX_GPU_CUDA
540     pmeGpu->maxGridWidthX = pmeGpu->deviceInfo->prop.maxGridSize[0];
541 #elif GMX_GPU == GMX_GPU_OPENCL
542     pmeGpu->maxGridWidthX = INT32_MAX / 2;
543     // TODO: is there no really global work size limit in OpenCL?
544 #endif
545
546     /* Creating a PME GPU stream:
547      * - default high priority with CUDA
548      * - no priorities implemented yet with OpenCL; see #2532
549      */
550 #if GMX_GPU == GMX_GPU_CUDA
551     cudaError_t stat;
552     int         highest_priority, lowest_priority;
553     stat = cudaDeviceGetStreamPriorityRange(&lowest_priority, &highest_priority);
554     CU_RET_ERR(stat, "PME cudaDeviceGetStreamPriorityRange failed");
555     stat = cudaStreamCreateWithPriority(&pmeGpu->archSpecific->pmeStream,
556                                         cudaStreamDefault, // cudaStreamNonBlocking,
557                                         highest_priority);
558     CU_RET_ERR(stat, "cudaStreamCreateWithPriority on the PME stream failed");
559 #elif GMX_GPU == GMX_GPU_OPENCL
560     cl_command_queue_properties queueProperties =
561             pmeGpu->archSpecific->useTiming ? CL_QUEUE_PROFILING_ENABLE : 0;
562     cl_device_id device_id = pmeGpu->deviceInfo->ocl_gpu_id.ocl_device_id;
563     cl_int       clError;
564     pmeGpu->archSpecific->pmeStream =
565             clCreateCommandQueue(pmeGpu->archSpecific->context, device_id, queueProperties, &clError);
566     if (clError != CL_SUCCESS)
567     {
568         GMX_THROW(gmx::InternalError("Failed to create PME command queue"));
569     }
570 #endif
571 }
572
573 void pme_gpu_destroy_specific(const PmeGpu* pmeGpu)
574 {
575 #if GMX_GPU == GMX_GPU_CUDA
576     /* Destroy the CUDA stream */
577     cudaError_t stat = cudaStreamDestroy(pmeGpu->archSpecific->pmeStream);
578     CU_RET_ERR(stat, "PME cudaStreamDestroy error");
579 #elif GMX_GPU == GMX_GPU_OPENCL
580     cl_int clError = clReleaseCommandQueue(pmeGpu->archSpecific->pmeStream);
581     if (clError != CL_SUCCESS)
582     {
583         gmx_warning("Failed to destroy PME command queue");
584     }
585 #endif
586 }
587
588 void pme_gpu_reinit_3dfft(const PmeGpu* pmeGpu)
589 {
590     if (pme_gpu_performs_FFT(pmeGpu))
591     {
592         pmeGpu->archSpecific->fftSetup.resize(0);
593         for (int i = 0; i < pmeGpu->common->ngrids; i++)
594         {
595             pmeGpu->archSpecific->fftSetup.push_back(std::make_unique<GpuParallel3dFft>(pmeGpu));
596         }
597     }
598 }
599
600 void pme_gpu_destroy_3dfft(const PmeGpu* pmeGpu)
601 {
602     pmeGpu->archSpecific->fftSetup.resize(0);
603 }
604
605 int getSplineParamFullIndex(int order, int splineIndex, int dimIndex, int atomIndex, int atomsPerWarp)
606 {
607     if (order != c_pmeGpuOrder)
608     {
609         throw order;
610     }
611     constexpr int fixedOrder = c_pmeGpuOrder;
612     GMX_UNUSED_VALUE(fixedOrder);
613
614     const int atomWarpIndex = atomIndex % atomsPerWarp;
615     const int warpIndex     = atomIndex / atomsPerWarp;
616     int       indexBase, result;
617     switch (atomsPerWarp)
618     {
619         case 1:
620             indexBase = getSplineParamIndexBase<fixedOrder, 1>(warpIndex, atomWarpIndex);
621             result    = getSplineParamIndex<fixedOrder, 1>(indexBase, dimIndex, splineIndex);
622             break;
623
624         case 2:
625             indexBase = getSplineParamIndexBase<fixedOrder, 2>(warpIndex, atomWarpIndex);
626             result    = getSplineParamIndex<fixedOrder, 2>(indexBase, dimIndex, splineIndex);
627             break;
628
629         case 4:
630             indexBase = getSplineParamIndexBase<fixedOrder, 4>(warpIndex, atomWarpIndex);
631             result    = getSplineParamIndex<fixedOrder, 4>(indexBase, dimIndex, splineIndex);
632             break;
633
634         case 8:
635             indexBase = getSplineParamIndexBase<fixedOrder, 8>(warpIndex, atomWarpIndex);
636             result    = getSplineParamIndex<fixedOrder, 8>(indexBase, dimIndex, splineIndex);
637             break;
638
639         default:
640             GMX_THROW(gmx::NotImplementedError(
641                     gmx::formatString("Test function call not unrolled for atomsPerWarp = %d in "
642                                       "getSplineParamFullIndex",
643                                       atomsPerWarp)));
644     }
645     return result;
646 }
647
648 void pme_gpu_getEnergyAndVirial(const gmx_pme_t& pme, PmeOutput* output)
649 {
650     const PmeGpu* pmeGpu = pme.gpu;
651     for (int j = 0; j < c_virialAndEnergyCount; j++)
652     {
653         GMX_ASSERT(std::isfinite(pmeGpu->staging.h_virialAndEnergy[j]),
654                    "PME GPU produces incorrect energy/virial.");
655     }
656
657     unsigned int j                 = 0;
658     output->coulombVirial_[XX][XX] = 0.25F * pmeGpu->staging.h_virialAndEnergy[j++];
659     output->coulombVirial_[YY][YY] = 0.25F * pmeGpu->staging.h_virialAndEnergy[j++];
660     output->coulombVirial_[ZZ][ZZ] = 0.25F * pmeGpu->staging.h_virialAndEnergy[j++];
661     output->coulombVirial_[XX][YY] = output->coulombVirial_[YY][XX] =
662             0.25F * pmeGpu->staging.h_virialAndEnergy[j++];
663     output->coulombVirial_[XX][ZZ] = output->coulombVirial_[ZZ][XX] =
664             0.25F * pmeGpu->staging.h_virialAndEnergy[j++];
665     output->coulombVirial_[YY][ZZ] = output->coulombVirial_[ZZ][YY] =
666             0.25F * pmeGpu->staging.h_virialAndEnergy[j++];
667     output->coulombEnergy_ = 0.5F * pmeGpu->staging.h_virialAndEnergy[j++];
668 }
669
670 /*! \brief Sets the force-related members in \p output
671  *
672  * \param[in]   pmeGpu      PME GPU data structure
673  * \param[out]  output      Pointer to PME output data structure
674  */
675 static void pme_gpu_getForceOutput(PmeGpu* pmeGpu, PmeOutput* output)
676 {
677     output->haveForceOutput_ = !pmeGpu->settings.useGpuForceReduction;
678     if (output->haveForceOutput_)
679     {
680         output->forces_ = pmeGpu->staging.h_forces;
681     }
682 }
683
684 PmeOutput pme_gpu_getOutput(const gmx_pme_t& pme, const int flags)
685 {
686     PmeGpu*    pmeGpu                      = pme.gpu;
687     const bool haveComputedEnergyAndVirial = (flags & GMX_PME_CALC_ENER_VIR) != 0;
688
689     PmeOutput output;
690
691     pme_gpu_getForceOutput(pmeGpu, &output);
692
693     // The caller knows from the flags that the energy and the virial are not usable
694     // on the else branch
695     if (haveComputedEnergyAndVirial)
696     {
697         if (pme_gpu_performs_solve(pmeGpu))
698         {
699             pme_gpu_getEnergyAndVirial(pme, &output);
700         }
701         else
702         {
703             get_pme_ener_vir_q(pme.solve_work, pme.nthread, &output);
704         }
705     }
706     return output;
707 }
708
709 void pme_gpu_update_input_box(PmeGpu gmx_unused* pmeGpu, const matrix gmx_unused box)
710 {
711 #if GMX_DOUBLE
712     GMX_THROW(gmx::NotImplementedError("PME is implemented for single-precision only on GPU"));
713 #else
714     matrix scaledBox;
715     pmeGpu->common->boxScaler->scaleBox(box, scaledBox);
716     auto* kernelParamsPtr              = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
717     kernelParamsPtr->current.boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
718     GMX_ASSERT(kernelParamsPtr->current.boxVolume != 0.0F, "Zero volume of the unit cell");
719     matrix recipBox;
720     gmx::invertBoxMatrix(scaledBox, recipBox);
721
722     /* The GPU recipBox is transposed as compared to the CPU recipBox.
723      * Spread uses matrix columns (while solve and gather use rows).
724      * There is no particular reason for this; it might be further rethought/optimized for better access patterns.
725      */
726     const real newRecipBox[DIM][DIM] = { { recipBox[XX][XX], recipBox[YY][XX], recipBox[ZZ][XX] },
727                                          { 0.0, recipBox[YY][YY], recipBox[ZZ][YY] },
728                                          { 0.0, 0.0, recipBox[ZZ][ZZ] } };
729     memcpy(kernelParamsPtr->current.recipBox, newRecipBox, sizeof(matrix));
730 #endif
731 }
732
733 /*! \brief \libinternal
734  * (Re-)initializes all the PME GPU data related to the grid size and cut-off.
735  *
736  * \param[in] pmeGpu            The PME GPU structure.
737  */
738 static void pme_gpu_reinit_grids(PmeGpu* pmeGpu)
739 {
740     auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
741     kernelParamsPtr->grid.ewaldFactor =
742             (M_PI * M_PI) / (pmeGpu->common->ewaldcoeff_q * pmeGpu->common->ewaldcoeff_q);
743
744     /* The grid size variants */
745     for (int i = 0; i < DIM; i++)
746     {
747         kernelParamsPtr->grid.realGridSize[i] = pmeGpu->common->nk[i];
748         kernelParamsPtr->grid.realGridSizeFP[i] =
749                 static_cast<float>(kernelParamsPtr->grid.realGridSize[i]);
750         kernelParamsPtr->grid.realGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
751
752         // The complex grid currently uses no padding;
753         // if it starts to do so, then another test should be added for that
754         kernelParamsPtr->grid.complexGridSize[i]       = kernelParamsPtr->grid.realGridSize[i];
755         kernelParamsPtr->grid.complexGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
756     }
757     /* FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
758     if (!pme_gpu_performs_FFT(pmeGpu))
759     {
760         // This allows for GPU spreading grid and CPU fftgrid to have the same layout, so that we can copy the data directly
761         kernelParamsPtr->grid.realGridSizePadded[ZZ] =
762                 (kernelParamsPtr->grid.realGridSize[ZZ] / 2 + 1) * 2;
763     }
764
765     /* GPU FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
766     kernelParamsPtr->grid.complexGridSize[ZZ] /= 2;
767     kernelParamsPtr->grid.complexGridSize[ZZ]++;
768     kernelParamsPtr->grid.complexGridSizePadded[ZZ] = kernelParamsPtr->grid.complexGridSize[ZZ];
769
770     pme_gpu_realloc_and_copy_fract_shifts(pmeGpu);
771     pme_gpu_realloc_and_copy_bspline_values(pmeGpu);
772     pme_gpu_realloc_grids(pmeGpu);
773     pme_gpu_reinit_3dfft(pmeGpu);
774 }
775
776 /* Several GPU functions that refer to the CPU PME data live here.
777  * We would like to keep these away from the GPU-framework specific code for clarity,
778  * as well as compilation issues with MPI.
779  */
780
781 /*! \brief \libinternal
782  * Copies everything useful from the PME CPU to the PME GPU structure.
783  * The goal is to minimize interaction with the PME CPU structure in the GPU code.
784  *
785  * \param[in] pme         The PME structure.
786  */
787 static void pme_gpu_copy_common_data_from(const gmx_pme_t* pme)
788 {
789     /* TODO: Consider refactoring the CPU PME code to use the same structure,
790      * so that this function becomes 2 lines */
791     PmeGpu* pmeGpu               = pme->gpu;
792     pmeGpu->common->ngrids       = pme->ngrids;
793     pmeGpu->common->epsilon_r    = pme->epsilon_r;
794     pmeGpu->common->ewaldcoeff_q = pme->ewaldcoeff_q;
795     pmeGpu->common->nk[XX]       = pme->nkx;
796     pmeGpu->common->nk[YY]       = pme->nky;
797     pmeGpu->common->nk[ZZ]       = pme->nkz;
798     pmeGpu->common->pme_order    = pme->pme_order;
799     if (pmeGpu->common->pme_order != c_pmeGpuOrder)
800     {
801         GMX_THROW(gmx::NotImplementedError("pme_order != 4 is not implemented!"));
802     }
803     for (int i = 0; i < DIM; i++)
804     {
805         pmeGpu->common->bsp_mod[i].assign(pme->bsp_mod[i], pme->bsp_mod[i] + pmeGpu->common->nk[i]);
806     }
807     const int cellCount = c_pmeNeighborUnitcellCount;
808     pmeGpu->common->fsh.resize(0);
809     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshx, pme->fshx + cellCount * pme->nkx);
810     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshy, pme->fshy + cellCount * pme->nky);
811     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshz, pme->fshz + cellCount * pme->nkz);
812     pmeGpu->common->nn.resize(0);
813     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnx, pme->nnx + cellCount * pme->nkx);
814     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nny, pme->nny + cellCount * pme->nky);
815     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnz, pme->nnz + cellCount * pme->nkz);
816     pmeGpu->common->runMode       = pme->runMode;
817     pmeGpu->common->isRankPmeOnly = !pme->bPPnode;
818     pmeGpu->common->boxScaler     = pme->boxScaler;
819 }
820
821 /*! \libinternal \brief
822  * uses heuristics to select the best performing PME gather and scatter kernels
823  *
824  * \param[in,out] pmeGpu         The PME GPU structure.
825  */
826 static void pme_gpu_select_best_performing_pme_spreadgather_kernels(PmeGpu* pmeGpu)
827 {
828     if (pmeGpu->kernelParams->atoms.nAtoms > c_pmeGpuPerformanceAtomLimit && (GMX_GPU == GMX_GPU_CUDA))
829     {
830         pmeGpu->settings.useOrderThreadsPerAtom = true;
831         pmeGpu->settings.recalculateSplines     = true;
832     }
833     else
834     {
835         pmeGpu->settings.useOrderThreadsPerAtom = false;
836         pmeGpu->settings.recalculateSplines     = false;
837     }
838 }
839
840
841 /*! \libinternal \brief
842  * Initializes the PME GPU data at the beginning of the run.
843  * TODO: this should become PmeGpu::PmeGpu()
844  *
845  * \param[in,out] pme            The PME structure.
846  * \param[in,out] gpuInfo        The GPU information structure.
847  * \param[in]     pmeGpuProgram  The handle to the program/kernel data created outside (e.g. in unit tests/runner)
848  */
849 static void pme_gpu_init(gmx_pme_t* pme, const gmx_device_info_t* gpuInfo, PmeGpuProgramHandle pmeGpuProgram)
850 {
851     pme->gpu       = new PmeGpu();
852     PmeGpu* pmeGpu = pme->gpu;
853     changePinningPolicy(&pmeGpu->staging.h_forces, pme_get_pinning_policy());
854     pmeGpu->common = std::make_shared<PmeShared>();
855
856     /* These settings are set here for the whole run; dynamic ones are set in pme_gpu_reinit() */
857     /* A convenience variable. */
858     pmeGpu->settings.useDecomposition = (pme->nnodes == 1);
859     /* TODO: CPU gather with GPU spread is broken due to different theta/dtheta layout. */
860     pmeGpu->settings.performGPUGather = true;
861     // By default GPU-side reduction is off (explicitly set here for tests, otherwise reset per-step)
862     pmeGpu->settings.useGpuForceReduction = false;
863
864     pme_gpu_set_testing(pmeGpu, false);
865
866     pmeGpu->deviceInfo = gpuInfo;
867     GMX_ASSERT(pmeGpuProgram != nullptr, "GPU kernels must be already compiled");
868     pmeGpu->programHandle_ = pmeGpuProgram;
869
870     pmeGpu->initializedClfftLibrary_ = std::make_unique<gmx::ClfftInitializer>();
871
872     pme_gpu_init_internal(pmeGpu);
873     pme_gpu_alloc_energy_virial(pmeGpu);
874
875     pme_gpu_copy_common_data_from(pme);
876
877     GMX_ASSERT(pmeGpu->common->epsilon_r != 0.0F, "PME GPU: bad electrostatic coefficient");
878
879     auto* kernelParamsPtr               = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
880     kernelParamsPtr->constants.elFactor = ONE_4PI_EPS0 / pmeGpu->common->epsilon_r;
881 }
882
883 void pme_gpu_transform_spline_atom_data(const PmeGpu*      pmeGpu,
884                                         const PmeAtomComm* atc,
885                                         PmeSplineDataType  type,
886                                         int                dimIndex,
887                                         PmeLayoutTransform transform)
888 {
889     // The GPU atom spline data is laid out in a different way currently than the CPU one.
890     // This function converts the data from GPU to CPU layout (in the host memory).
891     // It is only intended for testing purposes so far.
892     // Ideally we should use similar layouts on CPU and GPU if we care about mixed modes and their
893     // performance (e.g. spreading on GPU, gathering on CPU).
894     GMX_RELEASE_ASSERT(atc->nthread == 1, "Only the serial PME data layout is supported");
895     const uintmax_t threadIndex  = 0;
896     const auto      atomCount    = pme_gpu_get_kernel_params_base_ptr(pmeGpu)->atoms.nAtoms;
897     const auto      atomsPerWarp = pme_gpu_get_atoms_per_warp(pmeGpu);
898     const auto      pmeOrder     = pmeGpu->common->pme_order;
899     GMX_ASSERT(pmeOrder == c_pmeGpuOrder, "Only PME order 4 is implemented");
900
901     real*  cpuSplineBuffer;
902     float* h_splineBuffer;
903     switch (type)
904     {
905         case PmeSplineDataType::Values:
906             cpuSplineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
907             h_splineBuffer  = pmeGpu->staging.h_theta;
908             break;
909
910         case PmeSplineDataType::Derivatives:
911             cpuSplineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
912             h_splineBuffer  = pmeGpu->staging.h_dtheta;
913             break;
914
915         default: GMX_THROW(gmx::InternalError("Unknown spline data type"));
916     }
917
918     for (auto atomIndex = 0; atomIndex < atomCount; atomIndex++)
919     {
920         for (auto orderIndex = 0; orderIndex < pmeOrder; orderIndex++)
921         {
922             const auto gpuValueIndex =
923                     getSplineParamFullIndex(pmeOrder, orderIndex, dimIndex, atomIndex, atomsPerWarp);
924             const auto cpuValueIndex = atomIndex * pmeOrder + orderIndex;
925             GMX_ASSERT(cpuValueIndex < atomCount * pmeOrder,
926                        "Atom spline data index out of bounds (while transforming GPU data layout "
927                        "for host)");
928             switch (transform)
929             {
930                 case PmeLayoutTransform::GpuToHost:
931                     cpuSplineBuffer[cpuValueIndex] = h_splineBuffer[gpuValueIndex];
932                     break;
933
934                 case PmeLayoutTransform::HostToGpu:
935                     h_splineBuffer[gpuValueIndex] = cpuSplineBuffer[cpuValueIndex];
936                     break;
937
938                 default: GMX_THROW(gmx::InternalError("Unknown layout transform"));
939             }
940         }
941     }
942 }
943
944 void pme_gpu_get_real_grid_sizes(const PmeGpu* pmeGpu, gmx::IVec* gridSize, gmx::IVec* paddedGridSize)
945 {
946     GMX_ASSERT(gridSize != nullptr, "");
947     GMX_ASSERT(paddedGridSize != nullptr, "");
948     GMX_ASSERT(pmeGpu != nullptr, "");
949     auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
950     for (int i = 0; i < DIM; i++)
951     {
952         (*gridSize)[i]       = kernelParamsPtr->grid.realGridSize[i];
953         (*paddedGridSize)[i] = kernelParamsPtr->grid.realGridSizePadded[i];
954     }
955 }
956
957 void pme_gpu_reinit(gmx_pme_t* pme, const gmx_device_info_t* gpuInfo, PmeGpuProgramHandle pmeGpuProgram)
958 {
959     if (!pme_gpu_active(pme))
960     {
961         return;
962     }
963
964     if (!pme->gpu)
965     {
966         /* First-time initialization */
967         pme_gpu_init(pme, gpuInfo, pmeGpuProgram);
968     }
969     else
970     {
971         /* After this call nothing in the GPU code should refer to the gmx_pme_t *pme itself - until the next pme_gpu_reinit */
972         pme_gpu_copy_common_data_from(pme);
973     }
974     /* GPU FFT will only get used for a single rank.*/
975     pme->gpu->settings.performGPUFFT =
976             (pme->gpu->common->runMode == PmeRunMode::GPU) && !pme_gpu_uses_dd(pme->gpu);
977     pme->gpu->settings.performGPUSolve = (pme->gpu->common->runMode == PmeRunMode::GPU);
978
979     /* Reinit active timers */
980     pme_gpu_reinit_timings(pme->gpu);
981
982     pme_gpu_reinit_grids(pme->gpu);
983     // Note: if timing the reinit launch overhead becomes more relevant
984     // (e.g. with regulat PP-PME re-balancing), we should pass wcycle here.
985     pme_gpu_reinit_computation(pme, nullptr);
986     /* Clear the previous box - doesn't hurt, and forces the PME CPU recipbox
987      * update for mixed mode on grid switch. TODO: use shared recipbox field.
988      */
989     std::memset(pme->gpu->common->previousBox, 0, sizeof(pme->gpu->common->previousBox));
990     pme_gpu_select_best_performing_pme_spreadgather_kernels(pme->gpu);
991 }
992
993 void pme_gpu_destroy(PmeGpu* pmeGpu)
994 {
995     /* Free lots of data */
996     pme_gpu_free_energy_virial(pmeGpu);
997     pme_gpu_free_bspline_values(pmeGpu);
998     pme_gpu_free_forces(pmeGpu);
999     pme_gpu_free_coefficients(pmeGpu);
1000     pme_gpu_free_spline_data(pmeGpu);
1001     pme_gpu_free_grid_indices(pmeGpu);
1002     pme_gpu_free_fract_shifts(pmeGpu);
1003     pme_gpu_free_grids(pmeGpu);
1004
1005     pme_gpu_destroy_3dfft(pmeGpu);
1006
1007     /* Free the GPU-framework specific data last */
1008     pme_gpu_destroy_specific(pmeGpu);
1009
1010     delete pmeGpu;
1011 }
1012
1013 void pme_gpu_reinit_atoms(PmeGpu* pmeGpu, const int nAtoms, const real* charges)
1014 {
1015     auto* kernelParamsPtr         = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
1016     kernelParamsPtr->atoms.nAtoms = nAtoms;
1017     const int alignment           = pme_gpu_get_atom_data_alignment(pmeGpu);
1018     pmeGpu->nAtomsPadded          = ((nAtoms + alignment - 1) / alignment) * alignment;
1019     const int  nAtomsAlloc        = c_usePadding ? pmeGpu->nAtomsPadded : nAtoms;
1020     const bool haveToRealloc =
1021             (pmeGpu->nAtomsAlloc < nAtomsAlloc); /* This check might be redundant, but is logical */
1022     pmeGpu->nAtomsAlloc = nAtomsAlloc;
1023
1024 #if GMX_DOUBLE
1025     GMX_RELEASE_ASSERT(false, "Only single precision supported");
1026     GMX_UNUSED_VALUE(charges);
1027 #else
1028     pme_gpu_realloc_and_copy_input_coefficients(pmeGpu, reinterpret_cast<const float*>(charges));
1029     /* Could also be checked for haveToRealloc, but the copy always needs to be performed */
1030 #endif
1031
1032     if (haveToRealloc)
1033     {
1034         pme_gpu_realloc_forces(pmeGpu);
1035         pme_gpu_realloc_spline_data(pmeGpu);
1036         pme_gpu_realloc_grid_indices(pmeGpu);
1037     }
1038 }
1039
1040 void pme_gpu_3dfft(const PmeGpu* pmeGpu, gmx_fft_direction dir, int grid_index)
1041 {
1042     int timerId = (dir == GMX_FFT_REAL_TO_COMPLEX) ? gtPME_FFT_R2C : gtPME_FFT_C2R;
1043
1044     pme_gpu_start_timing(pmeGpu, timerId);
1045     pmeGpu->archSpecific->fftSetup[grid_index]->perform3dFft(
1046             dir, pme_gpu_fetch_timing_event(pmeGpu, timerId));
1047     pme_gpu_stop_timing(pmeGpu, timerId);
1048 }
1049
1050 /*! \brief
1051  * Given possibly large \p blockCount, returns a compact 1D or 2D grid for kernel scheduling,
1052  * to minimize number of unused blocks.
1053  */
1054 std::pair<int, int> inline pmeGpuCreateGrid(const PmeGpu* pmeGpu, int blockCount)
1055 {
1056     // How many maximum widths in X do we need (hopefully just one)
1057     const int minRowCount = (blockCount + pmeGpu->maxGridWidthX - 1) / pmeGpu->maxGridWidthX;
1058     // Trying to make things even
1059     const int colCount = (blockCount + minRowCount - 1) / minRowCount;
1060     GMX_ASSERT((colCount * minRowCount - blockCount) >= 0, "pmeGpuCreateGrid: totally wrong");
1061     GMX_ASSERT((colCount * minRowCount - blockCount) < minRowCount,
1062                "pmeGpuCreateGrid: excessive blocks");
1063     return std::pair<int, int>(colCount, minRowCount);
1064 }
1065
1066 /*! \brief
1067  * Returns a pointer to appropriate spline and spread kernel based on the input bool values
1068  *
1069  * \param[in]  pmeGpu                   The PME GPU structure.
1070  * \param[in]  useOrderThreadsPerAtom   bool controlling if we should use order or order*order threads per atom
1071  * \param[in]  writeSplinesToGlobal     bool controlling if we should write spline data to global memory
1072  *
1073  * \return Pointer to CUDA kernel
1074  */
1075 static auto selectSplineAndSpreadKernelPtr(const PmeGpu* pmeGpu, bool useOrderThreadsPerAtom, bool writeSplinesToGlobal)
1076 {
1077     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1078     if (writeSplinesToGlobal)
1079     {
1080         if (useOrderThreadsPerAtom)
1081         {
1082             kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesThPerAtom4;
1083         }
1084         else
1085         {
1086             kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplines;
1087         }
1088     }
1089     else
1090     {
1091         if (useOrderThreadsPerAtom)
1092         {
1093             kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4;
1094         }
1095         else
1096         {
1097             kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernel;
1098         }
1099     }
1100
1101     return kernelPtr;
1102 }
1103
1104 /*! \brief
1105  * Returns a pointer to appropriate spline kernel based on the input bool values
1106  *
1107  * \param[in]  pmeGpu                   The PME GPU structure.
1108  * \param[in]  useOrderThreadsPerAtom   bool controlling if we should use order or order*order threads per atom
1109  * \param[in]  writeSplinesToGlobal     bool controlling if we should write spline data to global memory
1110  *
1111  * \return Pointer to CUDA kernel
1112  */
1113 static auto selectSplineKernelPtr(const PmeGpu* pmeGpu, bool useOrderThreadsPerAtom, bool gmx_unused writeSplinesToGlobal)
1114 {
1115     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1116     GMX_ASSERT(
1117             writeSplinesToGlobal,
1118             "Spline data should always be written to global memory when just calculating splines");
1119
1120     if (useOrderThreadsPerAtom)
1121     {
1122         kernelPtr = pmeGpu->programHandle_->impl_->splineKernelThPerAtom4;
1123     }
1124     else
1125     {
1126         kernelPtr = pmeGpu->programHandle_->impl_->splineKernel;
1127     }
1128     return kernelPtr;
1129 }
1130
1131 /*! \brief
1132  * Returns a pointer to appropriate spread kernel based on the input bool values
1133  *
1134  * \param[in]  pmeGpu                   The PME GPU structure.
1135  * \param[in]  useOrderThreadsPerAtom   bool controlling if we should use order or order*order threads per atom
1136  * \param[in]  writeSplinesToGlobal     bool controlling if we should write spline data to global memory
1137  *
1138  * \return Pointer to CUDA kernel
1139  */
1140 static auto selectSpreadKernelPtr(const PmeGpu* pmeGpu, bool useOrderThreadsPerAtom, bool writeSplinesToGlobal)
1141 {
1142     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1143     if (writeSplinesToGlobal)
1144     {
1145         if (useOrderThreadsPerAtom)
1146         {
1147             kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelThPerAtom4;
1148         }
1149         else
1150         {
1151             kernelPtr = pmeGpu->programHandle_->impl_->spreadKernel;
1152         }
1153     }
1154     else
1155     {
1156         /* if we are not saving the spline data we need to recalculate it
1157            using the spline and spread Kernel */
1158         if (useOrderThreadsPerAtom)
1159         {
1160             kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4;
1161         }
1162         else
1163         {
1164             kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernel;
1165         }
1166     }
1167     return kernelPtr;
1168 }
1169
1170 void pme_gpu_spread(const PmeGpu*         pmeGpu,
1171                     GpuEventSynchronizer* xReadyOnDevice,
1172                     int gmx_unused gridIndex,
1173                     real*          h_grid,
1174                     bool           computeSplines,
1175                     bool           spreadCharges)
1176 {
1177     GMX_ASSERT(computeSplines || spreadCharges,
1178                "PME spline/spread kernel has invalid input (nothing to do)");
1179     const auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1180     GMX_ASSERT(kernelParamsPtr->atoms.nAtoms > 0, "No atom data in PME GPU spread");
1181
1182     const size_t blockSize = pmeGpu->programHandle_->impl_->spreadWorkGroupSize;
1183
1184     const int order = pmeGpu->common->pme_order;
1185     GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1186     const bool writeGlobal            = pmeGpu->settings.copyAllOutputs;
1187     const bool useOrderThreadsPerAtom = pmeGpu->settings.useOrderThreadsPerAtom;
1188     const bool recalculateSplines     = pmeGpu->settings.recalculateSplines;
1189 #if GMX_GPU == GMX_GPU_OPENCL
1190     GMX_ASSERT(!useOrderThreadsPerAtom, "Only 16 threads per atom supported in OpenCL");
1191     GMX_ASSERT(!recalculateSplines, "Recalculating splines not supported in OpenCL");
1192 #endif
1193     const int atomsPerBlock = useOrderThreadsPerAtom ? blockSize / c_pmeSpreadGatherThreadsPerAtom4ThPerAtom
1194                                                      : blockSize / c_pmeSpreadGatherThreadsPerAtom;
1195
1196     // TODO: pick smaller block size in runtime if needed
1197     // (e.g. on 660 Ti where 50% occupancy is ~25% faster than 100% occupancy with RNAse (~17.8k atoms))
1198     // If doing so, change atomsPerBlock in the kernels as well.
1199     // TODO: test varying block sizes on modern arch-s as well
1200     // TODO: also consider using cudaFuncSetCacheConfig() for preferring shared memory on older architectures
1201     //(for spline data mostly)
1202     GMX_ASSERT(!c_usePadding || !(c_pmeAtomDataAlignment % atomsPerBlock),
1203                "inconsistent atom data padding vs. spreading block size");
1204
1205     // Ensure that coordinates are ready on the device before launching spread;
1206     // only needed with CUDA on PP+PME ranks, not on separate PME ranks, in unit tests
1207     // nor in OpenCL as these cases use a single stream (hence xReadyOnDevice == nullptr).
1208     GMX_ASSERT(xReadyOnDevice != nullptr || (GMX_GPU != GMX_GPU_CUDA)
1209                        || pmeGpu->common->isRankPmeOnly || pme_gpu_is_testing(pmeGpu),
1210                "Need a valid coordinate synchronizer on PP+PME ranks with CUDA.");
1211     if (xReadyOnDevice)
1212     {
1213         xReadyOnDevice->enqueueWaitEvent(pmeGpu->archSpecific->pmeStream);
1214     }
1215
1216     const int blockCount = pmeGpu->nAtomsPadded / atomsPerBlock;
1217     auto      dimGrid    = pmeGpuCreateGrid(pmeGpu, blockCount);
1218
1219     KernelLaunchConfig config;
1220     config.blockSize[0] = order;
1221     config.blockSize[1] = useOrderThreadsPerAtom ? 1 : order;
1222     config.blockSize[2] = atomsPerBlock;
1223     config.gridSize[0]  = dimGrid.first;
1224     config.gridSize[1]  = dimGrid.second;
1225     config.stream       = pmeGpu->archSpecific->pmeStream;
1226
1227     int                                timingId;
1228     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1229     if (computeSplines)
1230     {
1231         if (spreadCharges)
1232         {
1233             timingId  = gtPME_SPLINEANDSPREAD;
1234             kernelPtr = selectSplineAndSpreadKernelPtr(pmeGpu, useOrderThreadsPerAtom,
1235                                                        writeGlobal || (!recalculateSplines));
1236         }
1237         else
1238         {
1239             timingId  = gtPME_SPLINE;
1240             kernelPtr = selectSplineKernelPtr(pmeGpu, useOrderThreadsPerAtom,
1241                                               writeGlobal || (!recalculateSplines));
1242         }
1243     }
1244     else
1245     {
1246         timingId  = gtPME_SPREAD;
1247         kernelPtr = selectSpreadKernelPtr(pmeGpu, useOrderThreadsPerAtom,
1248                                           writeGlobal || (!recalculateSplines));
1249     }
1250
1251
1252     pme_gpu_start_timing(pmeGpu, timingId);
1253     auto* timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1254 #if c_canEmbedBuffers
1255     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1256 #else
1257     const auto kernelArgs = prepareGpuKernelArguments(
1258             kernelPtr, config, kernelParamsPtr, &kernelParamsPtr->atoms.d_theta,
1259             &kernelParamsPtr->atoms.d_dtheta, &kernelParamsPtr->atoms.d_gridlineIndices,
1260             &kernelParamsPtr->grid.d_realGrid, &kernelParamsPtr->grid.d_fractShiftsTable,
1261             &kernelParamsPtr->grid.d_gridlineIndicesTable, &kernelParamsPtr->atoms.d_coefficients,
1262             &kernelParamsPtr->atoms.d_coordinates);
1263 #endif
1264
1265     launchGpuKernel(kernelPtr, config, timingEvent, "PME spline/spread", kernelArgs);
1266     pme_gpu_stop_timing(pmeGpu, timingId);
1267
1268     const bool copyBackGrid =
1269             spreadCharges && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_FFT(pmeGpu));
1270     if (copyBackGrid)
1271     {
1272         pme_gpu_copy_output_spread_grid(pmeGpu, h_grid);
1273     }
1274     const bool copyBackAtomData =
1275             computeSplines && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_gather(pmeGpu));
1276     if (copyBackAtomData)
1277     {
1278         pme_gpu_copy_output_spread_atom_data(pmeGpu);
1279     }
1280 }
1281
1282 void pme_gpu_solve(const PmeGpu* pmeGpu, t_complex* h_grid, GridOrdering gridOrdering, bool computeEnergyAndVirial)
1283 {
1284     const bool copyInputAndOutputGrid = pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_FFT(pmeGpu);
1285
1286     auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1287
1288     float* h_gridFloat = reinterpret_cast<float*>(h_grid);
1289     if (copyInputAndOutputGrid)
1290     {
1291         copyToDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid, h_gridFloat, 0,
1292                            pmeGpu->archSpecific->complexGridSize, pmeGpu->archSpecific->pmeStream,
1293                            pmeGpu->settings.transferKind, nullptr);
1294     }
1295
1296     int majorDim = -1, middleDim = -1, minorDim = -1;
1297     switch (gridOrdering)
1298     {
1299         case GridOrdering::YZX:
1300             majorDim  = YY;
1301             middleDim = ZZ;
1302             minorDim  = XX;
1303             break;
1304
1305         case GridOrdering::XYZ:
1306             majorDim  = XX;
1307             middleDim = YY;
1308             minorDim  = ZZ;
1309             break;
1310
1311         default: GMX_ASSERT(false, "Implement grid ordering here and below for the kernel launch");
1312     }
1313
1314     const int maxBlockSize = pmeGpu->programHandle_->impl_->solveMaxWorkGroupSize;
1315
1316     const int gridLineSize      = pmeGpu->kernelParams->grid.complexGridSize[minorDim];
1317     const int gridLinesPerBlock = std::max(maxBlockSize / gridLineSize, 1);
1318     const int blocksPerGridLine = (gridLineSize + maxBlockSize - 1) / maxBlockSize;
1319     int       cellsPerBlock;
1320     if (blocksPerGridLine == 1)
1321     {
1322         cellsPerBlock = gridLineSize * gridLinesPerBlock;
1323     }
1324     else
1325     {
1326         cellsPerBlock = (gridLineSize + blocksPerGridLine - 1) / blocksPerGridLine;
1327     }
1328     const int warpSize  = pmeGpu->programHandle_->impl_->warpSize;
1329     const int blockSize = (cellsPerBlock + warpSize - 1) / warpSize * warpSize;
1330
1331     static_assert(GMX_GPU != GMX_GPU_CUDA || c_solveMaxWarpsPerBlock / 2 >= 4,
1332                   "The CUDA solve energy kernels needs at least 4 warps. "
1333                   "Here we launch at least half of the max warps.");
1334
1335     KernelLaunchConfig config;
1336     config.blockSize[0] = blockSize;
1337     config.gridSize[0]  = blocksPerGridLine;
1338     // rounding up to full warps so that shuffle operations produce defined results
1339     config.gridSize[1] = (pmeGpu->kernelParams->grid.complexGridSize[middleDim] + gridLinesPerBlock - 1)
1340                          / gridLinesPerBlock;
1341     config.gridSize[2] = pmeGpu->kernelParams->grid.complexGridSize[majorDim];
1342     config.stream      = pmeGpu->archSpecific->pmeStream;
1343
1344     int                                timingId  = gtPME_SOLVE;
1345     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1346     if (gridOrdering == GridOrdering::YZX)
1347     {
1348         kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveYZXEnergyKernel
1349                                            : pmeGpu->programHandle_->impl_->solveYZXKernel;
1350     }
1351     else if (gridOrdering == GridOrdering::XYZ)
1352     {
1353         kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveXYZEnergyKernel
1354                                            : pmeGpu->programHandle_->impl_->solveXYZKernel;
1355     }
1356
1357     pme_gpu_start_timing(pmeGpu, timingId);
1358     auto* timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1359 #if c_canEmbedBuffers
1360     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1361 #else
1362     const auto kernelArgs = prepareGpuKernelArguments(
1363             kernelPtr, config, kernelParamsPtr, &kernelParamsPtr->grid.d_splineModuli,
1364             &kernelParamsPtr->constants.d_virialAndEnergy, &kernelParamsPtr->grid.d_fourierGrid);
1365 #endif
1366     launchGpuKernel(kernelPtr, config, timingEvent, "PME solve", kernelArgs);
1367     pme_gpu_stop_timing(pmeGpu, timingId);
1368
1369     if (computeEnergyAndVirial)
1370     {
1371         copyFromDeviceBuffer(pmeGpu->staging.h_virialAndEnergy,
1372                              &kernelParamsPtr->constants.d_virialAndEnergy, 0, c_virialAndEnergyCount,
1373                              pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
1374     }
1375
1376     if (copyInputAndOutputGrid)
1377     {
1378         copyFromDeviceBuffer(h_gridFloat, &kernelParamsPtr->grid.d_fourierGrid, 0,
1379                              pmeGpu->archSpecific->complexGridSize, pmeGpu->archSpecific->pmeStream,
1380                              pmeGpu->settings.transferKind, nullptr);
1381     }
1382 }
1383
1384 /*! \brief
1385  * Returns a pointer to appropriate gather kernel based on the inputvalues
1386  *
1387  * \param[in]  pmeGpu                   The PME GPU structure.
1388  * \param[in]  useOrderThreadsPerAtom   bool controlling if we should use order or order*order threads per atom
1389  * \param[in]  readSplinesFromGlobal    bool controlling if we should write spline data to global memory
1390  * \param[in]  forceTreatment           Controls if the forces from the gather should increment or replace the input forces.
1391  *
1392  * \return Pointer to CUDA kernel
1393  */
1394 inline auto selectGatherKernelPtr(const PmeGpu*          pmeGpu,
1395                                   bool                   useOrderThreadsPerAtom,
1396                                   bool                   readSplinesFromGlobal,
1397                                   PmeForceOutputHandling forceTreatment)
1398
1399 {
1400     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1401
1402     if (readSplinesFromGlobal)
1403     {
1404         if (useOrderThreadsPerAtom)
1405         {
1406             kernelPtr = (forceTreatment == PmeForceOutputHandling::Set)
1407                                 ? pmeGpu->programHandle_->impl_->gatherKernelReadSplinesThPerAtom4
1408                                 : pmeGpu->programHandle_->impl_->gatherReduceWithInputKernelReadSplinesThPerAtom4;
1409         }
1410         else
1411         {
1412             kernelPtr = (forceTreatment == PmeForceOutputHandling::Set)
1413                                 ? pmeGpu->programHandle_->impl_->gatherKernelReadSplines
1414                                 : pmeGpu->programHandle_->impl_->gatherReduceWithInputKernelReadSplines;
1415         }
1416     }
1417     else
1418     {
1419         if (useOrderThreadsPerAtom)
1420         {
1421             kernelPtr = (forceTreatment == PmeForceOutputHandling::Set)
1422                                 ? pmeGpu->programHandle_->impl_->gatherKernelThPerAtom4
1423                                 : pmeGpu->programHandle_->impl_->gatherReduceWithInputKernelThPerAtom4;
1424         }
1425         else
1426         {
1427             kernelPtr = (forceTreatment == PmeForceOutputHandling::Set)
1428                                 ? pmeGpu->programHandle_->impl_->gatherKernel
1429                                 : pmeGpu->programHandle_->impl_->gatherReduceWithInputKernel;
1430         }
1431     }
1432     return kernelPtr;
1433 }
1434
1435
1436 void pme_gpu_gather(PmeGpu* pmeGpu, PmeForceOutputHandling forceTreatment, const float* h_grid)
1437 {
1438     /* Copying the input CPU forces for reduction */
1439     if (forceTreatment != PmeForceOutputHandling::Set)
1440     {
1441         pme_gpu_copy_input_forces(pmeGpu);
1442     }
1443
1444     if (!pme_gpu_performs_FFT(pmeGpu) || pme_gpu_is_testing(pmeGpu))
1445     {
1446         pme_gpu_copy_input_gather_grid(pmeGpu, const_cast<float*>(h_grid));
1447     }
1448
1449     if (pme_gpu_is_testing(pmeGpu))
1450     {
1451         pme_gpu_copy_input_gather_atom_data(pmeGpu);
1452     }
1453
1454     /* Set if we have unit tests */
1455     const bool   readGlobal             = pmeGpu->settings.copyAllOutputs;
1456     const size_t blockSize              = pmeGpu->programHandle_->impl_->gatherWorkGroupSize;
1457     const bool   useOrderThreadsPerAtom = pmeGpu->settings.useOrderThreadsPerAtom;
1458     const bool   recalculateSplines     = pmeGpu->settings.recalculateSplines;
1459 #if GMX_GPU == GMX_GPU_OPENCL
1460     GMX_ASSERT(!useOrderThreadsPerAtom, "Only 16 threads per atom supported in OpenCL");
1461     GMX_ASSERT(!recalculateSplines, "Recalculating splines not supported in OpenCL");
1462 #endif
1463     const int atomsPerBlock = useOrderThreadsPerAtom ? blockSize / c_pmeSpreadGatherThreadsPerAtom4ThPerAtom
1464                                                      : blockSize / c_pmeSpreadGatherThreadsPerAtom;
1465
1466     GMX_ASSERT(!c_usePadding || !(c_pmeAtomDataAlignment % atomsPerBlock),
1467                "inconsistent atom data padding vs. gathering block size");
1468
1469     const int blockCount = pmeGpu->nAtomsPadded / atomsPerBlock;
1470     auto      dimGrid    = pmeGpuCreateGrid(pmeGpu, blockCount);
1471
1472     const int order = pmeGpu->common->pme_order;
1473     GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1474
1475     KernelLaunchConfig config;
1476     config.blockSize[0] = order;
1477     config.blockSize[1] = useOrderThreadsPerAtom ? 1 : order;
1478     config.blockSize[2] = atomsPerBlock;
1479     config.gridSize[0]  = dimGrid.first;
1480     config.gridSize[1]  = dimGrid.second;
1481     config.stream       = pmeGpu->archSpecific->pmeStream;
1482
1483     // TODO test different cache configs
1484
1485     int                                timingId  = gtPME_GATHER;
1486     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = selectGatherKernelPtr(
1487             pmeGpu, useOrderThreadsPerAtom, readGlobal || (!recalculateSplines), forceTreatment);
1488     // TODO design kernel selection getters and make PmeGpu a friend of PmeGpuProgramImpl
1489
1490     pme_gpu_start_timing(pmeGpu, timingId);
1491     auto*       timingEvent     = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1492     const auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1493 #if c_canEmbedBuffers
1494     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1495 #else
1496     const auto kernelArgs = prepareGpuKernelArguments(
1497             kernelPtr, config, kernelParamsPtr, &kernelParamsPtr->atoms.d_coefficients,
1498             &kernelParamsPtr->grid.d_realGrid, &kernelParamsPtr->atoms.d_theta,
1499             &kernelParamsPtr->atoms.d_dtheta, &kernelParamsPtr->atoms.d_gridlineIndices,
1500             &kernelParamsPtr->atoms.d_forces);
1501 #endif
1502     launchGpuKernel(kernelPtr, config, timingEvent, "PME gather", kernelArgs);
1503     pme_gpu_stop_timing(pmeGpu, timingId);
1504
1505     if (pmeGpu->settings.useGpuForceReduction)
1506     {
1507         pmeGpu->archSpecific->pmeForcesReady.markEvent(pmeGpu->archSpecific->pmeStream);
1508     }
1509     else
1510     {
1511         pme_gpu_copy_output_forces(pmeGpu);
1512     }
1513 }
1514
1515 DeviceBuffer<float> pme_gpu_get_kernelparam_coordinates(const PmeGpu* pmeGpu)
1516 {
1517     GMX_ASSERT(pmeGpu && pmeGpu->kernelParams,
1518                "PME GPU device buffer was requested in non-GPU build or before the GPU PME was "
1519                "initialized.");
1520
1521     return pmeGpu->kernelParams->atoms.d_coordinates;
1522 }
1523
1524 void* pme_gpu_get_kernelparam_forces(const PmeGpu* pmeGpu)
1525 {
1526     if (pmeGpu && pmeGpu->kernelParams)
1527     {
1528         return pmeGpu->kernelParams->atoms.d_forces;
1529     }
1530     else
1531     {
1532         return nullptr;
1533     }
1534 }
1535
1536 /*! \brief Check the validity of the device buffer.
1537  *
1538  * Checks if the buffer is not nullptr and, when possible, if it is big enough.
1539  *
1540  * \todo Split and move this function to gpu_utils.
1541  *
1542  * \param[in] buffer        Device buffer to be checked.
1543  * \param[in] requiredSize  Number of elements that the buffer will have to accommodate.
1544  *
1545  * \returns If the device buffer can be set.
1546  */
1547 template<typename T>
1548 static bool checkDeviceBuffer(gmx_unused DeviceBuffer<T> buffer, gmx_unused int requiredSize)
1549 {
1550 #if GMX_GPU == GMX_GPU_CUDA
1551     GMX_ASSERT(buffer != nullptr, "The device pointer is nullptr");
1552     return buffer != nullptr;
1553 #elif GMX_GPU == GMX_GPU_OPENCL
1554     size_t size;
1555     int    retval = clGetMemObjectInfo(buffer, CL_MEM_SIZE, sizeof(size), &size, nullptr);
1556     GMX_ASSERT(retval == CL_SUCCESS,
1557                gmx::formatString("clGetMemObjectInfo failed with error code #%d", retval).c_str());
1558     GMX_ASSERT(static_cast<int>(size) >= requiredSize,
1559                "Number of atoms in device buffer is smaller then required size.");
1560     return retval == CL_SUCCESS && static_cast<int>(size) >= requiredSize;
1561 #elif GMX_GPU == GMX_GPU_NONE
1562     GMX_ASSERT(false, "Setter for device-side coordinates was called in non-GPU build.");
1563     return false;
1564 #endif
1565 }
1566
1567 void pme_gpu_set_kernelparam_coordinates(const PmeGpu* pmeGpu, DeviceBuffer<float> d_x)
1568 {
1569     GMX_ASSERT(pmeGpu && pmeGpu->kernelParams,
1570                "PME GPU device buffer can not be set in non-GPU builds or before the GPU PME was "
1571                "initialized.");
1572
1573     GMX_ASSERT(checkDeviceBuffer(d_x, pmeGpu->kernelParams->atoms.nAtoms),
1574                "The device-side buffer can not be set.");
1575
1576     pmeGpu->kernelParams->atoms.d_coordinates = d_x;
1577 }
1578
1579 void* pme_gpu_get_stream(const PmeGpu* pmeGpu)
1580 {
1581     if (pmeGpu)
1582     {
1583         return static_cast<void*>(&pmeGpu->archSpecific->pmeStream);
1584     }
1585     else
1586     {
1587         return nullptr;
1588     }
1589 }
1590
1591 void* pme_gpu_get_context(const PmeGpu* pmeGpu)
1592 {
1593     if (pmeGpu)
1594     {
1595         return static_cast<void*>(&pmeGpu->archSpecific->context);
1596     }
1597     else
1598     {
1599         return nullptr;
1600     }
1601 }
1602
1603 GpuEventSynchronizer* pme_gpu_get_forces_ready_synchronizer(const PmeGpu* pmeGpu)
1604 {
1605     if (pmeGpu && pmeGpu->kernelParams)
1606     {
1607         return &pmeGpu->archSpecific->pmeForcesReady;
1608     }
1609     else
1610     {
1611         return nullptr;
1612     }
1613 }