2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 * \brief This file contains internal function implementations
39 * for performing the PME calculations on GPU.
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.
44 * \author Aleksei Iupinov <a.yupinov@gmail.com>
45 * \ingroup module_ewald
50 #include "pme-gpu-internal.h"
57 #include "gromacs/compat/make_unique.h"
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"
69 #if GMX_GPU == GMX_GPU_CUDA
70 #include "gromacs/gpu_utils/pmalloc_cuda.h"
73 #elif GMX_GPU == GMX_GPU_OPENCL
74 #include "gromacs/gpu_utils/gmxopencl.h"
77 #include "gromacs/ewald/pme.h"
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"
88 #include "pme-internal.h"
91 * Wrapper for getting a pointer to the plain C++ part of the GPU kernel parameters structure.
93 * \param[in] pmeGpu The PME GPU structure.
94 * \returns The pointer to the kernel parameters.
96 static PmeGpuKernelParamsBase *pme_gpu_get_kernel_params_base_ptr(const PmeGpu *pmeGpu)
98 // reinterpret_cast is needed because the derived CUDA structure is not known in this file
99 auto *kernelParamsPtr = reinterpret_cast<PmeGpuKernelParamsBase *>(pmeGpu->kernelParams.get());
100 return kernelParamsPtr;
103 int pme_gpu_get_atom_data_alignment(const PmeGpu * /*unused*/)
105 //TODO: this can be simplified, as PME_ATOM_DATA_ALIGNMENT is now constant
106 return PME_ATOM_DATA_ALIGNMENT;
109 int pme_gpu_get_atoms_per_warp(const PmeGpu *pmeGpu)
111 const int order = pmeGpu->common->pme_order;
112 GMX_ASSERT(order > 0, "Invalid PME order");
113 return pmeGpu->programHandle_->impl_->warpSize / PME_SPREADGATHER_THREADS_PER_ATOM;
116 void pme_gpu_synchronize(const PmeGpu *pmeGpu)
118 gpuStreamSynchronize(pmeGpu->archSpecific->pmeStream);
121 void pme_gpu_alloc_energy_virial(PmeGpu *pmeGpu)
123 const size_t energyAndVirialSize = c_virialAndEnergyCount * sizeof(float);
124 allocateDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy, c_virialAndEnergyCount, pmeGpu->archSpecific->context);
125 pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_virialAndEnergy), energyAndVirialSize);
128 void pme_gpu_free_energy_virial(PmeGpu *pmeGpu)
130 freeDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy);
131 pfree(pmeGpu->staging.h_virialAndEnergy);
132 pmeGpu->staging.h_virialAndEnergy = nullptr;
135 void pme_gpu_clear_energy_virial(const PmeGpu *pmeGpu)
137 clearDeviceBufferAsync(&pmeGpu->kernelParams->constants.d_virialAndEnergy, 0,
138 c_virialAndEnergyCount, pmeGpu->archSpecific->pmeStream);
141 void pme_gpu_realloc_and_copy_bspline_values(PmeGpu *pmeGpu)
143 const int splineValuesOffset[DIM] = {
145 pmeGpu->kernelParams->grid.realGridSize[XX],
146 pmeGpu->kernelParams->grid.realGridSize[XX] + pmeGpu->kernelParams->grid.realGridSize[YY]
148 memcpy(&pmeGpu->kernelParams->grid.splineValuesOffset, &splineValuesOffset, sizeof(splineValuesOffset));
150 const int newSplineValuesSize = pmeGpu->kernelParams->grid.realGridSize[XX] +
151 pmeGpu->kernelParams->grid.realGridSize[YY] +
152 pmeGpu->kernelParams->grid.realGridSize[ZZ];
153 const bool shouldRealloc = (newSplineValuesSize > pmeGpu->archSpecific->splineValuesSize);
154 reallocateDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli, newSplineValuesSize,
155 &pmeGpu->archSpecific->splineValuesSize, &pmeGpu->archSpecific->splineValuesSizeAlloc, pmeGpu->archSpecific->context);
158 /* Reallocate the host buffer */
159 pfree(pmeGpu->staging.h_splineModuli);
160 pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_splineModuli), newSplineValuesSize * sizeof(float));
162 for (int i = 0; i < DIM; i++)
164 memcpy(pmeGpu->staging.h_splineModuli + splineValuesOffset[i], pmeGpu->common->bsp_mod[i].data(), pmeGpu->common->bsp_mod[i].size() * sizeof(float));
166 /* TODO: pin original buffer instead! */
167 copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli, pmeGpu->staging.h_splineModuli,
168 0, newSplineValuesSize,
169 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
172 void pme_gpu_free_bspline_values(const PmeGpu *pmeGpu)
174 pfree(pmeGpu->staging.h_splineModuli);
175 freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli);
178 void pme_gpu_realloc_forces(PmeGpu *pmeGpu)
180 const size_t newForcesSize = pmeGpu->nAtomsAlloc * DIM;
181 GMX_ASSERT(newForcesSize > 0, "Bad number of atoms in PME GPU");
182 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces, newForcesSize,
183 &pmeGpu->archSpecific->forcesSize, &pmeGpu->archSpecific->forcesSizeAlloc, pmeGpu->archSpecific->context);
184 pmeGpu->staging.h_forces.reserve(pmeGpu->nAtomsAlloc);
185 pmeGpu->staging.h_forces.resize(pmeGpu->kernelParams->atoms.nAtoms);
188 void pme_gpu_free_forces(const PmeGpu *pmeGpu)
190 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces);
193 void pme_gpu_copy_input_forces(PmeGpu *pmeGpu)
195 GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
196 float *h_forcesFloat = reinterpret_cast<float *>(pmeGpu->staging.h_forces.data());
197 copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces, h_forcesFloat,
198 0, DIM * pmeGpu->kernelParams->atoms.nAtoms,
199 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
202 void pme_gpu_copy_output_forces(PmeGpu *pmeGpu)
204 GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
205 float *h_forcesFloat = reinterpret_cast<float *>(pmeGpu->staging.h_forces.data());
206 copyFromDeviceBuffer(h_forcesFloat, &pmeGpu->kernelParams->atoms.d_forces,
207 0, DIM * pmeGpu->kernelParams->atoms.nAtoms,
208 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
211 void pme_gpu_realloc_coordinates(const PmeGpu *pmeGpu)
213 const size_t newCoordinatesSize = pmeGpu->nAtomsAlloc * DIM;
214 GMX_ASSERT(newCoordinatesSize > 0, "Bad number of atoms in PME GPU");
215 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates, newCoordinatesSize,
216 &pmeGpu->archSpecific->coordinatesSize, &pmeGpu->archSpecific->coordinatesSizeAlloc, pmeGpu->archSpecific->context);
219 const size_t paddingIndex = DIM * pmeGpu->kernelParams->atoms.nAtoms;
220 const size_t paddingCount = DIM * pmeGpu->nAtomsAlloc - paddingIndex;
221 if (paddingCount > 0)
223 clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coordinates, paddingIndex,
224 paddingCount, pmeGpu->archSpecific->pmeStream);
229 void pme_gpu_copy_input_coordinates(const PmeGpu *pmeGpu, const rvec *h_coordinates)
231 GMX_ASSERT(h_coordinates, "Bad host-side coordinate buffer in PME GPU");
233 GMX_RELEASE_ASSERT(false, "Only single precision is supported");
234 GMX_UNUSED_VALUE(h_coordinates);
236 const float *h_coordinatesFloat = reinterpret_cast<const float *>(h_coordinates);
237 copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates, h_coordinatesFloat,
238 0, pmeGpu->kernelParams->atoms.nAtoms * DIM,
239 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
243 void pme_gpu_free_coordinates(const PmeGpu *pmeGpu)
245 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coordinates);
248 void pme_gpu_realloc_and_copy_input_coefficients(const PmeGpu *pmeGpu, const float *h_coefficients)
250 GMX_ASSERT(h_coefficients, "Bad host-side charge buffer in PME GPU");
251 const size_t newCoefficientsSize = pmeGpu->nAtomsAlloc;
252 GMX_ASSERT(newCoefficientsSize > 0, "Bad number of atoms in PME GPU");
253 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients, newCoefficientsSize,
254 &pmeGpu->archSpecific->coefficientsSize, &pmeGpu->archSpecific->coefficientsSizeAlloc, pmeGpu->archSpecific->context);
255 copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients, const_cast<float *>(h_coefficients),
256 0, pmeGpu->kernelParams->atoms.nAtoms,
257 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
260 const size_t paddingIndex = pmeGpu->kernelParams->atoms.nAtoms;
261 const size_t paddingCount = pmeGpu->nAtomsAlloc - paddingIndex;
262 if (paddingCount > 0)
264 clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coefficients, paddingIndex,
265 paddingCount, pmeGpu->archSpecific->pmeStream);
270 void pme_gpu_free_coefficients(const PmeGpu *pmeGpu)
272 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients);
275 void pme_gpu_realloc_spline_data(PmeGpu *pmeGpu)
277 const int order = pmeGpu->common->pme_order;
278 const int alignment = pme_gpu_get_atoms_per_warp(pmeGpu);
279 const size_t nAtomsPadded = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
280 const int newSplineDataSize = DIM * order * nAtomsPadded;
281 GMX_ASSERT(newSplineDataSize > 0, "Bad number of atoms in PME GPU");
282 /* Two arrays of the same size */
283 const bool shouldRealloc = (newSplineDataSize > pmeGpu->archSpecific->splineDataSize);
284 int currentSizeTemp = pmeGpu->archSpecific->splineDataSize;
285 int currentSizeTempAlloc = pmeGpu->archSpecific->splineDataSizeAlloc;
286 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta, newSplineDataSize,
287 ¤tSizeTemp, ¤tSizeTempAlloc, pmeGpu->archSpecific->context);
288 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta, newSplineDataSize,
289 &pmeGpu->archSpecific->splineDataSize, &pmeGpu->archSpecific->splineDataSizeAlloc, pmeGpu->archSpecific->context);
290 // the host side reallocation
293 pfree(pmeGpu->staging.h_theta);
294 pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_theta), newSplineDataSize * sizeof(float));
295 pfree(pmeGpu->staging.h_dtheta);
296 pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_dtheta), newSplineDataSize * sizeof(float));
300 void pme_gpu_free_spline_data(const PmeGpu *pmeGpu)
302 /* Two arrays of the same size */
303 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta);
304 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta);
305 pfree(pmeGpu->staging.h_theta);
306 pfree(pmeGpu->staging.h_dtheta);
309 void pme_gpu_realloc_grid_indices(PmeGpu *pmeGpu)
311 const size_t newIndicesSize = DIM * pmeGpu->nAtomsAlloc;
312 GMX_ASSERT(newIndicesSize > 0, "Bad number of atoms in PME GPU");
313 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices, newIndicesSize,
314 &pmeGpu->archSpecific->gridlineIndicesSize, &pmeGpu->archSpecific->gridlineIndicesSizeAlloc, pmeGpu->archSpecific->context);
315 pfree(pmeGpu->staging.h_gridlineIndices);
316 pmalloc(reinterpret_cast<void **>(&pmeGpu->staging.h_gridlineIndices), newIndicesSize * sizeof(int));
319 void pme_gpu_free_grid_indices(const PmeGpu *pmeGpu)
321 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices);
322 pfree(pmeGpu->staging.h_gridlineIndices);
325 void pme_gpu_realloc_grids(PmeGpu *pmeGpu)
327 auto *kernelParamsPtr = pmeGpu->kernelParams.get();
328 const int newRealGridSize = kernelParamsPtr->grid.realGridSizePadded[XX] *
329 kernelParamsPtr->grid.realGridSizePadded[YY] *
330 kernelParamsPtr->grid.realGridSizePadded[ZZ];
331 const int newComplexGridSize = kernelParamsPtr->grid.complexGridSizePadded[XX] *
332 kernelParamsPtr->grid.complexGridSizePadded[YY] *
333 kernelParamsPtr->grid.complexGridSizePadded[ZZ] * 2;
334 // Multiplied by 2 because we count complex grid size for complex numbers, but all allocations/pointers are float
335 if (pmeGpu->archSpecific->performOutOfPlaceFFT)
337 /* 2 separate grids */
338 reallocateDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid, newComplexGridSize,
339 &pmeGpu->archSpecific->complexGridSize, &pmeGpu->archSpecific->complexGridSizeAlloc, pmeGpu->archSpecific->context);
340 reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid, newRealGridSize,
341 &pmeGpu->archSpecific->realGridSize, &pmeGpu->archSpecific->realGridSizeAlloc, pmeGpu->archSpecific->context);
345 /* A single buffer so that any grid will fit */
346 const int newGridsSize = std::max(newRealGridSize, newComplexGridSize);
347 reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid, newGridsSize,
348 &pmeGpu->archSpecific->realGridSize, &pmeGpu->archSpecific->realGridSizeAlloc, pmeGpu->archSpecific->context);
349 kernelParamsPtr->grid.d_fourierGrid = kernelParamsPtr->grid.d_realGrid;
350 pmeGpu->archSpecific->complexGridSize = pmeGpu->archSpecific->realGridSize;
351 // the size might get used later for copying the grid
355 void pme_gpu_free_grids(const PmeGpu *pmeGpu)
357 if (pmeGpu->archSpecific->performOutOfPlaceFFT)
359 freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_fourierGrid);
361 freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid);
364 void pme_gpu_clear_grids(const PmeGpu *pmeGpu)
366 clearDeviceBufferAsync(&pmeGpu->kernelParams->grid.d_realGrid, 0,
367 pmeGpu->archSpecific->realGridSize, pmeGpu->archSpecific->pmeStream);
370 void pme_gpu_realloc_and_copy_fract_shifts(PmeGpu *pmeGpu)
372 pme_gpu_free_fract_shifts(pmeGpu);
374 auto *kernelParamsPtr = pmeGpu->kernelParams.get();
376 const int nx = kernelParamsPtr->grid.realGridSize[XX];
377 const int ny = kernelParamsPtr->grid.realGridSize[YY];
378 const int nz = kernelParamsPtr->grid.realGridSize[ZZ];
379 const int cellCount = c_pmeNeighborUnitcellCount;
380 const int gridDataOffset[DIM] = {0, cellCount * nx, cellCount * (nx + ny)};
382 memcpy(kernelParamsPtr->grid.tablesOffsets, &gridDataOffset, sizeof(gridDataOffset));
384 const int newFractShiftsSize = cellCount * (nx + ny + nz);
386 #if GMX_GPU == GMX_GPU_CUDA
387 initParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable,
388 kernelParamsPtr->fractShiftsTableTexture,
389 pmeGpu->common->fsh.data(),
393 initParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
394 kernelParamsPtr->gridlineIndicesTableTexture,
395 pmeGpu->common->nn.data(),
398 #elif GMX_GPU == GMX_GPU_OPENCL
399 // No dedicated texture routines....
400 allocateDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable, newFractShiftsSize, pmeGpu->archSpecific->context);
401 allocateDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable, newFractShiftsSize, pmeGpu->archSpecific->context);
402 copyToDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable, pmeGpu->common->fsh.data(),
403 0, newFractShiftsSize,
404 pmeGpu->archSpecific->pmeStream, GpuApiCallBehavior::Async, nullptr);
405 copyToDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable, pmeGpu->common->nn.data(),
406 0, newFractShiftsSize,
407 pmeGpu->archSpecific->pmeStream, GpuApiCallBehavior::Async, nullptr);
411 void pme_gpu_free_fract_shifts(const PmeGpu *pmeGpu)
413 auto *kernelParamsPtr = pmeGpu->kernelParams.get();
414 #if GMX_GPU == GMX_GPU_CUDA
415 destroyParamLookupTable(kernelParamsPtr->grid.d_fractShiftsTable,
416 kernelParamsPtr->fractShiftsTableTexture,
418 destroyParamLookupTable(kernelParamsPtr->grid.d_gridlineIndicesTable,
419 kernelParamsPtr->gridlineIndicesTableTexture,
421 #elif GMX_GPU == GMX_GPU_OPENCL
422 freeDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable);
423 freeDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable);
427 bool pme_gpu_stream_query(const PmeGpu *pmeGpu)
429 return haveStreamTasksCompleted(pmeGpu->archSpecific->pmeStream);
432 void pme_gpu_copy_input_gather_grid(const PmeGpu *pmeGpu, float *h_grid)
434 copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid, h_grid,
435 0, pmeGpu->archSpecific->realGridSize,
436 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
439 void pme_gpu_copy_output_spread_grid(const PmeGpu *pmeGpu, float *h_grid)
441 copyFromDeviceBuffer(h_grid, &pmeGpu->kernelParams->grid.d_realGrid,
442 0, pmeGpu->archSpecific->realGridSize,
443 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
444 pmeGpu->archSpecific->syncSpreadGridD2H.markEvent(pmeGpu->archSpecific->pmeStream);
447 void pme_gpu_copy_output_spread_atom_data(const PmeGpu *pmeGpu)
449 const int alignment = pme_gpu_get_atoms_per_warp(pmeGpu);
450 const size_t nAtomsPadded = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
451 const size_t splinesCount = DIM * nAtomsPadded * pmeGpu->common->pme_order;
452 auto *kernelParamsPtr = pmeGpu->kernelParams.get();
453 copyFromDeviceBuffer(pmeGpu->staging.h_dtheta, &kernelParamsPtr->atoms.d_dtheta,
455 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
456 copyFromDeviceBuffer(pmeGpu->staging.h_theta, &kernelParamsPtr->atoms.d_theta,
458 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
459 copyFromDeviceBuffer(pmeGpu->staging.h_gridlineIndices, &kernelParamsPtr->atoms.d_gridlineIndices,
460 0, kernelParamsPtr->atoms.nAtoms * DIM,
461 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
464 void pme_gpu_copy_input_gather_atom_data(const PmeGpu *pmeGpu)
466 const int alignment = pme_gpu_get_atoms_per_warp(pmeGpu);
467 const size_t nAtomsPadded = ((pmeGpu->nAtomsAlloc + alignment - 1) / alignment) * alignment;
468 const size_t splinesCount = DIM * nAtomsPadded * pmeGpu->common->pme_order;
469 auto *kernelParamsPtr = pmeGpu->kernelParams.get();
472 // TODO: could clear only the padding and not the whole thing, but this is a test-exclusive code anyway
473 clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_gridlineIndices, 0,
474 pmeGpu->nAtomsAlloc * DIM, pmeGpu->archSpecific->pmeStream);
475 clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_dtheta, 0,
476 pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM, pmeGpu->archSpecific->pmeStream);
477 clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_theta, 0,
478 pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM, pmeGpu->archSpecific->pmeStream);
480 copyToDeviceBuffer(&kernelParamsPtr->atoms.d_dtheta, pmeGpu->staging.h_dtheta,
482 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
483 copyToDeviceBuffer(&kernelParamsPtr->atoms.d_theta, pmeGpu->staging.h_theta,
485 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
486 copyToDeviceBuffer(&kernelParamsPtr->atoms.d_gridlineIndices, pmeGpu->staging.h_gridlineIndices,
487 0, kernelParamsPtr->atoms.nAtoms * DIM,
488 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
491 void pme_gpu_sync_spread_grid(const PmeGpu *pmeGpu)
493 pmeGpu->archSpecific->syncSpreadGridD2H.waitForEvent();
496 void pme_gpu_init_internal(PmeGpu *pmeGpu)
498 /* Allocate the target-specific structures */
499 pmeGpu->archSpecific.reset(new PmeGpuSpecific());
500 pmeGpu->kernelParams.reset(new PmeGpuKernelParams());
502 pmeGpu->archSpecific->performOutOfPlaceFFT = true;
503 /* This should give better performance, according to the cuFFT documentation.
504 * The performance seems to be the same though.
505 * TODO: PME could also try to pick up nice grid sizes (with factors of 2, 3, 5, 7).
508 // timing enabling - TODO put this in gpu_utils (even though generally this is just option handling?) and reuse in NB
509 if (GMX_GPU == GMX_GPU_CUDA)
511 /* WARNING: CUDA timings are incorrect with multiple streams.
512 * This is the main reason why they are disabled by default.
514 // TODO: Consider turning on by default when we can detect nr of streams.
515 pmeGpu->archSpecific->useTiming = (getenv("GMX_ENABLE_GPU_TIMING") != nullptr);
517 else if (GMX_GPU == GMX_GPU_OPENCL)
519 pmeGpu->archSpecific->useTiming = (getenv("GMX_DISABLE_GPU_TIMING") == nullptr);
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;
525 #if GMX_GPU == GMX_GPU_CUDA
526 // Prepare to use the device that this PME task was assigned earlier.
527 CU_RET_ERR(cudaSetDevice(pmeGpu->deviceInfo->id), "Switching to PME CUDA device");
530 #if GMX_GPU == GMX_GPU_CUDA
531 pmeGpu->maxGridWidthX = pmeGpu->deviceInfo->prop.maxGridSize[0];
532 #elif GMX_GPU == GMX_GPU_OPENCL
533 //TODO we'll need work size checks for OpenCL too
536 /* Creating a PME GPU stream:
537 * - default high priority with CUDA
538 * - no priorities implemented yet with OpenCL; see #2532
540 #if GMX_GPU == GMX_GPU_CUDA
542 int highest_priority, lowest_priority;
543 stat = cudaDeviceGetStreamPriorityRange(&lowest_priority, &highest_priority);
544 CU_RET_ERR(stat, "PME cudaDeviceGetStreamPriorityRange failed");
545 stat = cudaStreamCreateWithPriority(&pmeGpu->archSpecific->pmeStream,
546 cudaStreamDefault, //cudaStreamNonBlocking,
548 CU_RET_ERR(stat, "cudaStreamCreateWithPriority on the PME stream failed");
549 #elif GMX_GPU == GMX_GPU_OPENCL
550 cl_command_queue_properties queueProperties = pmeGpu->archSpecific->useTiming ? CL_QUEUE_PROFILING_ENABLE : 0;
551 cl_device_id device_id = pmeGpu->deviceInfo->ocl_gpu_id.ocl_device_id;
553 pmeGpu->archSpecific->pmeStream = clCreateCommandQueue(pmeGpu->archSpecific->context,
554 device_id, queueProperties, &clError);
555 if (clError != CL_SUCCESS)
557 GMX_THROW(gmx::InternalError("Failed to create PME command queue"));
562 void pme_gpu_destroy_specific(const PmeGpu *pmeGpu)
564 #if GMX_GPU == GMX_GPU_CUDA
565 /* Destroy the CUDA stream */
566 cudaError_t stat = cudaStreamDestroy(pmeGpu->archSpecific->pmeStream);
567 CU_RET_ERR(stat, "PME cudaStreamDestroy error");
568 #elif GMX_GPU == GMX_GPU_OPENCL
569 cl_int clError = clReleaseCommandQueue(pmeGpu->archSpecific->pmeStream);
570 if (clError != CL_SUCCESS)
572 gmx_warning("Failed to destroy PME command queue");
577 void pme_gpu_reinit_3dfft(const PmeGpu *pmeGpu)
579 if (pme_gpu_performs_FFT(pmeGpu))
581 pmeGpu->archSpecific->fftSetup.resize(0);
582 for (int i = 0; i < pmeGpu->common->ngrids; i++)
584 pmeGpu->archSpecific->fftSetup.push_back(gmx::compat::make_unique<GpuParallel3dFft>(pmeGpu));
589 void pme_gpu_destroy_3dfft(const PmeGpu *pmeGpu)
591 pmeGpu->archSpecific->fftSetup.resize(0);
594 int getSplineParamFullIndex(int order, int splineIndex, int dimIndex, int atomIndex, int atomsPerWarp)
600 constexpr int fixedOrder = 4;
601 GMX_UNUSED_VALUE(fixedOrder);
603 const int atomWarpIndex = atomIndex % atomsPerWarp;
604 const int warpIndex = atomIndex / atomsPerWarp;
605 int indexBase, result;
606 switch (atomsPerWarp)
609 indexBase = getSplineParamIndexBase<fixedOrder, 2>(warpIndex, atomWarpIndex);
610 result = getSplineParamIndex<fixedOrder, 2>(indexBase, dimIndex, splineIndex);
614 GMX_THROW(gmx::NotImplementedError(gmx::formatString("Test function call not unrolled for atomsPerWarp = %d in getSplineParamFullIndex", atomsPerWarp)));
619 gmx::ArrayRef<gmx::RVec> pme_gpu_get_forces(PmeGpu *pmeGpu)
621 return pmeGpu->staging.h_forces;
624 void pme_gpu_get_energy_virial(const PmeGpu *pmeGpu, real *energy, matrix virial)
626 for (int j = 0; j < c_virialAndEnergyCount; j++)
628 GMX_ASSERT(std::isfinite(pmeGpu->staging.h_virialAndEnergy[j]), "PME GPU produces incorrect energy/virial.");
631 GMX_ASSERT(energy, "Invalid energy output pointer in PME GPU");
633 virial[XX][XX] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
634 virial[YY][YY] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
635 virial[ZZ][ZZ] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
636 virial[XX][YY] = virial[YY][XX] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
637 virial[XX][ZZ] = virial[ZZ][XX] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
638 virial[YY][ZZ] = virial[ZZ][YY] = 0.25f * pmeGpu->staging.h_virialAndEnergy[j++];
639 *energy = 0.5f * pmeGpu->staging.h_virialAndEnergy[j++];
642 void pme_gpu_update_input_box(PmeGpu gmx_unused *pmeGpu,
643 const matrix gmx_unused box)
646 GMX_THROW(gmx::NotImplementedError("PME is implemented for single-precision only on GPU"));
649 pmeGpu->common->boxScaler->scaleBox(box, scaledBox);
650 auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
651 kernelParamsPtr->current.boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
652 GMX_ASSERT(kernelParamsPtr->current.boxVolume != 0.0f, "Zero volume of the unit cell");
654 gmx::invertBoxMatrix(scaledBox, recipBox);
656 /* The GPU recipBox is transposed as compared to the CPU recipBox.
657 * Spread uses matrix columns (while solve and gather use rows).
658 * There is no particular reason for this; it might be further rethought/optimized for better access patterns.
660 const real newRecipBox[DIM][DIM] =
662 {recipBox[XX][XX], recipBox[YY][XX], recipBox[ZZ][XX]},
663 { 0.0, recipBox[YY][YY], recipBox[ZZ][YY]},
664 { 0.0, 0.0, recipBox[ZZ][ZZ]}
666 memcpy(kernelParamsPtr->current.recipBox, newRecipBox, sizeof(matrix));
670 /*! \brief \libinternal
671 * (Re-)initializes all the PME GPU data related to the grid size and cut-off.
673 * \param[in] pmeGpu The PME GPU structure.
675 static void pme_gpu_reinit_grids(PmeGpu *pmeGpu)
677 auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
678 kernelParamsPtr->grid.ewaldFactor = (M_PI * M_PI) / (pmeGpu->common->ewaldcoeff_q * pmeGpu->common->ewaldcoeff_q);
680 /* The grid size variants */
681 for (int i = 0; i < DIM; i++)
683 kernelParamsPtr->grid.realGridSize[i] = pmeGpu->common->nk[i];
684 kernelParamsPtr->grid.realGridSizeFP[i] = static_cast<float>(kernelParamsPtr->grid.realGridSize[i]);
685 kernelParamsPtr->grid.realGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
687 // The complex grid currently uses no padding;
688 // if it starts to do so, then another test should be added for that
689 kernelParamsPtr->grid.complexGridSize[i] = kernelParamsPtr->grid.realGridSize[i];
690 kernelParamsPtr->grid.complexGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
692 /* FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
693 if (!pme_gpu_performs_FFT(pmeGpu))
695 // This allows for GPU spreading grid and CPU fftgrid to have the same layout, so that we can copy the data directly
696 kernelParamsPtr->grid.realGridSizePadded[ZZ] = (kernelParamsPtr->grid.realGridSize[ZZ] / 2 + 1) * 2;
699 /* GPU FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
700 kernelParamsPtr->grid.complexGridSize[ZZ] /= 2;
701 kernelParamsPtr->grid.complexGridSize[ZZ]++;
702 kernelParamsPtr->grid.complexGridSizePadded[ZZ] = kernelParamsPtr->grid.complexGridSize[ZZ];
704 pme_gpu_realloc_and_copy_fract_shifts(pmeGpu);
705 pme_gpu_realloc_and_copy_bspline_values(pmeGpu);
706 pme_gpu_realloc_grids(pmeGpu);
707 pme_gpu_reinit_3dfft(pmeGpu);
710 /* Several GPU functions that refer to the CPU PME data live here.
711 * We would like to keep these away from the GPU-framework specific code for clarity,
712 * as well as compilation issues with MPI.
715 /*! \brief \libinternal
716 * Copies everything useful from the PME CPU to the PME GPU structure.
717 * The goal is to minimize interaction with the PME CPU structure in the GPU code.
719 * \param[in] pme The PME structure.
721 static void pme_gpu_copy_common_data_from(const gmx_pme_t *pme)
723 /* TODO: Consider refactoring the CPU PME code to use the same structure,
724 * so that this function becomes 2 lines */
725 PmeGpu *pmeGpu = pme->gpu;
726 pmeGpu->common->ngrids = pme->ngrids;
727 pmeGpu->common->epsilon_r = pme->epsilon_r;
728 pmeGpu->common->ewaldcoeff_q = pme->ewaldcoeff_q;
729 pmeGpu->common->nk[XX] = pme->nkx;
730 pmeGpu->common->nk[YY] = pme->nky;
731 pmeGpu->common->nk[ZZ] = pme->nkz;
732 pmeGpu->common->pme_order = pme->pme_order;
733 for (int i = 0; i < DIM; i++)
735 pmeGpu->common->bsp_mod[i].assign(pme->bsp_mod[i], pme->bsp_mod[i] + pmeGpu->common->nk[i]);
737 const int cellCount = c_pmeNeighborUnitcellCount;
738 pmeGpu->common->fsh.resize(0);
739 pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshx, pme->fshx + cellCount * pme->nkx);
740 pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshy, pme->fshy + cellCount * pme->nky);
741 pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshz, pme->fshz + cellCount * pme->nkz);
742 pmeGpu->common->nn.resize(0);
743 pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnx, pme->nnx + cellCount * pme->nkx);
744 pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nny, pme->nny + cellCount * pme->nky);
745 pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnz, pme->nnz + cellCount * pme->nkz);
746 pmeGpu->common->runMode = pme->runMode;
747 pmeGpu->common->boxScaler = pme->boxScaler;
750 /*! \libinternal \brief
751 * Initializes the PME GPU data at the beginning of the run.
752 * TODO: this should become PmeGpu::PmeGpu()
754 * \param[in,out] pme The PME structure.
755 * \param[in,out] gpuInfo The GPU information structure.
756 * \param[in] pmeGpuProgram The handle to the program/kernel data created outside (e.g. in unit tests/runner)
758 static void pme_gpu_init(gmx_pme_t *pme,
759 gmx_device_info_t *gpuInfo,
760 PmeGpuProgramHandle pmeGpuProgram)
762 pme->gpu = new PmeGpu();
763 PmeGpu *pmeGpu = pme->gpu;
764 changePinningPolicy(&pmeGpu->staging.h_forces, pme_get_pinning_policy());
765 pmeGpu->common = std::make_shared<PmeShared>();
767 /* These settings are set here for the whole run; dynamic ones are set in pme_gpu_reinit() */
768 /* A convenience variable. */
769 pmeGpu->settings.useDecomposition = (pme->nnodes == 1);
770 /* TODO: CPU gather with GPU spread is broken due to different theta/dtheta layout. */
771 pmeGpu->settings.performGPUGather = true;
773 pme_gpu_set_testing(pmeGpu, false);
775 pmeGpu->deviceInfo = gpuInfo;
776 GMX_ASSERT(pmeGpuProgram != nullptr, "GPU kernels must be already compiled");
777 pmeGpu->programHandle_ = pmeGpuProgram;
779 pme_gpu_init_internal(pmeGpu);
780 pme_gpu_alloc_energy_virial(pmeGpu);
782 pme_gpu_copy_common_data_from(pme);
784 GMX_ASSERT(pmeGpu->common->epsilon_r != 0.0f, "PME GPU: bad electrostatic coefficient");
786 auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
787 kernelParamsPtr->constants.elFactor = ONE_4PI_EPS0 / pmeGpu->common->epsilon_r;
790 void pme_gpu_transform_spline_atom_data(const PmeGpu *pmeGpu, const pme_atomcomm_t *atc,
791 PmeSplineDataType type, int dimIndex, PmeLayoutTransform transform)
793 // The GPU atom spline data is laid out in a different way currently than the CPU one.
794 // This function converts the data from GPU to CPU layout (in the host memory).
795 // It is only intended for testing purposes so far.
796 // Ideally we should use similar layouts on CPU and GPU if we care about mixed modes and their performance
797 // (e.g. spreading on GPU, gathering on CPU).
798 GMX_RELEASE_ASSERT(atc->nthread == 1, "Only the serial PME data layout is supported");
799 const uintmax_t threadIndex = 0;
800 const auto atomCount = pme_gpu_get_kernel_params_base_ptr(pmeGpu)->atoms.nAtoms;
801 const auto atomsPerWarp = pme_gpu_get_atoms_per_warp(pmeGpu);
802 const auto pmeOrder = pmeGpu->common->pme_order;
804 real *cpuSplineBuffer;
805 float *h_splineBuffer;
808 case PmeSplineDataType::Values:
809 cpuSplineBuffer = atc->spline[threadIndex].theta[dimIndex];
810 h_splineBuffer = pmeGpu->staging.h_theta;
813 case PmeSplineDataType::Derivatives:
814 cpuSplineBuffer = atc->spline[threadIndex].dtheta[dimIndex];
815 h_splineBuffer = pmeGpu->staging.h_dtheta;
819 GMX_THROW(gmx::InternalError("Unknown spline data type"));
822 for (auto atomIndex = 0; atomIndex < atomCount; atomIndex++)
824 for (auto orderIndex = 0; orderIndex < pmeOrder; orderIndex++)
826 const auto gpuValueIndex = getSplineParamFullIndex(pmeOrder, orderIndex, dimIndex, atomIndex, atomsPerWarp);
827 const auto cpuValueIndex = atomIndex * pmeOrder + orderIndex;
828 GMX_ASSERT(cpuValueIndex < atomCount * pmeOrder, "Atom spline data index out of bounds (while transforming GPU data layout for host)");
831 case PmeLayoutTransform::GpuToHost:
832 cpuSplineBuffer[cpuValueIndex] = h_splineBuffer[gpuValueIndex];
835 case PmeLayoutTransform::HostToGpu:
836 h_splineBuffer[gpuValueIndex] = cpuSplineBuffer[cpuValueIndex];
840 GMX_THROW(gmx::InternalError("Unknown layout transform"));
846 void pme_gpu_get_real_grid_sizes(const PmeGpu *pmeGpu, gmx::IVec *gridSize, gmx::IVec *paddedGridSize)
848 GMX_ASSERT(gridSize != nullptr, "");
849 GMX_ASSERT(paddedGridSize != nullptr, "");
850 GMX_ASSERT(pmeGpu != nullptr, "");
851 auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
852 for (int i = 0; i < DIM; i++)
854 (*gridSize)[i] = kernelParamsPtr->grid.realGridSize[i];
855 (*paddedGridSize)[i] = kernelParamsPtr->grid.realGridSizePadded[i];
859 void pme_gpu_reinit(gmx_pme_t *pme,
860 gmx_device_info_t *gpuInfo,
861 PmeGpuProgramHandle pmeGpuProgram)
863 if (!pme_gpu_active(pme))
870 /* First-time initialization */
871 pme_gpu_init(pme, gpuInfo, pmeGpuProgram);
875 /* After this call nothing in the GPU code should refer to the gmx_pme_t *pme itself - until the next pme_gpu_reinit */
876 pme_gpu_copy_common_data_from(pme);
878 /* GPU FFT will only get used for a single rank.*/
879 pme->gpu->settings.performGPUFFT = (pme->gpu->common->runMode == PmeRunMode::GPU) && !pme_gpu_uses_dd(pme->gpu);
880 pme->gpu->settings.performGPUSolve = (pme->gpu->common->runMode == PmeRunMode::GPU);
882 /* Reinit active timers */
883 pme_gpu_reinit_timings(pme->gpu);
885 pme_gpu_reinit_grids(pme->gpu);
886 // Note: if timing the reinit launch overhead becomes more relevant
887 // (e.g. with regulat PP-PME re-balancing), we should pass wcycle here.
888 pme_gpu_reinit_computation(pme, nullptr);
889 /* Clear the previous box - doesn't hurt, and forces the PME CPU recipbox
890 * update for mixed mode on grid switch. TODO: use shared recipbox field.
892 std::memset(pme->gpu->common->previousBox, 0, sizeof(pme->gpu->common->previousBox));
895 void pme_gpu_destroy(PmeGpu *pmeGpu)
897 /* Free lots of data */
898 pme_gpu_free_energy_virial(pmeGpu);
899 pme_gpu_free_bspline_values(pmeGpu);
900 pme_gpu_free_forces(pmeGpu);
901 pme_gpu_free_coordinates(pmeGpu);
902 pme_gpu_free_coefficients(pmeGpu);
903 pme_gpu_free_spline_data(pmeGpu);
904 pme_gpu_free_grid_indices(pmeGpu);
905 pme_gpu_free_fract_shifts(pmeGpu);
906 pme_gpu_free_grids(pmeGpu);
908 pme_gpu_destroy_3dfft(pmeGpu);
910 /* Free the GPU-framework specific data last */
911 pme_gpu_destroy_specific(pmeGpu);
916 void pme_gpu_reinit_atoms(PmeGpu *pmeGpu, const int nAtoms, const real *charges)
918 auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
919 kernelParamsPtr->atoms.nAtoms = nAtoms;
920 const int alignment = pme_gpu_get_atom_data_alignment(pmeGpu);
921 pmeGpu->nAtomsPadded = ((nAtoms + alignment - 1) / alignment) * alignment;
922 const int nAtomsAlloc = c_usePadding ? pmeGpu->nAtomsPadded : nAtoms;
923 const bool haveToRealloc = (pmeGpu->nAtomsAlloc < nAtomsAlloc); /* This check might be redundant, but is logical */
924 pmeGpu->nAtomsAlloc = nAtomsAlloc;
927 GMX_RELEASE_ASSERT(false, "Only single precision supported");
928 GMX_UNUSED_VALUE(charges);
930 pme_gpu_realloc_and_copy_input_coefficients(pmeGpu, reinterpret_cast<const float *>(charges));
931 /* Could also be checked for haveToRealloc, but the copy always needs to be performed */
936 pme_gpu_realloc_coordinates(pmeGpu);
937 pme_gpu_realloc_forces(pmeGpu);
938 pme_gpu_realloc_spline_data(pmeGpu);
939 pme_gpu_realloc_grid_indices(pmeGpu);
943 void pme_gpu_3dfft(const PmeGpu *pmeGpu, gmx_fft_direction dir, int grid_index)
945 int timerId = (dir == GMX_FFT_REAL_TO_COMPLEX) ? gtPME_FFT_R2C : gtPME_FFT_C2R;
947 auto *timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timerId);
948 pme_gpu_start_timing(pmeGpu, timerId);
949 pmeGpu->archSpecific->fftSetup[grid_index]->perform3dFft(dir, timingEvent);
950 pme_gpu_stop_timing(pmeGpu, timerId);
954 * Given possibly large \p blockCount, returns a compact 1D or 2D grid for kernel scheduling,
955 * to minimize number of unused blocks.
957 std::pair<int, int> inline pmeGpuCreateGrid(const PmeGpu *pmeGpu, int blockCount)
959 // How many maximum widths in X do we need (hopefully just one)
960 const int minRowCount = (blockCount + pmeGpu->maxGridWidthX - 1) / pmeGpu->maxGridWidthX;
961 // Trying to make things even
962 const int colCount = (blockCount + minRowCount - 1) / minRowCount;
963 GMX_ASSERT((colCount * minRowCount - blockCount) >= 0, "pmeGpuCreateGrid: totally wrong");
964 GMX_ASSERT((colCount * minRowCount - blockCount) < minRowCount, "pmeGpuCreateGrid: excessive blocks");
965 return std::pair<int, int>(colCount, minRowCount);
968 void pme_gpu_spread(const PmeGpu *pmeGpu,
969 int gmx_unused gridIndex,
974 GMX_ASSERT(computeSplines || spreadCharges, "PME spline/spread kernel has invalid input (nothing to do)");
975 const auto *kernelParamsPtr = pmeGpu->kernelParams.get();
976 GMX_ASSERT(kernelParamsPtr->atoms.nAtoms > 0, "No atom data in PME GPU spread");
978 const int maxThreadsPerBlock = c_spreadMaxWarpsPerBlock * pmeGpu->programHandle_->impl_->warpSize;
979 // TODO take device limitations (CL_DEVICE_MAX_WORK_GROUP_SIZE) into account for all 3 kernels
980 const int order = pmeGpu->common->pme_order;
981 const int atomsPerBlock = maxThreadsPerBlock / PME_SPREADGATHER_THREADS_PER_ATOM;
982 // TODO: pick smaller block size in runtime if needed
983 // (e.g. on 660 Ti where 50% occupancy is ~25% faster than 100% occupancy with RNAse (~17.8k atoms))
984 // If doing so, change atomsPerBlock in the kernels as well.
985 // TODO: test varying block sizes on modern arch-s as well
986 // TODO: also consider using cudaFuncSetCacheConfig() for preferring shared memory on older architectures
987 //(for spline data mostly, together with varying PME_GPU_PARALLEL_SPLINE define)
988 GMX_ASSERT(!c_usePadding || !(PME_ATOM_DATA_ALIGNMENT % atomsPerBlock), "inconsistent atom data padding vs. spreading block size");
990 const int blockCount = pmeGpu->nAtomsPadded / atomsPerBlock;
991 auto dimGrid = pmeGpuCreateGrid(pmeGpu, blockCount);
993 KernelLaunchConfig config;
994 config.blockSize[0] = config.blockSize[1] = order;
995 config.blockSize[2] = atomsPerBlock;
996 config.gridSize[0] = dimGrid.first;
997 config.gridSize[1] = dimGrid.second;
998 config.stream = pmeGpu->archSpecific->pmeStream;
1002 GMX_THROW(gmx::NotImplementedError("The code for pme_order != 4 was not implemented!"));
1006 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1011 timingId = gtPME_SPLINEANDSPREAD;
1012 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernel;
1016 timingId = gtPME_SPLINE;
1017 kernelPtr = pmeGpu->programHandle_->impl_->splineKernel;
1022 timingId = gtPME_SPREAD;
1023 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernel;
1026 pme_gpu_start_timing(pmeGpu, timingId);
1027 auto *timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1028 const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1029 launchGpuKernel(kernelPtr, config, timingEvent, "PME spline/spread", kernelArgs);
1030 pme_gpu_stop_timing(pmeGpu, timingId);
1032 const bool copyBackGrid = spreadCharges && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_FFT(pmeGpu));
1035 pme_gpu_copy_output_spread_grid(pmeGpu, h_grid);
1037 const bool copyBackAtomData = computeSplines && (pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_gather(pmeGpu));
1038 if (copyBackAtomData)
1040 pme_gpu_copy_output_spread_atom_data(pmeGpu);
1044 void pme_gpu_solve(const PmeGpu *pmeGpu, t_complex *h_grid,
1045 GridOrdering gridOrdering, bool computeEnergyAndVirial)
1047 const bool copyInputAndOutputGrid = pme_gpu_is_testing(pmeGpu) || !pme_gpu_performs_FFT(pmeGpu);
1049 auto *kernelParamsPtr = pmeGpu->kernelParams.get();
1051 float *h_gridFloat = reinterpret_cast<float *>(h_grid);
1052 if (copyInputAndOutputGrid)
1054 copyToDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid, h_gridFloat,
1055 0, pmeGpu->archSpecific->complexGridSize,
1056 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
1059 int majorDim = -1, middleDim = -1, minorDim = -1;
1060 switch (gridOrdering)
1062 case GridOrdering::YZX:
1068 case GridOrdering::XYZ:
1075 GMX_ASSERT(false, "Implement grid ordering here and below for the kernel launch");
1078 const int maxThreadsPerBlock = c_solveMaxWarpsPerBlock * pmeGpu->programHandle_->impl_->warpSize;
1079 const int maxBlockSize = maxThreadsPerBlock;
1080 const int gridLineSize = pmeGpu->kernelParams->grid.complexGridSize[minorDim];
1081 const int gridLinesPerBlock = std::max(maxBlockSize / gridLineSize, 1);
1082 const int blocksPerGridLine = (gridLineSize + maxBlockSize - 1) / maxBlockSize;
1083 const int cellsPerBlock = gridLineSize * gridLinesPerBlock;
1084 const size_t warpSize = pmeGpu->programHandle_->impl_->warpSize;
1085 const int blockSize = (cellsPerBlock + warpSize - 1) / warpSize * warpSize;
1088 KernelLaunchConfig config;
1089 config.blockSize[0] = blockSize;
1090 config.gridSize[0] = blocksPerGridLine;
1091 // rounding up to full warps so that shuffle operations produce defined results
1092 config.gridSize[1] = (pmeGpu->kernelParams->grid.complexGridSize[middleDim] + gridLinesPerBlock - 1) / gridLinesPerBlock;
1093 config.gridSize[2] = pmeGpu->kernelParams->grid.complexGridSize[majorDim];
1094 config.stream = pmeGpu->archSpecific->pmeStream;
1096 int timingId = gtPME_SOLVE;
1097 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1098 if (gridOrdering == GridOrdering::YZX)
1100 kernelPtr = computeEnergyAndVirial ?
1101 pmeGpu->programHandle_->impl_->solveYZXEnergyKernel :
1102 pmeGpu->programHandle_->impl_->solveYZXKernel;
1104 else if (gridOrdering == GridOrdering::XYZ)
1106 kernelPtr = computeEnergyAndVirial ?
1107 pmeGpu->programHandle_->impl_->solveXYZEnergyKernel :
1108 pmeGpu->programHandle_->impl_->solveXYZKernel;
1111 pme_gpu_start_timing(pmeGpu, timingId);
1112 auto *timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1113 const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1114 launchGpuKernel(kernelPtr, config, timingEvent, "PME solve", kernelArgs);
1115 pme_gpu_stop_timing(pmeGpu, timingId);
1117 if (computeEnergyAndVirial)
1119 copyFromDeviceBuffer(pmeGpu->staging.h_virialAndEnergy, &kernelParamsPtr->constants.d_virialAndEnergy,
1120 0, c_virialAndEnergyCount,
1121 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
1124 if (copyInputAndOutputGrid)
1126 copyFromDeviceBuffer(h_gridFloat, &kernelParamsPtr->grid.d_fourierGrid,
1127 0, pmeGpu->archSpecific->complexGridSize,
1128 pmeGpu->archSpecific->pmeStream, pmeGpu->settings.transferKind, nullptr);
1132 void pme_gpu_gather(PmeGpu *pmeGpu,
1133 PmeForceOutputHandling forceTreatment,
1137 /* Copying the input CPU forces for reduction */
1138 if (forceTreatment != PmeForceOutputHandling::Set)
1140 pme_gpu_copy_input_forces(pmeGpu);
1143 const int order = pmeGpu->common->pme_order;
1144 const auto *kernelParamsPtr = pmeGpu->kernelParams.get();
1146 if (!pme_gpu_performs_FFT(pmeGpu) || pme_gpu_is_testing(pmeGpu))
1148 pme_gpu_copy_input_gather_grid(pmeGpu, const_cast<float *>(h_grid));
1151 if (pme_gpu_is_testing(pmeGpu))
1153 pme_gpu_copy_input_gather_atom_data(pmeGpu);
1156 const int maxThreadsPerBlock = c_gatherMaxWarpsPerBlock * pmeGpu->programHandle_->impl_->warpSize;
1157 const int atomsPerBlock = maxThreadsPerBlock / PME_SPREADGATHER_THREADS_PER_ATOM;
1158 GMX_ASSERT(!c_usePadding || !(PME_ATOM_DATA_ALIGNMENT % atomsPerBlock), "inconsistent atom data padding vs. gathering block size");
1160 const int blockCount = pmeGpu->nAtomsPadded / atomsPerBlock;
1161 auto dimGrid = pmeGpuCreateGrid(pmeGpu, blockCount);
1163 KernelLaunchConfig config;
1164 config.blockSize[0] = config.blockSize[1] = order;
1165 config.blockSize[2] = atomsPerBlock;
1166 config.gridSize[0] = dimGrid.first;
1167 config.gridSize[1] = dimGrid.second;
1168 config.stream = pmeGpu->archSpecific->pmeStream;
1172 GMX_THROW(gmx::NotImplementedError("The code for pme_order != 4 was not implemented!"));
1175 // TODO test different cache configs
1177 int timingId = gtPME_GATHER;
1178 // TODO design kernel selection getters and make PmeGpu a friend of PmeGpuProgramImpl
1179 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = (forceTreatment == PmeForceOutputHandling::Set) ?
1180 pmeGpu->programHandle_->impl_->gatherKernel :
1181 pmeGpu->programHandle_->impl_->gatherReduceWithInputKernel;
1183 pme_gpu_start_timing(pmeGpu, timingId);
1184 auto *timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1185 const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1186 launchGpuKernel(kernelPtr, config, timingEvent, "PME gather", kernelArgs);
1187 pme_gpu_stop_timing(pmeGpu, timingId);
1189 pme_gpu_copy_output_forces(pmeGpu);