a8d98e96cf11f71a92630bc24d76d2a2b4bf22fc
[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 /*! \internal \brief
92  * Wrapper for getting a pointer to the plain C++ part of the GPU kernel parameters structure.
93  *
94  * \param[in] pmeGpu  The PME GPU structure.
95  * \returns The pointer to the kernel parameters.
96  */
97 static PmeGpuKernelParamsBase *pme_gpu_get_kernel_params_base_ptr(const PmeGpu *pmeGpu)
98 {
99     // reinterpret_cast is needed because the derived CUDA structure is not known in this file
100     auto *kernelParamsPtr = reinterpret_cast<PmeGpuKernelParamsBase *>(pmeGpu->kernelParams.get());
101     return kernelParamsPtr;
102 }
103
104 int pme_gpu_get_atom_data_alignment(const PmeGpu * /*unused*/)
105 {
106     //TODO: this can be simplified, as c_pmeAtomDataAlignment is now constant
107     return c_pmeAtomDataAlignment;
108 }
109
110 int pme_gpu_get_atoms_per_warp(const PmeGpu *pmeGpu)
111 {
112     return pmeGpu->programHandle_->impl_->warpSize / c_pmeSpreadGatherThreadsPerAtom;
113 }
114
115 void pme_gpu_synchronize(const PmeGpu *pmeGpu)
116 {
117     gpuStreamSynchronize(pmeGpu->archSpecific->pmeStream);
118 }
119
120 void pme_gpu_alloc_energy_virial(PmeGpu *pmeGpu)
121 {
122     const size_t energyAndVirialSize = c_virialAndEnergyCount * sizeof(float);
123     allocateDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy, c_virialAndEnergyCount, pmeGpu->archSpecific->context);
124     pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_virialAndEnergy), energyAndVirialSize);
125 }
126
127 void pme_gpu_free_energy_virial(PmeGpu *pmeGpu)
128 {
129     freeDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy);
130     pfree(pmeGpu->staging.h_virialAndEnergy);
131     pmeGpu->staging.h_virialAndEnergy = nullptr;
132 }
133
134 void pme_gpu_clear_energy_virial(const PmeGpu *pmeGpu)
135 {
136     clearDeviceBufferAsync(&pmeGpu->kernelParams->constants.d_virialAndEnergy, 0,
137                            c_virialAndEnergyCount, pmeGpu->archSpecific->pmeStream);
138 }
139
140 void pme_gpu_realloc_and_copy_bspline_values(PmeGpu *pmeGpu)
141 {
142     const int splineValuesOffset[DIM] = {
143         0,
144         pmeGpu->kernelParams->grid.realGridSize[XX],
145         pmeGpu->kernelParams->grid.realGridSize[XX] + pmeGpu->kernelParams->grid.realGridSize[YY]
146     };
147     memcpy(&pmeGpu->kernelParams->grid.splineValuesOffset, &splineValuesOffset, sizeof(splineValuesOffset));
148
149     const int newSplineValuesSize = pmeGpu->kernelParams->grid.realGridSize[XX] +
150         pmeGpu->kernelParams->grid.realGridSize[YY] +
151         pmeGpu->kernelParams->grid.realGridSize[ZZ];
152     const bool shouldRealloc = (newSplineValuesSize > pmeGpu->archSpecific->splineValuesSize);
153     reallocateDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli, newSplineValuesSize,
154                            &pmeGpu->archSpecific->splineValuesSize, &pmeGpu->archSpecific->splineValuesSizeAlloc, pmeGpu->archSpecific->context);
155     if (shouldRealloc)
156     {
157         /* Reallocate the host buffer */
158         pfree(pmeGpu->staging.h_splineModuli);
159         pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_splineModuli), newSplineValuesSize * sizeof(float));
160     }
161     for (int i = 0; i < DIM; i++)
162     {
163         memcpy(pmeGpu->staging.h_splineModuli + splineValuesOffset[i], pmeGpu->common->bsp_mod[i].data(), pmeGpu->common->bsp_mod[i].size() * sizeof(float));
164     }
165     /* TODO: pin original buffer instead! */
166     copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli, pmeGpu->staging.h_splineModuli,
167                        0, newSplineValuesSize,
168                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
169 }
170
171 void pme_gpu_free_bspline_values(const PmeGpu *pmeGpu)
172 {
173     pfree(pmeGpu->staging.h_splineModuli);
174     freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli);
175 }
176
177 void pme_gpu_realloc_forces(PmeGpu *pmeGpu)
178 {
179     const size_t newForcesSize = pmeGpu->nAtomsAlloc * DIM;
180     GMX_ASSERT(newForcesSize > 0, "Bad number of atoms in PME GPU");
181     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces, newForcesSize,
182                            &pmeGpu->archSpecific->forcesSize, &pmeGpu->archSpecific->forcesSizeAlloc, pmeGpu->archSpecific->context);
183     pmeGpu->staging.h_forces.reserveWithPadding(pmeGpu->nAtomsAlloc);
184     pmeGpu->staging.h_forces.resizeWithPadding(pmeGpu->kernelParams->atoms.nAtoms);
185 }
186
187 void pme_gpu_free_forces(const PmeGpu *pmeGpu)
188 {
189     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces);
190 }
191
192 void pme_gpu_copy_input_forces(PmeGpu *pmeGpu)
193 {
194     GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
195     float *h_forcesFloat = reinterpret_cast<float *>(pmeGpu->staging.h_forces.data());
196     copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces, h_forcesFloat,
197                        0, DIM * pmeGpu->kernelParams->atoms.nAtoms,
198                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
199 }
200
201 void pme_gpu_copy_output_forces(PmeGpu *pmeGpu)
202 {
203     GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
204     float *h_forcesFloat = reinterpret_cast<float *>(pmeGpu->staging.h_forces.data());
205     copyFromDeviceBuffer(h_forcesFloat, &pmeGpu->kernelParams->atoms.d_forces,
206                          0, DIM * pmeGpu->kernelParams->atoms.nAtoms,
207                          pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
208 }
209
210 void pme_gpu_realloc_coordinates(const PmeGpu *pmeGpu)
211 {
212     const size_t newCoordinatesSize = pmeGpu->nAtomsAlloc * DIM;
213     GMX_ASSERT(newCoordinatesSize > 0, "Bad number of atoms in PME GPU");
214     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates, newCoordinatesSize,
215                            &pmeGpu->archSpecific->coordinatesSize, &pmeGpu->archSpecific->coordinatesSizeAlloc, pmeGpu->archSpecific->context);
216     if (c_usePadding)
217     {
218         const size_t paddingIndex = DIM * pmeGpu->kernelParams->atoms.nAtoms;
219         const size_t paddingCount = DIM * pmeGpu->nAtomsAlloc - paddingIndex;
220         if (paddingCount > 0)
221         {
222             clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coordinates, paddingIndex,
223                                    paddingCount, pmeGpu->archSpecific->pmeStream);
224         }
225     }
226 }
227
228 void pme_gpu_copy_input_coordinates(const PmeGpu *pmeGpu, const rvec *h_coordinates)
229 {
230     GMX_ASSERT(h_coordinates, "Bad host-side coordinate buffer in PME GPU");
231 #if GMX_DOUBLE
232     GMX_RELEASE_ASSERT(false, "Only single precision is supported");
233     GMX_UNUSED_VALUE(h_coordinates);
234 #else
235     const float *h_coordinatesFloat = reinterpret_cast<const float *>(h_coordinates);
236     copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates, h_coordinatesFloat,
237                        0, pmeGpu->kernelParams->atoms.nAtoms * DIM,
238                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
239 #endif
240 }
241
242 void pme_gpu_free_coordinates(const PmeGpu *pmeGpu)
243 {
244     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates);
245 }
246
247 void pme_gpu_realloc_and_copy_input_coefficients(const PmeGpu *pmeGpu, const float *h_coefficients)
248 {
249     GMX_ASSERT(h_coefficients, "Bad host-side charge buffer in PME GPU");
250     const size_t newCoefficientsSize = pmeGpu->nAtomsAlloc;
251     GMX_ASSERT(newCoefficientsSize > 0, "Bad number of atoms in PME GPU");
252     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients, newCoefficientsSize,
253                            &pmeGpu->archSpecific->coefficientsSize, &pmeGpu->archSpecific->coefficientsSizeAlloc, pmeGpu->archSpecific->context);
254     copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients, const_cast<float *>(h_coefficients),
255                        0, pmeGpu->kernelParams->atoms.nAtoms,
256                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
257     if (c_usePadding)
258     {
259         const size_t paddingIndex = pmeGpu->kernelParams->atoms.nAtoms;
260         const size_t paddingCount = pmeGpu->nAtomsAlloc - paddingIndex;
261         if (paddingCount > 0)
262         {
263             clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coefficients, paddingIndex,
264                                    paddingCount, pmeGpu->archSpecific->pmeStream);
265         }
266     }
267 }
268
269 void pme_gpu_free_coefficients(const PmeGpu *pmeGpu)
270 {
271     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients);
272 }
273
274 void pme_gpu_realloc_spline_data(PmeGpu *pmeGpu)
275 {
276     const int    order             = pmeGpu->common->pme_order;
277     const int    alignment         = pme_gpu_get_atoms_per_warp(pmeGpu);
278     const size_t nAtomsPadded      = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
279     const int    newSplineDataSize = DIM * order * nAtomsPadded;
280     GMX_ASSERT(newSplineDataSize > 0, "Bad number of atoms in PME GPU");
281     /* Two arrays of the same size */
282     const bool shouldRealloc        = (newSplineDataSize > pmeGpu->archSpecific->splineDataSize);
283     int        currentSizeTemp      = pmeGpu->archSpecific->splineDataSize;
284     int        currentSizeTempAlloc = pmeGpu->archSpecific->splineDataSizeAlloc;
285     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta, newSplineDataSize,
286                            &currentSizeTemp, &currentSizeTempAlloc, pmeGpu->archSpecific->context);
287     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta, newSplineDataSize,
288                            &pmeGpu->archSpecific->splineDataSize, &pmeGpu->archSpecific->splineDataSizeAlloc, pmeGpu->archSpecific->context);
289     // the host side reallocation
290     if (shouldRealloc)
291     {
292         pfree(pmeGpu->staging.h_theta);
293         pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_theta), newSplineDataSize * sizeof(float));
294         pfree(pmeGpu->staging.h_dtheta);
295         pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_dtheta), newSplineDataSize * sizeof(float));
296     }
297 }
298
299 void pme_gpu_free_spline_data(const PmeGpu *pmeGpu)
300 {
301     /* Two arrays of the same size */
302     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta);
303     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta);
304     pfree(pmeGpu->staging.h_theta);
305     pfree(pmeGpu->staging.h_dtheta);
306 }
307
308 void pme_gpu_realloc_grid_indices(PmeGpu *pmeGpu)
309 {
310     const size_t newIndicesSize = DIM * pmeGpu->nAtomsAlloc;
311     GMX_ASSERT(newIndicesSize > 0, "Bad number of atoms in PME GPU");
312     reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices, newIndicesSize,
313                            &pmeGpu->archSpecific->gridlineIndicesSize, &pmeGpu->archSpecific->gridlineIndicesSizeAlloc, pmeGpu->archSpecific->context);
314     pfree(pmeGpu->staging.h_gridlineIndices);
315     pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_gridlineIndices), newIndicesSize * sizeof(int));
316 }
317
318 void pme_gpu_free_grid_indices(const PmeGpu *pmeGpu)
319 {
320     freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices);
321     pfree(pmeGpu->staging.h_gridlineIndices);
322 }
323
324 void pme_gpu_realloc_grids(PmeGpu *pmeGpu)
325 {
326     auto     *kernelParamsPtr = pmeGpu->kernelParams.get();
327     const int newRealGridSize = kernelParamsPtr->grid.realGridSizePadded[XX] *
328         kernelParamsPtr->grid.realGridSizePadded[YY] *
329         kernelParamsPtr->grid.realGridSizePadded[ZZ];
330     const int newComplexGridSize = kernelParamsPtr->grid.complexGridSizePadded[XX] *
331         kernelParamsPtr->grid.complexGridSizePadded[YY] *
332         kernelParamsPtr->grid.complexGridSizePadded[ZZ] * 2;
333     // Multiplied by 2 because we count complex grid size for complex numbers, but all allocations/pointers are float
334     if (pmeGpu->archSpecific->performOutOfPlaceFFT)
335     {
336         /* 2 separate grids */
337         reallocateDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid, newComplexGridSize,
338                                &pmeGpu->archSpecific->complexGridSize, &pmeGpu->archSpecific->complexGridSizeAlloc, pmeGpu->archSpecific->context);
339         reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid, newRealGridSize,
340                                &pmeGpu->archSpecific->realGridSize, &pmeGpu->archSpecific->realGridSizeAlloc, pmeGpu->archSpecific->context);
341     }
342     else
343     {
344         /* A single buffer so that any grid will fit */
345         const int newGridsSize = std::max(newRealGridSize, newComplexGridSize);
346         reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid, newGridsSize,
347                                &pmeGpu->archSpecific->realGridSize, &pmeGpu->archSpecific->realGridSizeAlloc, pmeGpu->archSpecific->context);
348         kernelParamsPtr->grid.d_fourierGrid   = kernelParamsPtr->grid.d_realGrid;
349         pmeGpu->archSpecific->complexGridSize = pmeGpu->archSpecific->realGridSize;
350         // the size might get used later for copying the grid
351     }
352 }
353
354 void pme_gpu_free_grids(const PmeGpu *pmeGpu)
355 {
356     if (pmeGpu->archSpecific->performOutOfPlaceFFT)
357     {
358         freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_fourierGrid);
359     }
360     freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid);
361 }
362
363 void pme_gpu_clear_grids(const PmeGpu *pmeGpu)
364 {
365     clearDeviceBufferAsync(&pmeGpu->kernelParams->grid.d_realGrid, 0,
366                            pmeGpu->archSpecific->realGridSize, pmeGpu->archSpecific->pmeStream);
367 }
368
369 void pme_gpu_realloc_and_copy_fract_shifts(PmeGpu *pmeGpu)
370 {
371     pme_gpu_free_fract_shifts(pmeGpu);
372
373     auto        *kernelParamsPtr = pmeGpu->kernelParams.get();
374
375     const int    nx                  = kernelParamsPtr->grid.realGridSize[XX];
376     const int    ny                  = kernelParamsPtr->grid.realGridSize[YY];
377     const int    nz                  = kernelParamsPtr->grid.realGridSize[ZZ];
378     const int    cellCount           = c_pmeNeighborUnitcellCount;
379     const int    gridDataOffset[DIM] = {0, cellCount * nx, cellCount * (nx + ny)};
380
381     memcpy(kernelParamsPtr->grid.tablesOffsets, &gridDataOffset, sizeof(gridDataOffset));
382
383     const int    newFractShiftsSize  = cellCount * (nx + ny + nz);
384
385 #if GMX_GPU == GMX_GPU_CUDA
386     initParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable,
387                          kernelParamsPtr->fractShiftsTableTexture,
388                          pmeGpu->common->fsh.data(),
389                          newFractShiftsSize);
390
391     initParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
392                          kernelParamsPtr->gridlineIndicesTableTexture,
393                          pmeGpu->common->nn.data(),
394                          newFractShiftsSize);
395 #elif GMX_GPU == GMX_GPU_OPENCL
396     // No dedicated texture routines....
397     allocateDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable, newFractShiftsSize, pmeGpu->archSpecific->context);
398     allocateDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable, newFractShiftsSize, pmeGpu->archSpecific->context);
399     copyToDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable, pmeGpu->common->fsh.data(),
400                        0, newFractShiftsSize,
401                        pmeGpu->archSpecific->pmeStream, GpuApiCallBehavior::Async, nullptr);
402     copyToDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable, pmeGpu->common->nn.data(),
403                        0, newFractShiftsSize,
404                        pmeGpu->archSpecific->pmeStream, GpuApiCallBehavior::Async, nullptr);
405 #endif
406 }
407
408 void pme_gpu_free_fract_shifts(const PmeGpu *pmeGpu)
409 {
410     auto *kernelParamsPtr = pmeGpu->kernelParams.get();
411 #if GMX_GPU == GMX_GPU_CUDA
412     destroyParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable,
413                             kernelParamsPtr->fractShiftsTableTexture);
414     destroyParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
415                             kernelParamsPtr->gridlineIndicesTableTexture);
416 #elif GMX_GPU == GMX_GPU_OPENCL
417     freeDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable);
418     freeDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable);
419 #endif
420 }
421
422 bool pme_gpu_stream_query(const PmeGpu *pmeGpu)
423 {
424     return haveStreamTasksCompleted(pmeGpu->archSpecific->pmeStream);
425 }
426
427 void pme_gpu_copy_input_gather_grid(const PmeGpu *pmeGpu, float *h_grid)
428 {
429     copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid, h_grid,
430                        0, pmeGpu->archSpecific->realGridSize,
431                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
432 }
433
434 void pme_gpu_copy_output_spread_grid(const PmeGpu *pmeGpu, float *h_grid)
435 {
436     copyFromDeviceBuffer(h_grid, &pmeGpu->kernelParams->grid.d_realGrid,
437                          0, pmeGpu->archSpecific->realGridSize,
438                          pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
439     pmeGpu->archSpecific->syncSpreadGridD2H.markEvent(pmeGpu->archSpecific->pmeStream);
440 }
441
442 void pme_gpu_copy_output_spread_atom_data(const PmeGpu *pmeGpu)
443 {
444     const int    alignment        = pme_gpu_get_atoms_per_warp(pmeGpu);
445     const size_t nAtomsPadded     = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
446     const size_t splinesCount     = DIM * nAtomsPadded * pmeGpu->common->pme_order;
447     auto        *kernelParamsPtr  = pmeGpu->kernelParams.get();
448     copyFromDeviceBuffer(pmeGpu->staging.h_dtheta, &kernelParamsPtr->atoms.d_dtheta,
449                          0, splinesCount,
450                          pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
451     copyFromDeviceBuffer(pmeGpu->staging.h_theta, &kernelParamsPtr->atoms.d_theta,
452                          0, splinesCount,
453                          pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
454     copyFromDeviceBuffer(pmeGpu->staging.h_gridlineIndices, &kernelParamsPtr->atoms.d_gridlineIndices,
455                          0, kernelParamsPtr->atoms.nAtoms * DIM,
456                          pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
457 }
458
459 void pme_gpu_copy_input_gather_atom_data(const PmeGpu *pmeGpu)
460 {
461     const int    alignment       = pme_gpu_get_atoms_per_warp(pmeGpu);
462     const size_t nAtomsPadded    = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
463     const size_t splinesCount    = DIM * nAtomsPadded * pmeGpu->common->pme_order;
464     auto        *kernelParamsPtr = pmeGpu->kernelParams.get();
465     if (c_usePadding)
466     {
467         // TODO: could clear only the padding and not the whole thing, but this is a test-exclusive code anyway
468         clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_gridlineIndices, 0,
469                                pmeGpu->nAtomsAlloc * DIM, pmeGpu->archSpecific->pmeStream);
470         clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_dtheta, 0,
471                                pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM, pmeGpu->archSpecific->pmeStream);
472         clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_theta, 0,
473                                pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM, pmeGpu->archSpecific->pmeStream);
474     }
475     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_dtheta, pmeGpu->staging.h_dtheta,
476                        0, splinesCount,
477                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
478     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_theta, pmeGpu->staging.h_theta,
479                        0, splinesCount,
480                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
481     copyToDeviceBuffer(&kernelParamsPtr->atoms.d_gridlineIndices, pmeGpu->staging.h_gridlineIndices,
482                        0, kernelParamsPtr->atoms.nAtoms * DIM,
483                        pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
484 }
485
486 void pme_gpu_sync_spread_grid(const PmeGpu *pmeGpu)
487 {
488     pmeGpu->archSpecific->syncSpreadGridD2H.waitForEvent();
489 }
490
491 void pme_gpu_init_internal(PmeGpu *pmeGpu)
492 {
493 #if GMX_GPU == GMX_GPU_CUDA
494     // Prepare to use the device that this PME task was assigned earlier.
495     // Other entities, such as CUDA timing events, are known to implicitly use the device context.
496     CU_RET_ERR(cudaSetDevice(pmeGpu->deviceInfo->id), "Switching to PME CUDA device");
497 #endif
498
499     /* Allocate the target-specific structures */
500     pmeGpu->archSpecific.reset(new PmeGpuSpecific());
501     pmeGpu->kernelParams.reset(new PmeGpuKernelParams());
502
503     pmeGpu->archSpecific->performOutOfPlaceFFT = true;
504     /* This should give better performance, according to the cuFFT documentation.
505      * The performance seems to be the same though.
506      * TODO: PME could also try to pick up nice grid sizes (with factors of 2, 3, 5, 7).
507      */
508
509     // TODO: this is just a convenient reuse because programHandle_ currently is in charge of creating context
510     pmeGpu->archSpecific->context = pmeGpu->programHandle_->impl_->context;
511
512     // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?) and reuse in NB
513     if (GMX_GPU == GMX_GPU_CUDA)
514     {
515         /* WARNING: CUDA timings are incorrect with multiple streams.
516          *          This is the main reason why they are disabled by default.
517          */
518         // TODO: Consider turning on by default when we can detect nr of streams.
519         pmeGpu->archSpecific->useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
520     }
521     else if (GMX_GPU == GMX_GPU_OPENCL)
522     {
523         pmeGpu->archSpecific->useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
524     }
525
526 #if GMX_GPU == GMX_GPU_CUDA
527     pmeGpu->maxGridWidthX = pmeGpu->deviceInfo->prop.maxGridSize[0];
528 #elif GMX_GPU == GMX_GPU_OPENCL
529     pmeGpu->maxGridWidthX = INT32_MAX / 2;
530     // TODO: is there no really global work size limit in OpenCL?
531 #endif
532
533     /* Creating a PME GPU stream:
534      * - default high priority with CUDA
535      * - no priorities implemented yet with OpenCL; see #2532
536      */
537 #if GMX_GPU == GMX_GPU_CUDA
538     cudaError_t stat;
539     int         highest_priority, lowest_priority;
540     stat = cudaDeviceGetStreamPriorityRange(&lowest_priority, &highest_priority);
541     CU_RET_ERR(stat, "PME cudaDeviceGetStreamPriorityRange failed");
542     stat = cudaStreamCreateWithPriority(&pmeGpu->archSpecific->pmeStream,
543                                         cudaStreamDefault, //cudaStreamNonBlocking,
544                                         highest_priority);
545     CU_RET_ERR(stat, "cudaStreamCreateWithPriority on the PME stream failed");
546 #elif GMX_GPU == GMX_GPU_OPENCL
547     cl_command_queue_properties queueProperties = pmeGpu->archSpecific->useTiming ? CL_QUEUE_PROFILING_ENABLE : 0;
548     cl_device_id                device_id       = pmeGpu->deviceInfo->ocl_gpu_id.ocl_device_id;
549     cl_int                      clError;
550     pmeGpu->archSpecific->pmeStream = clCreateCommandQueue(pmeGpu->archSpecific->context,
551                                                            device_id, queueProperties, &clError);
552     if (clError != CL_SUCCESS)
553     {
554         GMX_THROW(gmx::InternalError("Failed to create PME command queue"));
555     }
556 #endif
557 }
558
559 void pme_gpu_destroy_specific(const PmeGpu *pmeGpu)
560 {
561 #if GMX_GPU == GMX_GPU_CUDA
562     /* Destroy the CUDA stream */
563     cudaError_t stat = cudaStreamDestroy(pmeGpu->archSpecific->pmeStream);
564     CU_RET_ERR(stat, "PME cudaStreamDestroy error");
565 #elif GMX_GPU == GMX_GPU_OPENCL
566     cl_int clError = clReleaseCommandQueue(pmeGpu->archSpecific->pmeStream);
567     if (clError != CL_SUCCESS)
568     {
569         gmx_warning("Failed to destroy PME command queue");
570     }
571 #endif
572 }
573
574 void pme_gpu_reinit_3dfft(const PmeGpu *pmeGpu)
575 {
576     if (pme_gpu_performs_FFT(pmeGpu))
577     {
578         pmeGpu->archSpecific->fftSetup.resize(0);
579         for (int i = 0; i < pmeGpu->common->ngrids; i++)
580         {
581             pmeGpu->archSpecific->fftSetup.push_back(std::make_unique<GpuParallel3dFft>(pmeGpu));
582         }
583     }
584 }
585
586 void pme_gpu_destroy_3dfft(const PmeGpu *pmeGpu)
587 {
588     pmeGpu->archSpecific->fftSetup.resize(0);
589 }
590
591 int getSplineParamFullIndex(int order, int splineIndex, int dimIndex, int atomIndex, int atomsPerWarp)
592 {
593     if (order != c_pmeGpuOrder)
594     {
595         throw order;
596     }
597     constexpr int fixedOrder = c_pmeGpuOrder;
598     GMX_UNUSED_VALUE(fixedOrder);
599
600     const int atomWarpIndex = atomIndex % atomsPerWarp;
601     const int warpIndex     = atomIndex / atomsPerWarp;
602     int       indexBase, result;
603     switch (atomsPerWarp)
604     {
605         case 1:
606             indexBase = getSplineParamIndexBase<fixedOrder, 1>(warpIndex, atomWarpIndex);
607             result    = getSplineParamIndex<fixedOrder, 1>(indexBase, dimIndex, splineIndex);
608             break;
609
610         case 2:
611             indexBase = getSplineParamIndexBase<fixedOrder, 2>(warpIndex, atomWarpIndex);
612             result    = getSplineParamIndex<fixedOrder, 2>(indexBase, dimIndex, splineIndex);
613             break;
614
615         case 4:
616             indexBase = getSplineParamIndexBase<fixedOrder, 4>(warpIndex, atomWarpIndex);
617             result    = getSplineParamIndex<fixedOrder, 4>(indexBase, dimIndex, splineIndex);
618             break;
619
620         case 8:
621             indexBase = getSplineParamIndexBase<fixedOrder, 8>(warpIndex, atomWarpIndex);
622             result    = getSplineParamIndex<fixedOrder, 8>(indexBase, dimIndex, splineIndex);
623             break;
624
625         default:
626             GMX_THROW(gmx::NotImplementedError(gmx::formatString("Test function call not unrolled for atomsPerWarp = %d in getSplineParamFullIndex", atomsPerWarp)));
627     }
628     return result;
629 }
630
631 PmeOutput pme_gpu_getEnergyAndVirial(const gmx_pme_t &pme)
632 {
633     PmeOutput     output;
634
635     const PmeGpu *pmeGpu = pme.gpu;
636     for (int j = 0; j < c_virialAndEnergyCount; j++)
637     {
638         GMX_ASSERT(std::isfinite(pmeGpu->staging.h_virialAndEnergy[j]), "PME GPU produces incorrect energy/virial.");
639     }
640
641     unsigned int j = 0;
642     output.coulombVirial_[XX][XX] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
643     output.coulombVirial_[YY][YY] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
644     output.coulombVirial_[ZZ][ZZ] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
645     output.coulombVirial_[XX][YY] = output.coulombVirial_[YY][XX] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
646     output.coulombVirial_[XX][ZZ] = output.coulombVirial_[ZZ][XX] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
647     output.coulombVirial_[YY][ZZ] = output.coulombVirial_[ZZ][YY] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
648     output.coulombEnergy_         = 0.5f * pmeGpu->staging.h_virialAndEnergy[j++];
649
650     return output;
651 }
652
653 PmeOutput pme_gpu_getOutput(const gmx_pme_t &pme,
654                             const int        flags)
655 {
656     PmeGpu    *pmeGpu = pme.gpu;
657     const bool haveComputedEnergyAndVirial = (flags & GMX_PME_CALC_ENER_VIR) != 0;
658     if (!haveComputedEnergyAndVirial)
659     {
660         // The caller knows from the flags that the energy and the virial are not usable
661         PmeOutput output;
662         output.forces_ = pmeGpu->staging.h_forces;
663         return output;
664     }
665
666     if (pme_gpu_performs_solve(pmeGpu))
667     {
668         PmeOutput output = pme_gpu_getEnergyAndVirial(pme);
669         output.forces_ = pmeGpu->staging.h_forces;
670         return output;
671     }
672     else
673     {
674         PmeOutput output;
675         get_pme_ener_vir_q(pme.solve_work, pme.nthread, &output);
676         output.forces_ = pmeGpu->staging.h_forces;
677         return output;
678     }
679
680 }
681
682 void pme_gpu_update_input_box(PmeGpu gmx_unused       *pmeGpu,
683                               const matrix gmx_unused  box)
684 {
685 #if GMX_DOUBLE
686     GMX_THROW(gmx::NotImplementedError("PME is implemented for single-precision only on GPU"));
687 #else
688     matrix  scaledBox;
689     pmeGpu->common->boxScaler->scaleBox(box, scaledBox);
690     auto   *kernelParamsPtr      = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
691     kernelParamsPtr->current.boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
692     GMX_ASSERT(kernelParamsPtr->current.boxVolume != 0.0f, "Zero volume of the unit cell");
693     matrix recipBox;
694     gmx::invertBoxMatrix(scaledBox, recipBox);
695
696     /* The GPU recipBox is transposed as compared to the CPU recipBox.
697      * Spread uses matrix columns (while solve and gather use rows).
698      * There is no particular reason for this; it might be further rethought/optimized for better access patterns.
699      */
700     const real newRecipBox[DIM][DIM] =
701     {
702         {recipBox[XX][XX], recipBox[YY][XX], recipBox[ZZ][XX]},
703         {             0.0, recipBox[YY][YY], recipBox[ZZ][YY]},
704         {             0.0,              0.0, recipBox[ZZ][ZZ]}
705     };
706     memcpy(kernelParamsPtr->current.recipBox, newRecipBox, sizeof(matrix));
707 #endif
708 }
709
710 /*! \brief \libinternal
711  * (Re-)initializes all the PME GPU data related to the grid size and cut-off.
712  *
713  * \param[in] pmeGpu            The PME GPU structure.
714  */
715 static void pme_gpu_reinit_grids(PmeGpu *pmeGpu)
716 {
717     auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
718     kernelParamsPtr->grid.ewaldFactor = (M_PI * M_PI) / (pmeGpu->common->ewaldcoeff_q * pmeGpu->common->ewaldcoeff_q);
719
720     /* The grid size variants */
721     for (int i = 0; i < DIM; i++)
722     {
723         kernelParamsPtr->grid.realGridSize[i]       = pmeGpu->common->nk[i];
724         kernelParamsPtr->grid.realGridSizeFP[i]     = static_cast<float>(kernelParamsPtr->grid.realGridSize[i]);
725         kernelParamsPtr->grid.realGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
726
727         // The complex grid currently uses no padding;
728         // if it starts to do so, then another test should be added for that
729         kernelParamsPtr->grid.complexGridSize[i]       = kernelParamsPtr->grid.realGridSize[i];
730         kernelParamsPtr->grid.complexGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
731     }
732     /* FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
733     if (!pme_gpu_performs_FFT(pmeGpu))
734     {
735         // This allows for GPU spreading grid and CPU fftgrid to have the same layout, so that we can copy the data directly
736         kernelParamsPtr->grid.realGridSizePadded[ZZ] = (kernelParamsPtr->grid.realGridSize[ZZ] / 2 + 1) * 2;
737     }
738
739     /* GPU FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
740     kernelParamsPtr->grid.complexGridSize[ZZ] /= 2;
741     kernelParamsPtr->grid.complexGridSize[ZZ]++;
742     kernelParamsPtr->grid.complexGridSizePadded[ZZ] = kernelParamsPtr->grid.complexGridSize[ZZ];
743
744     pme_gpu_realloc_and_copy_fract_shifts(pmeGpu);
745     pme_gpu_realloc_and_copy_bspline_values(pmeGpu);
746     pme_gpu_realloc_grids(pmeGpu);
747     pme_gpu_reinit_3dfft(pmeGpu);
748 }
749
750 /* Several GPU functions that refer to the CPU PME data live here.
751  * We would like to keep these away from the GPU-framework specific code for clarity,
752  * as well as compilation issues with MPI.
753  */
754
755 /*! \brief \libinternal
756  * Copies everything useful from the PME CPU to the PME GPU structure.
757  * The goal is to minimize interaction with the PME CPU structure in the GPU code.
758  *
759  * \param[in] pme         The PME structure.
760  */
761 static void pme_gpu_copy_common_data_from(const gmx_pme_t *pme)
762 {
763     /* TODO: Consider refactoring the CPU PME code to use the same structure,
764      * so that this function becomes 2 lines */
765     PmeGpu *pmeGpu             = pme->gpu;
766     pmeGpu->common->ngrids        = pme->ngrids;
767     pmeGpu->common->epsilon_r     = pme->epsilon_r;
768     pmeGpu->common->ewaldcoeff_q  = pme->ewaldcoeff_q;
769     pmeGpu->common->nk[XX]        = pme->nkx;
770     pmeGpu->common->nk[YY]        = pme->nky;
771     pmeGpu->common->nk[ZZ]        = pme->nkz;
772     pmeGpu->common->pme_order     = pme->pme_order;
773     if (pmeGpu->common->pme_order != c_pmeGpuOrder)
774     {
775         GMX_THROW(gmx::NotImplementedError("pme_order != 4 is not implemented!"));
776     }
777     for (int i = 0; i < DIM; i++)
778     {
779         pmeGpu->common->bsp_mod[i].assign(pme->bsp_mod[i], pme->bsp_mod[i] + pmeGpu->common->nk[i]);
780     }
781     const int cellCount = c_pmeNeighborUnitcellCount;
782     pmeGpu->common->fsh.resize(0);
783     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshx, pme->fshx + cellCount * pme->nkx);
784     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshy, pme->fshy + cellCount * pme->nky);
785     pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshz, pme->fshz + cellCount * pme->nkz);
786     pmeGpu->common->nn.resize(0);
787     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnx, pme->nnx + cellCount * pme->nkx);
788     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nny, pme->nny + cellCount * pme->nky);
789     pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnz, pme->nnz + cellCount * pme->nkz);
790     pmeGpu->common->runMode   = pme->runMode;
791     pmeGpu->common->boxScaler = pme->boxScaler;
792 }
793
794 /*! \libinternal \brief
795  * Initializes the PME GPU data at the beginning of the run.
796  * TODO: this should become PmeGpu::PmeGpu()
797  *
798  * \param[in,out] pme            The PME structure.
799  * \param[in,out] gpuInfo        The GPU information structure.
800  * \param[in]     pmeGpuProgram  The handle to the program/kernel data created outside (e.g. in unit tests/runner)
801  */
802 static void pme_gpu_init(gmx_pme_t               *pme,
803                          const gmx_device_info_t *gpuInfo,
804                          PmeGpuProgramHandle      pmeGpuProgram)
805 {
806     pme->gpu          = new PmeGpu();
807     PmeGpu *pmeGpu = pme->gpu;
808     changePinningPolicy(&pmeGpu->staging.h_forces, pme_get_pinning_policy());
809     pmeGpu->common = std::make_shared<PmeShared>();
810
811     /* These settings are set here for the whole run; dynamic ones are set in pme_gpu_reinit() */
812     /* A convenience variable. */
813     pmeGpu->settings.useDecomposition = (pme->nnodes == 1);
814     /* TODO: CPU gather with GPU spread is broken due to different theta/dtheta layout. */
815     pmeGpu->settings.performGPUGather = true;
816
817     pme_gpu_set_testing(pmeGpu, false);
818
819     pmeGpu->deviceInfo = gpuInfo;
820     GMX_ASSERT(pmeGpuProgram != nullptr, "GPU kernels must be already compiled");
821     pmeGpu->programHandle_ = pmeGpuProgram;
822
823     pme_gpu_init_internal(pmeGpu);
824     pme_gpu_alloc_energy_virial(pmeGpu);
825
826     pme_gpu_copy_common_data_from(pme);
827
828     GMX_ASSERT(pmeGpu->common->epsilon_r != 0.0f, "PME GPU: bad electrostatic coefficient");
829
830     auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
831     kernelParamsPtr->constants.elFactor = ONE_4PI_EPS0 / pmeGpu->common->epsilon_r;
832 }
833
834 void pme_gpu_transform_spline_atom_data(const PmeGpu *pmeGpu, const pme_atomcomm_t *atc,
835                                         PmeSplineDataType type, int dimIndex, PmeLayoutTransform transform)
836 {
837     // The GPU atom spline data is laid out in a different way currently than the CPU one.
838     // This function converts the data from GPU to CPU layout (in the host memory).
839     // It is only intended for testing purposes so far.
840     // Ideally we should use similar layouts on CPU and GPU if we care about mixed modes and their performance
841     // (e.g. spreading on GPU, gathering on CPU).
842     GMX_RELEASE_ASSERT(atc->nthread == 1, "Only the serial PME data layout is supported");
843     const uintmax_t threadIndex  = 0;
844     const auto      atomCount    = pme_gpu_get_kernel_params_base_ptr(pmeGpu)->atoms.nAtoms;
845     const auto      atomsPerWarp = pme_gpu_get_atoms_per_warp(pmeGpu);
846     const auto      pmeOrder     = pmeGpu->common->pme_order;
847     GMX_ASSERT(pmeOrder == c_pmeGpuOrder, "Only PME order 4 is implemented");
848
849     real           *cpuSplineBuffer;
850     float          *h_splineBuffer;
851     switch (type)
852     {
853         case PmeSplineDataType::Values:
854             cpuSplineBuffer = atc->spline[threadIndex].theta[dimIndex];
855             h_splineBuffer  = pmeGpu->staging.h_theta;
856             break;
857
858         case PmeSplineDataType::Derivatives:
859             cpuSplineBuffer = atc->spline[threadIndex].dtheta[dimIndex];
860             h_splineBuffer  = pmeGpu->staging.h_dtheta;
861             break;
862
863         default:
864             GMX_THROW(gmx::InternalError("Unknown spline data type"));
865     }
866
867     for (auto atomIndex = 0; atomIndex < atomCount; atomIndex++)
868     {
869         for (auto orderIndex = 0; orderIndex < pmeOrder; orderIndex++)
870         {
871             const auto gpuValueIndex = getSplineParamFullIndex(pmeOrder, orderIndex, dimIndex, atomIndex, atomsPerWarp);
872             const auto cpuValueIndex = atomIndex * pmeOrder + orderIndex;
873             GMX_ASSERT(cpuValueIndex < atomCount * pmeOrder, "Atom spline data index out of bounds (while transforming GPU data layout for host)");
874             switch (transform)
875             {
876                 case PmeLayoutTransform::GpuToHost:
877                     cpuSplineBuffer[cpuValueIndex] = h_splineBuffer[gpuValueIndex];
878                     break;
879
880                 case PmeLayoutTransform::HostToGpu:
881                     h_splineBuffer[gpuValueIndex] = cpuSplineBuffer[cpuValueIndex];
882                     break;
883
884                 default:
885                     GMX_THROW(gmx::InternalError("Unknown layout transform"));
886             }
887         }
888     }
889 }
890
891 void pme_gpu_get_real_grid_sizes(const PmeGpu *pmeGpu, gmx::IVec *gridSize, gmx::IVec *paddedGridSize)
892 {
893     GMX_ASSERT(gridSize != nullptr, "");
894     GMX_ASSERT(paddedGridSize != nullptr, "");
895     GMX_ASSERT(pmeGpu != nullptr, "");
896     auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
897     for (int i = 0; i < DIM; i++)
898     {
899         (*gridSize)[i]       = kernelParamsPtr->grid.realGridSize[i];
900         (*paddedGridSize)[i] = kernelParamsPtr->grid.realGridSizePadded[i];
901     }
902 }
903
904 void pme_gpu_reinit(gmx_pme_t               *pme,
905                     const gmx_device_info_t *gpuInfo,
906                     PmeGpuProgramHandle      pmeGpuProgram)
907 {
908     if (!pme_gpu_active(pme))
909     {
910         return;
911     }
912
913     if (!pme->gpu)
914     {
915         /* First-time initialization */
916         pme_gpu_init(pme, gpuInfo, pmeGpuProgram);
917     }
918     else
919     {
920         /* After this call nothing in the GPU code should refer to the gmx_pme_t *pme itself - until the next pme_gpu_reinit */
921         pme_gpu_copy_common_data_from(pme);
922     }
923     /* GPU FFT will only get used for a single rank.*/
924     pme->gpu->settings.performGPUFFT   = (pme->gpu->common->runMode == PmeRunMode::GPU) && !pme_gpu_uses_dd(pme->gpu);
925     pme->gpu->settings.performGPUSolve = (pme->gpu->common->runMode == PmeRunMode::GPU);
926
927     /* Reinit active timers */
928     pme_gpu_reinit_timings(pme->gpu);
929
930     pme_gpu_reinit_grids(pme->gpu);
931     // Note: if timing the reinit launch overhead becomes more relevant
932     // (e.g. with regulat PP-PME re-balancing), we should pass wcycle here.
933     pme_gpu_reinit_computation(pme, nullptr);
934     /* Clear the previous box - doesn't hurt, and forces the PME CPU recipbox
935      * update for mixed mode on grid switch. TODO: use shared recipbox field.
936      */
937     std::memset(pme->gpu->common->previousBox, 0, sizeof(pme->gpu->common->previousBox));
938 }
939
940 void pme_gpu_destroy(PmeGpu *pmeGpu)
941 {
942     /* Free lots of data */
943     pme_gpu_free_energy_virial(pmeGpu);
944     pme_gpu_free_bspline_values(pmeGpu);
945     pme_gpu_free_forces(pmeGpu);
946     pme_gpu_free_coordinates(pmeGpu);
947     pme_gpu_free_coefficients(pmeGpu);
948     pme_gpu_free_spline_data(pmeGpu);
949     pme_gpu_free_grid_indices(pmeGpu);
950     pme_gpu_free_fract_shifts(pmeGpu);
951     pme_gpu_free_grids(pmeGpu);
952
953     pme_gpu_destroy_3dfft(pmeGpu);
954
955     /* Free the GPU-framework specific data last */
956     pme_gpu_destroy_specific(pmeGpu);
957
958     delete pmeGpu;
959 }
960
961 void pme_gpu_reinit_atoms(PmeGpu *pmeGpu, const int nAtoms, const real *charges)
962 {
963     auto      *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
964     kernelParamsPtr->atoms.nAtoms = nAtoms;
965     const int  alignment = pme_gpu_get_atom_data_alignment(pmeGpu);
966     pmeGpu->nAtomsPadded = ((nAtoms + alignment - 1) / alignment) * alignment;
967     const int  nAtomsAlloc   = c_usePadding ? pmeGpu->nAtomsPadded : nAtoms;
968     const bool haveToRealloc = (pmeGpu->nAtomsAlloc < nAtomsAlloc); /* This check might be redundant, but is logical */
969     pmeGpu->nAtomsAlloc = nAtomsAlloc;
970
971 #if GMX_DOUBLE
972     GMX_RELEASE_ASSERT(false, "Only single precision supported");
973     GMX_UNUSED_VALUE(charges);
974 #else
975     pme_gpu_realloc_and_copy_input_coefficients(pmeGpu, reinterpret_cast<const float *>(charges));
976     /* Could also be checked for haveToRealloc, but the copy always needs to be performed */
977 #endif
978
979     if (haveToRealloc)
980     {
981         pme_gpu_realloc_coordinates(pmeGpu);
982         pme_gpu_realloc_forces(pmeGpu);
983         pme_gpu_realloc_spline_data(pmeGpu);
984         pme_gpu_realloc_grid_indices(pmeGpu);
985     }
986 }
987
988 void pme_gpu_3dfft(const PmeGpu *pmeGpu, gmx_fft_direction dir, int grid_index)
989 {
990     int   timerId = (dir == GMX_FFT_REAL_TO_COMPLEX) ? gtPME_FFT_R2C : gtPME_FFT_C2R;
991
992     pme_gpu_start_timing(pmeGpu, timerId);
993     pmeGpu->archSpecific->fftSetup[grid_index]->perform3dFft(dir, pme_gpu_fetch_timing_event(pmeGpu, timerId));
994     pme_gpu_stop_timing(pmeGpu, timerId);
995 }
996
997 /*! \brief
998  * Given possibly large \p blockCount, returns a compact 1D or 2D grid for kernel scheduling,
999  * to minimize number of unused blocks.
1000  */
1001 std::pair<int, int> inline pmeGpuCreateGrid(const PmeGpu *pmeGpu, int blockCount)
1002 {
1003     // How many maximum widths in X do we need (hopefully just one)
1004     const int minRowCount = (blockCount + pmeGpu->maxGridWidthX - 1) / pmeGpu->maxGridWidthX;
1005     // Trying to make things even
1006     const int colCount = (blockCount + minRowCount - 1) / minRowCount;
1007     GMX_ASSERT((colCount * minRowCount - blockCount) >= 0, "pmeGpuCreateGrid: totally wrong");
1008     GMX_ASSERT((colCount * minRowCount - blockCount) < minRowCount, "pmeGpuCreateGrid: excessive blocks");
1009     return std::pair<int, int>(colCount, minRowCount);
1010 }
1011
1012 void pme_gpu_spread(const PmeGpu    *pmeGpu,
1013                     int gmx_unused   gridIndex,
1014                     real            *h_grid,
1015                     bool             computeSplines,
1016                     bool             spreadCharges)
1017 {
1018     GMX_ASSERT(computeSplines || spreadCharges, "PME spline/spread kernel has invalid input (nothing to do)");
1019     const auto   *kernelParamsPtr = pmeGpu->kernelParams.get();
1020     GMX_ASSERT(kernelParamsPtr->atoms.nAtoms > 0, "No atom data in PME GPU spread");
1021
1022     const size_t blockSize  = pmeGpu->programHandle_->impl_->spreadWorkGroupSize;
1023
1024     const int    order         = pmeGpu->common->pme_order;
1025     GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1026     const int    atomsPerBlock = blockSize / c_pmeSpreadGatherThreadsPerAtom;
1027     // TODO: pick smaller block size in runtime if needed
1028     // (e.g. on 660 Ti where 50% occupancy is ~25% faster than 100% occupancy with RNAse (~17.8k atoms))
1029     // If doing so, change atomsPerBlock in the kernels as well.
1030     // TODO: test varying block sizes on modern arch-s as well
1031     // TODO: also consider using cudaFuncSetCacheConfig() for preferring shared memory on older architectures
1032     //(for spline data mostly, together with varying PME_GPU_PARALLEL_SPLINE define)
1033     GMX_ASSERT(!c_usePadding || !(c_pmeAtomDataAlignment % atomsPerBlock), "inconsistent atom data padding vs. spreading block size");
1034
1035     const int          blockCount = pmeGpu->nAtomsPadded / atomsPerBlock;
1036     auto               dimGrid    = pmeGpuCreateGrid(pmeGpu, blockCount);
1037
1038     KernelLaunchConfig config;
1039     config.blockSize[0] = config.blockSize[1] = order;
1040     config.blockSize[2] = atomsPerBlock;
1041     config.gridSize[0]  = dimGrid.first;
1042     config.gridSize[1]  = dimGrid.second;
1043     config.stream       = pmeGpu->archSpecific->pmeStream;
1044
1045     int  timingId;
1046     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1047     if (computeSplines)
1048     {
1049         if (spreadCharges)
1050         {
1051             timingId  = gtPME_SPLINEANDSPREAD;
1052             kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernel;
1053         }
1054         else
1055         {
1056             timingId  = gtPME_SPLINE;
1057             kernelPtr = pmeGpu->programHandle_->impl_->splineKernel;
1058         }
1059     }
1060     else
1061     {
1062         timingId  = gtPME_SPREAD;
1063         kernelPtr = pmeGpu->programHandle_->impl_->spreadKernel;
1064     }
1065
1066     pme_gpu_start_timing(pmeGpu, timingId);
1067     auto      *timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1068 #if c_canEmbedBuffers
1069     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1070 #else
1071     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr,
1072                                                       &kernelParamsPtr->atoms.d_theta,
1073                                                       &kernelParamsPtr->atoms.d_dtheta,
1074                                                       &kernelParamsPtr->atoms.d_gridlineIndices,
1075                                                       &kernelParamsPtr->grid.d_realGrid,
1076                                                       &kernelParamsPtr->grid.d_fractShiftsTable,
1077                                                       &kernelParamsPtr->grid.d_gridlineIndicesTable,
1078                                                       &kernelParamsPtr->atoms.d_coefficients,
1079                                                       &kernelParamsPtr->atoms.d_coordinates
1080                                                       );
1081 #endif
1082
1083     launchGpuKernel(kernelPtr, config, timingEvent, "PME spline/spread", kernelArgs);
1084     pme_gpu_stop_timing(pmeGpu, timingId);
1085
1086     const bool copyBackGrid = spreadCharges && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_FFT(pmeGpu));
1087     if (copyBackGrid)
1088     {
1089         pme_gpu_copy_output_spread_grid(pmeGpu, h_grid);
1090     }
1091     const bool copyBackAtomData = computeSplines && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_gather(pmeGpu));
1092     if (copyBackAtomData)
1093     {
1094         pme_gpu_copy_output_spread_atom_data(pmeGpu);
1095     }
1096 }
1097
1098 void pme_gpu_solve(const PmeGpu *pmeGpu, t_complex *h_grid,
1099                    GridOrdering gridOrdering, bool computeEnergyAndVirial)
1100 {
1101     const bool   copyInputAndOutputGrid = pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_FFT(pmeGpu);
1102
1103     auto        *kernelParamsPtr = pmeGpu->kernelParams.get();
1104
1105     float       *h_gridFloat = reinterpret_cast<float *>(h_grid);
1106     if (copyInputAndOutputGrid)
1107     {
1108         copyToDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid, h_gridFloat,
1109                            0, pmeGpu->archSpecific->complexGridSize,
1110                            pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
1111     }
1112
1113     int majorDim = -1, middleDim = -1, minorDim = -1;
1114     switch (gridOrdering)
1115     {
1116         case GridOrdering::YZX:
1117             majorDim  = YY;
1118             middleDim = ZZ;
1119             minorDim  = XX;
1120             break;
1121
1122         case GridOrdering::XYZ:
1123             majorDim  = XX;
1124             middleDim = YY;
1125             minorDim  = ZZ;
1126             break;
1127
1128         default:
1129             GMX_ASSERT(false, "Implement grid ordering here and below for the kernel launch");
1130     }
1131
1132     const int    maxBlockSize       = pmeGpu->programHandle_->impl_->solveMaxWorkGroupSize;
1133
1134     const int    gridLineSize       = pmeGpu->kernelParams->grid.complexGridSize[minorDim];
1135     const int    gridLinesPerBlock  = std::max(maxBlockSize / gridLineSize, 1);
1136     const int    blocksPerGridLine  = (gridLineSize + maxBlockSize - 1) / maxBlockSize;
1137     int          cellsPerBlock;
1138     if (blocksPerGridLine == 1)
1139     {
1140         cellsPerBlock               = gridLineSize * gridLinesPerBlock;
1141     }
1142     else
1143     {
1144         cellsPerBlock               = (gridLineSize + blocksPerGridLine - 1) / blocksPerGridLine;
1145     }
1146     const int    warpSize           = pmeGpu->programHandle_->impl_->warpSize;
1147     const int    blockSize          = (cellsPerBlock + warpSize - 1) / warpSize * warpSize;
1148
1149     static_assert(GMX_GPU != GMX_GPU_CUDA || c_solveMaxWarpsPerBlock/2 >= 4,
1150                   "The CUDA solve energy kernels needs at least 4 warps. "
1151                   "Here we launch at least half of the max warps.");
1152
1153     KernelLaunchConfig config;
1154     config.blockSize[0] = blockSize;
1155     config.gridSize[0]  = blocksPerGridLine;
1156     // rounding up to full warps so that shuffle operations produce defined results
1157     config.gridSize[1]  = (pmeGpu->kernelParams->grid.complexGridSize[middleDim] + gridLinesPerBlock - 1) / gridLinesPerBlock;
1158     config.gridSize[2]  = pmeGpu->kernelParams->grid.complexGridSize[majorDim];
1159     config.stream       = pmeGpu->archSpecific->pmeStream;
1160
1161     int  timingId = gtPME_SOLVE;
1162     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1163     if (gridOrdering == GridOrdering::YZX)
1164     {
1165         kernelPtr = computeEnergyAndVirial ?
1166             pmeGpu->programHandle_->impl_->solveYZXEnergyKernel :
1167             pmeGpu->programHandle_->impl_->solveYZXKernel;
1168     }
1169     else if (gridOrdering == GridOrdering::XYZ)
1170     {
1171         kernelPtr = computeEnergyAndVirial ?
1172             pmeGpu->programHandle_->impl_->solveXYZEnergyKernel :
1173             pmeGpu->programHandle_->impl_->solveXYZKernel;
1174     }
1175
1176     pme_gpu_start_timing(pmeGpu, timingId);
1177     auto      *timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1178 #if c_canEmbedBuffers
1179     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1180 #else
1181     const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr,
1182                                                       &kernelParamsPtr->grid.d_splineModuli,
1183                                                       &kernelParamsPtr->constants.d_virialAndEnergy,
1184                                                       &kernelParamsPtr->grid.d_fourierGrid);
1185 #endif
1186     launchGpuKernel(kernelPtr, config, timingEvent, "PME solve", kernelArgs);
1187     pme_gpu_stop_timing(pmeGpu, timingId);
1188
1189     if (computeEnergyAndVirial)
1190     {
1191         copyFromDeviceBuffer(pmeGpu->staging.h_virialAndEnergy, &kernelParamsPtr->constants.d_virialAndEnergy,
1192                              0, c_virialAndEnergyCount,
1193                              pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
1194     }
1195
1196     if (copyInputAndOutputGrid)
1197     {
1198         copyFromDeviceBuffer(h_gridFloat, &kernelParamsPtr->grid.d_fourierGrid,
1199                              0, pmeGpu->archSpecific->complexGridSize,
1200                              pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
1201     }
1202 }
1203
1204 void pme_gpu_gather(PmeGpu                *pmeGpu,
1205                     PmeForceOutputHandling forceTreatment,
1206                     const float           *h_grid
1207                     )
1208 {
1209     /* Copying the input CPU forces for reduction */
1210     if (forceTreatment != PmeForceOutputHandling::Set)
1211     {
1212         pme_gpu_copy_input_forces(pmeGpu);
1213     }
1214
1215     if (!pme_gpu_performs_FFT(pmeGpu) || pme_gpu_is_testing(pmeGpu))
1216     {
1217         pme_gpu_copy_input_gather_grid(pmeGpu, const_cast<float *>(h_grid));
1218     }
1219
1220     if (pme_gpu_is_testing(pmeGpu))
1221     {
1222         pme_gpu_copy_input_gather_atom_data(pmeGpu);
1223     }
1224
1225     const size_t blockSize     = pmeGpu->programHandle_->impl_->gatherWorkGroupSize;
1226     const int    atomsPerBlock = blockSize / c_pmeSpreadGatherThreadsPerAtom;
1227     GMX_ASSERT(!c_usePadding || !(c_pmeAtomDataAlignment % atomsPerBlock), "inconsistent atom data padding vs. gathering block size");
1228
1229     const int          blockCount = pmeGpu->nAtomsPadded / atomsPerBlock;
1230     auto               dimGrid    = pmeGpuCreateGrid(pmeGpu, blockCount);
1231
1232
1233     const int order = pmeGpu->common->pme_order;
1234     GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1235
1236     KernelLaunchConfig config;
1237     config.blockSize[0] = config.blockSize[1] = order;
1238     config.blockSize[2] = atomsPerBlock;
1239     config.gridSize[0]  = dimGrid.first;
1240     config.gridSize[1]  = dimGrid.second;
1241     config.stream       = pmeGpu->archSpecific->pmeStream;
1242
1243     // TODO test different cache configs
1244
1245     int  timingId = gtPME_GATHER;
1246     // TODO design kernel selection getters and make PmeGpu a friend of PmeGpuProgramImpl
1247     PmeGpuProgramImpl::PmeKernelHandle kernelPtr = (forceTreatment == PmeForceOutputHandling::Set) ?
1248         pmeGpu->programHandle_->impl_->gatherKernel :
1249         pmeGpu->programHandle_->impl_->gatherReduceWithInputKernel;
1250
1251     pme_gpu_start_timing(pmeGpu, timingId);
1252     auto       *timingEvent     = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1253     const auto *kernelParamsPtr = pmeGpu->kernelParams.get();
1254 #if c_canEmbedBuffers
1255     const auto  kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1256 #else
1257     const auto  kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr,
1258                                                        &kernelParamsPtr->atoms.d_coefficients,
1259                                                        &kernelParamsPtr->grid.d_realGrid,
1260                                                        &kernelParamsPtr->atoms.d_theta,
1261                                                        &kernelParamsPtr->atoms.d_dtheta,
1262                                                        &kernelParamsPtr->atoms.d_gridlineIndices,
1263                                                        &kernelParamsPtr->atoms.d_forces);
1264 #endif
1265     launchGpuKernel(kernelPtr, config, timingEvent, "PME gather", kernelArgs);
1266     pme_gpu_stop_timing(pmeGpu, timingId);
1267
1268     pme_gpu_copy_output_forces(pmeGpu);
1269 }