Remove texture reference support in the CUDA
[alexxy/gromacs.git] / src / gromacs / ewald / pme.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 /*! \internal \file
37  * \brief This file contains internal CUDA function implementations
38  * for performing the PME calculations on GPU.
39  *
40  * \author Aleksei Iupinov <a.yupinov@gmail.com>
41  */
42
43 #include "gmxpre.h"
44
45 #include <cmath>
46
47 /* The rest */
48 #include "pme.h"
49
50 #include "gromacs/gpu_utils/cudautils.cuh"
51 #include "gromacs/gpu_utils/pmalloc_cuda.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/smalloc.h"
54
55 #include "pme.cuh"
56 #include "pme-3dfft.cuh"
57 #include "pme-grid.h"
58
59 int pme_gpu_get_atom_data_alignment(const PmeGpu *pmeGpu)
60 {
61     const int order = pmeGpu->common->pme_order;
62     GMX_ASSERT(order > 0, "Invalid PME order");
63     return PME_ATOM_DATA_ALIGNMENT;
64 }
65
66 int pme_gpu_get_atoms_per_warp(const PmeGpu *pmeGpu)
67 {
68     const int order = pmeGpu->common->pme_order;
69     GMX_ASSERT(order > 0, "Invalid PME order");
70     return PME_SPREADGATHER_ATOMS_PER_WARP;
71 }
72
73 void pme_gpu_synchronize(const PmeGpu *pmeGpu)
74 {
75     cudaError_t stat = cudaStreamSynchronize(pmeGpu->archSpecific->pmeStream);
76     CU_RET_ERR(stat, "Failed to synchronize the PME GPU stream!");
77 }
78
79 void pme_gpu_alloc_energy_virial(const PmeGpu *pmeGpu)
80 {
81     const size_t energyAndVirialSize = c_virialAndEnergyCount * sizeof(float);
82     cudaError_t  stat                = cudaMalloc((void **)&pmeGpu->kernelParams->constants.d_virialAndEnergy, energyAndVirialSize);
83     CU_RET_ERR(stat, "cudaMalloc failed on PME energy and virial");
84     pmalloc((void **)&pmeGpu->staging.h_virialAndEnergy, energyAndVirialSize);
85 }
86
87 void pme_gpu_free_energy_virial(PmeGpu *pmeGpu)
88 {
89     cudaError_t stat = cudaFree(pmeGpu->kernelParams->constants.d_virialAndEnergy);
90     CU_RET_ERR(stat, "cudaFree failed on PME energy and virial");
91     pmeGpu->kernelParams->constants.d_virialAndEnergy = nullptr;
92     pfree(pmeGpu->staging.h_virialAndEnergy);
93     pmeGpu->staging.h_virialAndEnergy = nullptr;
94 }
95
96 void pme_gpu_clear_energy_virial(const PmeGpu *pmeGpu)
97 {
98     cudaError_t stat = cudaMemsetAsync(pmeGpu->kernelParams->constants.d_virialAndEnergy, 0,
99                                        c_virialAndEnergyCount * sizeof(float), pmeGpu->archSpecific->pmeStream);
100     CU_RET_ERR(stat, "PME energy/virial cudaMemsetAsync error");
101 }
102
103 void pme_gpu_realloc_and_copy_bspline_values(const PmeGpu *pmeGpu)
104 {
105     const int splineValuesOffset[DIM] = {
106         0,
107         pmeGpu->kernelParams->grid.realGridSize[XX],
108         pmeGpu->kernelParams->grid.realGridSize[XX] + pmeGpu->kernelParams->grid.realGridSize[YY]
109     };
110     memcpy((void *)&pmeGpu->kernelParams->grid.splineValuesOffset, &splineValuesOffset, sizeof(splineValuesOffset));
111
112     const int newSplineValuesSize = pmeGpu->kernelParams->grid.realGridSize[XX] +
113         pmeGpu->kernelParams->grid.realGridSize[YY] +
114         pmeGpu->kernelParams->grid.realGridSize[ZZ];
115     const bool shouldRealloc = (newSplineValuesSize > pmeGpu->archSpecific->splineValuesSize);
116     cu_realloc_buffered((void **)&pmeGpu->kernelParams->grid.d_splineModuli, nullptr, sizeof(float),
117                         &pmeGpu->archSpecific->splineValuesSize, &pmeGpu->archSpecific->splineValuesSizeAlloc, newSplineValuesSize, pmeGpu->archSpecific->pmeStream, true);
118     if (shouldRealloc)
119     {
120         /* Reallocate the host buffer */
121         pfree(pmeGpu->staging.h_splineModuli);
122         pmalloc((void **)&pmeGpu->staging.h_splineModuli, newSplineValuesSize * sizeof(float));
123     }
124     for (int i = 0; i < DIM; i++)
125     {
126         memcpy(pmeGpu->staging.h_splineModuli + splineValuesOffset[i], pmeGpu->common->bsp_mod[i].data(), pmeGpu->common->bsp_mod[i].size() * sizeof(float));
127     }
128     /* TODO: pin original buffer instead! */
129     cu_copy_H2D(pmeGpu->kernelParams->grid.d_splineModuli, pmeGpu->staging.h_splineModuli,
130                 newSplineValuesSize * sizeof(float), pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
131 }
132
133 void pme_gpu_free_bspline_values(const PmeGpu *pmeGpu)
134 {
135     pfree(pmeGpu->staging.h_splineModuli);
136     cu_free_buffered(pmeGpu->kernelParams->grid.d_splineModuli, &pmeGpu->archSpecific->splineValuesSize,
137                      &pmeGpu->archSpecific->splineValuesSizeAlloc);
138 }
139
140 void pme_gpu_realloc_forces(PmeGpu *pmeGpu)
141 {
142     const size_t newForcesSize = pmeGpu->nAtomsAlloc * DIM;
143     GMX_ASSERT(newForcesSize > 0, "Bad number of atoms in PME GPU");
144     cu_realloc_buffered((void **)&pmeGpu->kernelParams->atoms.d_forces, nullptr, sizeof(float),
145                         &pmeGpu->archSpecific->forcesSize, &pmeGpu->archSpecific->forcesSizeAlloc, newForcesSize, pmeGpu->archSpecific->pmeStream, true);
146     pmeGpu->staging.h_forces.reserve(pmeGpu->nAtomsAlloc);
147     pmeGpu->staging.h_forces.resize(pmeGpu->kernelParams->atoms.nAtoms);
148 }
149
150 void pme_gpu_free_forces(const PmeGpu *pmeGpu)
151 {
152     cu_free_buffered(pmeGpu->kernelParams->atoms.d_forces, &pmeGpu->archSpecific->forcesSize, &pmeGpu->archSpecific->forcesSizeAlloc);
153 }
154
155 void pme_gpu_copy_input_forces(PmeGpu *pmeGpu)
156 {
157     const size_t forcesSize = DIM * pmeGpu->kernelParams->atoms.nAtoms * sizeof(float);
158     GMX_ASSERT(forcesSize > 0, "Bad number of atoms in PME GPU");
159     cu_copy_H2D(pmeGpu->kernelParams->atoms.d_forces, pmeGpu->staging.h_forces.data(), forcesSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
160 }
161
162 void pme_gpu_copy_output_forces(PmeGpu *pmeGpu)
163 {
164     const size_t forcesSize   = DIM * pmeGpu->kernelParams->atoms.nAtoms * sizeof(float);
165     GMX_ASSERT(forcesSize > 0, "Bad number of atoms in PME GPU");
166     cu_copy_D2H(pmeGpu->staging.h_forces.data(), pmeGpu->kernelParams->atoms.d_forces, forcesSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
167 }
168
169 void pme_gpu_realloc_coordinates(const PmeGpu *pmeGpu)
170 {
171     const size_t newCoordinatesSize = pmeGpu->nAtomsAlloc * DIM;
172     GMX_ASSERT(newCoordinatesSize > 0, "Bad number of atoms in PME GPU");
173     cu_realloc_buffered((void **)&pmeGpu->kernelParams->atoms.d_coordinates, nullptr, sizeof(float),
174                         &pmeGpu->archSpecific->coordinatesSize, &pmeGpu->archSpecific->coordinatesSizeAlloc, newCoordinatesSize, pmeGpu->archSpecific->pmeStream, true);
175     if (c_usePadding)
176     {
177         const size_t paddingIndex = DIM * pmeGpu->kernelParams->atoms.nAtoms;
178         const size_t paddingCount = DIM * pmeGpu->nAtomsAlloc - paddingIndex;
179         if (paddingCount > 0)
180         {
181             cudaError_t stat = cudaMemsetAsync(pmeGpu->kernelParams->atoms.d_coordinates + paddingIndex, 0, paddingCount * sizeof(float), pmeGpu->archSpecific->pmeStream);
182             CU_RET_ERR(stat, "PME failed to clear the padded coordinates");
183         }
184     }
185 }
186
187 void pme_gpu_copy_input_coordinates(const PmeGpu *pmeGpu, const rvec *h_coordinates)
188 {
189     GMX_ASSERT(h_coordinates, "Bad host-side coordinate buffer in PME GPU");
190 #if GMX_DOUBLE
191     GMX_RELEASE_ASSERT(false, "Only single precision is supported");
192     GMX_UNUSED_VALUE(h_coordinates);
193 #else
194     cu_copy_H2D(pmeGpu->kernelParams->atoms.d_coordinates, const_cast<rvec *>(h_coordinates),
195                 pmeGpu->kernelParams->atoms.nAtoms * sizeof(rvec), pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
196 #endif
197 }
198
199 void pme_gpu_free_coordinates(const PmeGpu *pmeGpu)
200 {
201     cu_free_buffered(pmeGpu->kernelParams->atoms.d_coordinates, &pmeGpu->archSpecific->coordinatesSize, &pmeGpu->archSpecific->coordinatesSizeAlloc);
202 }
203
204 void pme_gpu_realloc_and_copy_input_coefficients(const PmeGpu *pmeGpu, const float *h_coefficients)
205 {
206     GMX_ASSERT(h_coefficients, "Bad host-side charge buffer in PME GPU");
207     const size_t newCoefficientsSize = pmeGpu->nAtomsAlloc;
208     GMX_ASSERT(newCoefficientsSize > 0, "Bad number of atoms in PME GPU");
209     cu_realloc_buffered((void **)&pmeGpu->kernelParams->atoms.d_coefficients, nullptr, sizeof(float),
210                         &pmeGpu->archSpecific->coefficientsSize, &pmeGpu->archSpecific->coefficientsSizeAlloc,
211                         newCoefficientsSize, pmeGpu->archSpecific->pmeStream, true);
212     cu_copy_H2D(pmeGpu->kernelParams->atoms.d_coefficients, const_cast<float *>(h_coefficients),
213                 pmeGpu->kernelParams->atoms.nAtoms * sizeof(float), pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
214     if (c_usePadding)
215     {
216         const size_t paddingIndex = pmeGpu->kernelParams->atoms.nAtoms;
217         const size_t paddingCount = pmeGpu->nAtomsAlloc - paddingIndex;
218         if (paddingCount > 0)
219         {
220             cudaError_t stat = cudaMemsetAsync(pmeGpu->kernelParams->atoms.d_coefficients + paddingIndex, 0, paddingCount * sizeof(float), pmeGpu->archSpecific->pmeStream);
221             CU_RET_ERR(stat, "PME failed to clear the padded charges");
222         }
223     }
224 }
225
226 void pme_gpu_free_coefficients(const PmeGpu *pmeGpu)
227 {
228     cu_free_buffered(pmeGpu->kernelParams->atoms.d_coefficients, &pmeGpu->archSpecific->coefficientsSize, &pmeGpu->archSpecific->coefficientsSizeAlloc);
229 }
230
231 void pme_gpu_realloc_spline_data(const PmeGpu *pmeGpu)
232 {
233     const int    order             = pmeGpu->common->pme_order;
234     const int    alignment         = pme_gpu_get_atoms_per_warp(pmeGpu);
235     const size_t nAtomsPadded      = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
236     const int    newSplineDataSize = DIM * order * nAtomsPadded;
237     GMX_ASSERT(newSplineDataSize > 0, "Bad number of atoms in PME GPU");
238     /* Two arrays of the same size */
239     const bool shouldRealloc        = (newSplineDataSize > pmeGpu->archSpecific->splineDataSize);
240     int        currentSizeTemp      = pmeGpu->archSpecific->splineDataSize;
241     int        currentSizeTempAlloc = pmeGpu->archSpecific->splineDataSizeAlloc;
242     cu_realloc_buffered((void **)&pmeGpu->kernelParams->atoms.d_theta, nullptr, sizeof(float),
243                         &currentSizeTemp, &currentSizeTempAlloc, newSplineDataSize, pmeGpu->archSpecific->pmeStream, true);
244     cu_realloc_buffered((void **)&pmeGpu->kernelParams->atoms.d_dtheta, nullptr, sizeof(float),
245                         &pmeGpu->archSpecific->splineDataSize, &pmeGpu->archSpecific->splineDataSizeAlloc, newSplineDataSize, pmeGpu->archSpecific->pmeStream, true);
246     // the host side reallocation
247     if (shouldRealloc)
248     {
249         pfree(pmeGpu->staging.h_theta);
250         pmalloc((void **)&pmeGpu->staging.h_theta, newSplineDataSize * sizeof(float));
251         pfree(pmeGpu->staging.h_dtheta);
252         pmalloc((void **)&pmeGpu->staging.h_dtheta, newSplineDataSize * sizeof(float));
253     }
254 }
255
256 void pme_gpu_free_spline_data(const PmeGpu *pmeGpu)
257 {
258     /* Two arrays of the same size */
259     cu_free_buffered(pmeGpu->kernelParams->atoms.d_theta);
260     cu_free_buffered(pmeGpu->kernelParams->atoms.d_dtheta, &pmeGpu->archSpecific->splineDataSize, &pmeGpu->archSpecific->splineDataSizeAlloc);
261     pfree(pmeGpu->staging.h_theta);
262     pfree(pmeGpu->staging.h_dtheta);
263 }
264
265 void pme_gpu_realloc_grid_indices(const PmeGpu *pmeGpu)
266 {
267     const size_t newIndicesSize = DIM * pmeGpu->nAtomsAlloc;
268     GMX_ASSERT(newIndicesSize > 0, "Bad number of atoms in PME GPU");
269     cu_realloc_buffered((void **)&pmeGpu->kernelParams->atoms.d_gridlineIndices, nullptr, sizeof(int),
270                         &pmeGpu->archSpecific->gridlineIndicesSize, &pmeGpu->archSpecific->gridlineIndicesSizeAlloc, newIndicesSize, pmeGpu->archSpecific->pmeStream, true);
271     pfree(pmeGpu->staging.h_gridlineIndices);
272     pmalloc((void **)&pmeGpu->staging.h_gridlineIndices, newIndicesSize * sizeof(int));
273 }
274
275 void pme_gpu_free_grid_indices(const PmeGpu *pmeGpu)
276 {
277     cu_free_buffered(pmeGpu->kernelParams->atoms.d_gridlineIndices, &pmeGpu->archSpecific->gridlineIndicesSize, &pmeGpu->archSpecific->gridlineIndicesSizeAlloc);
278     pfree(pmeGpu->staging.h_gridlineIndices);
279 }
280
281 void pme_gpu_realloc_grids(PmeGpu *pmeGpu)
282 {
283     auto     *kernelParamsPtr = pmeGpu->kernelParams.get();
284     const int newRealGridSize = kernelParamsPtr->grid.realGridSizePadded[XX] *
285         kernelParamsPtr->grid.realGridSizePadded[YY] *
286         kernelParamsPtr->grid.realGridSizePadded[ZZ];
287     const int newComplexGridSize = kernelParamsPtr->grid.complexGridSizePadded[XX] *
288         kernelParamsPtr->grid.complexGridSizePadded[YY] *
289         kernelParamsPtr->grid.complexGridSizePadded[ZZ] * 2;
290     // Multiplied by 2 because we count complex grid size for complex numbers, but all allocations/pointers are float
291     if (pmeGpu->archSpecific->performOutOfPlaceFFT)
292     {
293         /* 2 separate grids */
294         cu_realloc_buffered((void **)&kernelParamsPtr->grid.d_fourierGrid, nullptr, sizeof(float),
295                             &pmeGpu->archSpecific->complexGridSize, &pmeGpu->archSpecific->complexGridSizeAlloc,
296                             newComplexGridSize, pmeGpu->archSpecific->pmeStream, true);
297         cu_realloc_buffered((void **)&kernelParamsPtr->grid.d_realGrid, nullptr, sizeof(float),
298                             &pmeGpu->archSpecific->realGridSize, &pmeGpu->archSpecific->realGridSizeAlloc,
299                             newRealGridSize, pmeGpu->archSpecific->pmeStream, true);
300     }
301     else
302     {
303         /* A single buffer so that any grid will fit */
304         const int newGridsSize = std::max(newRealGridSize, newComplexGridSize);
305         cu_realloc_buffered((void **)&kernelParamsPtr->grid.d_realGrid, nullptr, sizeof(float),
306                             &pmeGpu->archSpecific->realGridSize, &pmeGpu->archSpecific->realGridSizeAlloc,
307                             newGridsSize, pmeGpu->archSpecific->pmeStream, true);
308         kernelParamsPtr->grid.d_fourierGrid   = kernelParamsPtr->grid.d_realGrid;
309         pmeGpu->archSpecific->complexGridSize = pmeGpu->archSpecific->realGridSize;
310         // the size might get used later for copying the grid
311     }
312 }
313
314 void pme_gpu_free_grids(const PmeGpu *pmeGpu)
315 {
316     if (pmeGpu->archSpecific->performOutOfPlaceFFT)
317     {
318         cu_free_buffered(pmeGpu->kernelParams->grid.d_fourierGrid);
319     }
320     cu_free_buffered(pmeGpu->kernelParams->grid.d_realGrid,
321                      &pmeGpu->archSpecific->realGridSize, &pmeGpu->archSpecific->realGridSizeAlloc);
322 }
323
324 void pme_gpu_clear_grids(const PmeGpu *pmeGpu)
325 {
326     cudaError_t stat = cudaMemsetAsync(pmeGpu->kernelParams->grid.d_realGrid, 0,
327                                        pmeGpu->archSpecific->realGridSize * sizeof(float), pmeGpu->archSpecific->pmeStream);
328     /* Should the complex grid be cleared in some weird case? */
329     CU_RET_ERR(stat, "cudaMemsetAsync on the PME grid error");
330 }
331
332 void pme_gpu_realloc_and_copy_fract_shifts(PmeGpu *pmeGpu)
333 {
334     pme_gpu_free_fract_shifts(pmeGpu);
335
336     auto        *kernelParamsPtr = pmeGpu->kernelParams.get();
337
338     const int    nx                  = kernelParamsPtr->grid.realGridSize[XX];
339     const int    ny                  = kernelParamsPtr->grid.realGridSize[YY];
340     const int    nz                  = kernelParamsPtr->grid.realGridSize[ZZ];
341     const int    cellCount           = c_pmeNeighborUnitcellCount;
342     const int    gridDataOffset[DIM] = {0, cellCount * nx, cellCount * (nx + ny)};
343
344     memcpy(kernelParamsPtr->grid.tablesOffsets, &gridDataOffset, sizeof(gridDataOffset));
345
346     const int    newFractShiftsSize  = cellCount * (nx + ny + nz);
347
348     initParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable,
349                          kernelParamsPtr->fractShiftsTableTexture,
350                          pmeGpu->common->fsh.data(),
351                          newFractShiftsSize,
352                          pmeGpu->deviceInfo);
353
354     initParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
355                          kernelParamsPtr->gridlineIndicesTableTexture,
356                          pmeGpu->common->nn.data(),
357                          newFractShiftsSize,
358                          pmeGpu->deviceInfo);
359 }
360
361 void pme_gpu_free_fract_shifts(const PmeGpu *pmeGpu)
362 {
363     auto *kernelParamsPtr = pmeGpu->kernelParams.get();
364     destroyParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable,
365                             kernelParamsPtr->fractShiftsTableTexture,
366                             pmeGpu->deviceInfo);
367     destroyParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
368                             kernelParamsPtr->gridlineIndicesTableTexture,
369                             pmeGpu->deviceInfo);
370 }
371
372 bool pme_gpu_stream_query(const PmeGpu *pmeGpu)
373 {
374     return haveStreamTasksCompleted(pmeGpu->archSpecific->pmeStream);
375 }
376
377 void pme_gpu_copy_input_gather_grid(const PmeGpu *pmeGpu, float *h_grid)
378 {
379     const size_t gridSize = pmeGpu->archSpecific->realGridSize * sizeof(float);
380     cu_copy_H2D(pmeGpu->kernelParams->grid.d_realGrid, h_grid, gridSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
381 }
382
383 void pme_gpu_copy_output_spread_grid(const PmeGpu *pmeGpu, float *h_grid)
384 {
385     const size_t gridSize = pmeGpu->archSpecific->realGridSize * sizeof(float);
386     cu_copy_D2H(h_grid, pmeGpu->kernelParams->grid.d_realGrid, gridSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
387     cudaError_t  stat = cudaEventRecord(pmeGpu->archSpecific->syncSpreadGridD2H, pmeGpu->archSpecific->pmeStream);
388     CU_RET_ERR(stat, "PME spread grid sync event record failure");
389 }
390
391 void pme_gpu_copy_output_spread_atom_data(const PmeGpu *pmeGpu)
392 {
393     const int    alignment       = pme_gpu_get_atoms_per_warp(pmeGpu);
394     const size_t nAtomsPadded    = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
395     const size_t splinesSize     = DIM * nAtomsPadded * pmeGpu->common->pme_order * sizeof(float);
396     auto        *kernelParamsPtr = pmeGpu->kernelParams.get();
397     cu_copy_D2H(pmeGpu->staging.h_dtheta, kernelParamsPtr->atoms.d_dtheta, splinesSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
398     cu_copy_D2H(pmeGpu->staging.h_theta, kernelParamsPtr->atoms.d_theta, splinesSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
399     cu_copy_D2H(pmeGpu->staging.h_gridlineIndices, kernelParamsPtr->atoms.d_gridlineIndices,
400                 kernelParamsPtr->atoms.nAtoms * DIM * sizeof(int), pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
401 }
402
403 void pme_gpu_copy_input_gather_atom_data(const PmeGpu *pmeGpu)
404 {
405     const int    alignment       = pme_gpu_get_atoms_per_warp(pmeGpu);
406     const size_t nAtomsPadded    = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
407     const size_t splinesSize     = DIM * nAtomsPadded * pmeGpu->common->pme_order * sizeof(float);
408     auto        *kernelParamsPtr = pmeGpu->kernelParams.get();
409     if (c_usePadding)
410     {
411         const size_t gridlineIndicesSizePerAtom = DIM * sizeof(int);
412         const size_t splineDataSizePerAtom      = pmeGpu->common->pme_order * DIM * sizeof(float);
413         // TODO: could clear only the padding and not the whole thing, but this is a test-exclusive code anyway
414         CU_RET_ERR(cudaMemsetAsync(kernelParamsPtr->atoms.d_gridlineIndices, 0, pmeGpu->nAtomsAlloc * gridlineIndicesSizePerAtom, pmeGpu->archSpecific->pmeStream),
415                    "PME failed to clear the gridline indices");
416         CU_RET_ERR(cudaMemsetAsync(kernelParamsPtr->atoms.d_dtheta, 0, pmeGpu->nAtomsAlloc * splineDataSizePerAtom, pmeGpu->archSpecific->pmeStream),
417                    "PME failed to clear the spline derivatives");
418         CU_RET_ERR(cudaMemsetAsync(kernelParamsPtr->atoms.d_theta, 0, pmeGpu->nAtomsAlloc * splineDataSizePerAtom, pmeGpu->archSpecific->pmeStream),
419                    "PME failed to clear the spline values");
420     }
421     cu_copy_H2D(kernelParamsPtr->atoms.d_dtheta, pmeGpu->staging.h_dtheta, splinesSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
422     cu_copy_H2D(kernelParamsPtr->atoms.d_theta, pmeGpu->staging.h_theta, splinesSize, pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
423     cu_copy_H2D(kernelParamsPtr->atoms.d_gridlineIndices, pmeGpu->staging.h_gridlineIndices,
424                 kernelParamsPtr->atoms.nAtoms * DIM * sizeof(int), pmeGpu->settings.transferKind, pmeGpu->archSpecific->pmeStream);
425 }
426
427 void pme_gpu_sync_spread_grid(const PmeGpu *pmeGpu)
428 {
429     cudaError_t stat = cudaEventSynchronize(pmeGpu->archSpecific->syncSpreadGridD2H);
430     CU_RET_ERR(stat, "Error while waiting for the PME GPU spread grid to be copied to the host");
431 }
432
433 void pme_gpu_init_internal(PmeGpu *pmeGpu)
434 {
435     /* Allocate the target-specific structures */
436     pmeGpu->archSpecific.reset(new PmeGpuSpecific());
437     pmeGpu->kernelParams.reset(new PmeGpuKernelParams());
438
439     pmeGpu->archSpecific->performOutOfPlaceFFT = true;
440     /* This should give better performance, according to the cuFFT documentation.
441      * The performance seems to be the same though.
442      * TODO: PME could also try to pick up nice grid sizes (with factors of 2, 3, 5, 7).
443      */
444
445     /* WARNING: CUDA timings are incorrect with multiple streams.
446      *          This is the main reason why they are disabled by default.
447      */
448     // TODO: Consider turning on by default when we can detect nr of streams.
449     pmeGpu->archSpecific->useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
450
451     /* Creating a PME CUDA stream */
452     cudaError_t stat;
453     int         highest_priority, lowest_priority;
454     stat = cudaDeviceGetStreamPriorityRange(&lowest_priority, &highest_priority);
455     CU_RET_ERR(stat, "PME cudaDeviceGetStreamPriorityRange failed");
456     stat = cudaStreamCreateWithPriority(&pmeGpu->archSpecific->pmeStream,
457                                         cudaStreamDefault, //cudaStreamNonBlocking,
458                                         highest_priority);
459     CU_RET_ERR(stat, "cudaStreamCreateWithPriority on the PME stream failed");
460 }
461
462 void pme_gpu_destroy_specific(const PmeGpu *pmeGpu)
463 {
464     /* Destroy the CUDA stream */
465     cudaError_t stat = cudaStreamDestroy(pmeGpu->archSpecific->pmeStream);
466     CU_RET_ERR(stat, "PME cudaStreamDestroy error");
467 }
468
469 void pme_gpu_init_sync_events(const PmeGpu *pmeGpu)
470 {
471     const auto  eventFlags = cudaEventDisableTiming;
472     CU_RET_ERR(cudaEventCreateWithFlags(&pmeGpu->archSpecific->syncSpreadGridD2H, eventFlags), "cudaEventCreate on syncSpreadGridD2H failed");
473 }
474
475 void pme_gpu_destroy_sync_events(const PmeGpu *pmeGpu)
476 {
477     CU_RET_ERR(cudaEventDestroy(pmeGpu->archSpecific->syncSpreadGridD2H), "cudaEventDestroy failed on syncSpreadGridD2H");
478 }
479
480 void pme_gpu_reinit_3dfft(const PmeGpu *pmeGpu)
481 {
482     if (pme_gpu_performs_FFT(pmeGpu))
483     {
484         pmeGpu->archSpecific->fftSetup.resize(0);
485         for (int i = 0; i < pmeGpu->common->ngrids; i++)
486         {
487             pmeGpu->archSpecific->fftSetup.push_back(std::unique_ptr<GpuParallel3dFft>(new GpuParallel3dFft(pmeGpu)));
488         }
489     }
490 }
491
492 void pme_gpu_destroy_3dfft(const PmeGpu *pmeGpu)
493 {
494     pmeGpu->archSpecific->fftSetup.resize(0);
495 }