2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019,2020 by the GROMACS development team.
5 * Copyright (c) 2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
39 * \brief This file contains internal function implementations
40 * for performing the PME calculations on GPU.
42 * Note that this file is compiled as regular C++ source in OpenCL builds, but
43 * it is treated as CUDA source in CUDA-enabled GPU builds.
45 * \author Aleksei Iupinov <a.yupinov@gmail.com>
46 * \ingroup module_ewald
51 #include "pme_gpu_internal.h"
59 #include "gromacs/ewald/ewald_utils.h"
60 #include "gromacs/fft/gpu_3dfft.h"
61 #include "gromacs/gpu_utils/device_context.h"
62 #include "gromacs/gpu_utils/device_stream.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/gpu_utils/pmalloc.h"
66 # include "gromacs/gpu_utils/syclutils.h"
68 #include "gromacs/hardware/device_information.h"
69 #include "gromacs/math/invertmatrix.h"
70 #include "gromacs/math/units.h"
71 #include "gromacs/timing/gpu_timing.h"
72 #include "gromacs/utility/exceptions.h"
73 #include "gromacs/utility/fatalerror.h"
74 #include "gromacs/utility/gmxassert.h"
75 #include "gromacs/utility/logger.h"
76 #include "gromacs/utility/stringutil.h"
77 #include "gromacs/ewald/pme.h"
83 #include "pme_gpu_calculate_splines.h"
84 #include "pme_gpu_constants.h"
85 #include "pme_gpu_program_impl.h"
86 #include "pme_gpu_timings.h"
87 #include "pme_gpu_types.h"
88 #include "pme_gpu_types_host.h"
89 #include "pme_gpu_types_host_impl.h"
91 #include "pme_internal.h"
92 #include "pme_solve.h"
96 * Atom limit above which it is advantageous to turn on the
97 * recalculating of the splines in the gather and using less threads per atom in the spline and spread
99 constexpr int c_pmeGpuPerformanceAtomLimit = 23000;
102 * Wrapper for getting a pointer to the plain C++ part of the GPU kernel parameters structure.
104 * \param[in] pmeGpu The PME GPU structure.
105 * \returns The pointer to the kernel parameters.
107 static PmeGpuKernelParamsBase* pme_gpu_get_kernel_params_base_ptr(const PmeGpu* pmeGpu)
109 // reinterpret_cast is needed because the derived CUDA structure is not known in this file
110 auto* kernelParamsPtr = reinterpret_cast<PmeGpuKernelParamsBase*>(pmeGpu->kernelParams.get());
111 return kernelParamsPtr;
115 * Atom data block size (in terms of number of atoms).
116 * This is the least common multiple of number of atoms processed by
117 * a single block/workgroup of the spread and gather kernels.
118 * The GPU atom data buffers must be padded, which means that
119 * the numbers of atoms used for determining the size of the memory
120 * allocation must be divisible by this.
122 constexpr int c_pmeAtomDataBlockSize = 64;
124 int pme_gpu_get_atom_data_block_size()
126 return c_pmeAtomDataBlockSize;
129 void pme_gpu_synchronize(const PmeGpu* pmeGpu)
131 pmeGpu->archSpecific->pmeStream_.synchronize();
134 void pme_gpu_alloc_energy_virial(PmeGpu* pmeGpu)
136 const size_t energyAndVirialSize = c_virialAndEnergyCount * sizeof(float);
139 pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
140 "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
142 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
144 allocateDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy[gridIndex],
145 c_virialAndEnergyCount,
146 pmeGpu->archSpecific->deviceContext_);
147 pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_virialAndEnergy[gridIndex]), energyAndVirialSize);
151 void pme_gpu_free_energy_virial(PmeGpu* pmeGpu)
153 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
155 freeDeviceBuffer(&pmeGpu->kernelParams->constants.d_virialAndEnergy[gridIndex]);
156 pfree(pmeGpu->staging.h_virialAndEnergy[gridIndex]);
157 pmeGpu->staging.h_virialAndEnergy[gridIndex] = nullptr;
161 void pme_gpu_clear_energy_virial(const PmeGpu* pmeGpu)
163 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
165 clearDeviceBufferAsync(&pmeGpu->kernelParams->constants.d_virialAndEnergy[gridIndex],
167 c_virialAndEnergyCount,
168 pmeGpu->archSpecific->pmeStream_);
172 void pme_gpu_realloc_and_copy_bspline_values(PmeGpu* pmeGpu, const int gridIndex)
175 pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
176 "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
177 GMX_ASSERT(gridIndex < pmeGpu->common->ngrids,
178 "Invalid combination of gridIndex and number of grids");
180 const int splineValuesOffset[DIM] = { 0,
181 pmeGpu->kernelParams->grid.realGridSize[XX],
182 pmeGpu->kernelParams->grid.realGridSize[XX]
183 + pmeGpu->kernelParams->grid.realGridSize[YY] };
184 memcpy(&pmeGpu->kernelParams->grid.splineValuesOffset, &splineValuesOffset, sizeof(splineValuesOffset));
186 const int newSplineValuesSize = pmeGpu->kernelParams->grid.realGridSize[XX]
187 + pmeGpu->kernelParams->grid.realGridSize[YY]
188 + pmeGpu->kernelParams->grid.realGridSize[ZZ];
189 const bool shouldRealloc = (newSplineValuesSize > pmeGpu->archSpecific->splineValuesSize[gridIndex]);
190 reallocateDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli[gridIndex],
192 &pmeGpu->archSpecific->splineValuesSize[gridIndex],
193 &pmeGpu->archSpecific->splineValuesCapacity[gridIndex],
194 pmeGpu->archSpecific->deviceContext_);
197 /* Reallocate the host buffer */
198 pfree(pmeGpu->staging.h_splineModuli[gridIndex]);
199 pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_splineModuli[gridIndex]),
200 newSplineValuesSize * sizeof(float));
202 for (int i = 0; i < DIM; i++)
204 memcpy(pmeGpu->staging.h_splineModuli[gridIndex] + splineValuesOffset[i],
205 pmeGpu->common->bsp_mod[i].data(),
206 pmeGpu->common->bsp_mod[i].size() * sizeof(float));
208 /* TODO: pin original buffer instead! */
209 copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli[gridIndex],
210 pmeGpu->staging.h_splineModuli[gridIndex],
213 pmeGpu->archSpecific->pmeStream_,
214 pmeGpu->settings.transferKind,
218 void pme_gpu_free_bspline_values(const PmeGpu* pmeGpu)
220 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
222 pfree(pmeGpu->staging.h_splineModuli[gridIndex]);
223 freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_splineModuli[gridIndex]);
227 void pme_gpu_realloc_forces(PmeGpu* pmeGpu)
229 const size_t newForcesSize = pmeGpu->nAtomsAlloc;
230 GMX_ASSERT(newForcesSize > 0, "Bad number of atoms in PME GPU");
231 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces,
233 &pmeGpu->archSpecific->forcesSize,
234 &pmeGpu->archSpecific->forcesSizeAlloc,
235 pmeGpu->archSpecific->deviceContext_);
236 pmeGpu->staging.h_forces.reserveWithPadding(pmeGpu->nAtomsAlloc);
237 pmeGpu->staging.h_forces.resizeWithPadding(pmeGpu->kernelParams->atoms.nAtoms);
240 void pme_gpu_free_forces(const PmeGpu* pmeGpu)
242 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces);
245 void pme_gpu_copy_input_forces(PmeGpu* pmeGpu)
247 GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
248 copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_forces,
249 pmeGpu->staging.h_forces.data(),
251 pmeGpu->kernelParams->atoms.nAtoms,
252 pmeGpu->archSpecific->pmeStream_,
253 pmeGpu->settings.transferKind,
257 void pme_gpu_copy_output_forces(PmeGpu* pmeGpu)
259 GMX_ASSERT(pmeGpu->kernelParams->atoms.nAtoms > 0, "Bad number of atoms in PME GPU");
260 copyFromDeviceBuffer(pmeGpu->staging.h_forces.data(),
261 &pmeGpu->kernelParams->atoms.d_forces,
263 pmeGpu->kernelParams->atoms.nAtoms,
264 pmeGpu->archSpecific->pmeStream_,
265 pmeGpu->settings.transferKind,
269 void pme_gpu_realloc_and_copy_input_coefficients(const PmeGpu* pmeGpu,
270 const float* h_coefficients,
273 GMX_ASSERT(h_coefficients, "Bad host-side charge buffer in PME GPU");
274 const size_t newCoefficientsSize = pmeGpu->nAtomsAlloc;
275 GMX_ASSERT(newCoefficientsSize > 0, "Bad number of atoms in PME GPU");
276 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex],
278 &pmeGpu->archSpecific->coefficientsSize[gridIndex],
279 &pmeGpu->archSpecific->coefficientsCapacity[gridIndex],
280 pmeGpu->archSpecific->deviceContext_);
281 copyToDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex],
282 const_cast<float*>(h_coefficients),
284 pmeGpu->kernelParams->atoms.nAtoms,
285 pmeGpu->archSpecific->pmeStream_,
286 pmeGpu->settings.transferKind,
289 const size_t paddingIndex = pmeGpu->kernelParams->atoms.nAtoms;
290 const size_t paddingCount = pmeGpu->nAtomsAlloc - paddingIndex;
291 if (paddingCount > 0)
293 clearDeviceBufferAsync(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex],
296 pmeGpu->archSpecific->pmeStream_);
300 void pme_gpu_free_coefficients(const PmeGpu* pmeGpu)
302 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
304 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_coefficients[gridIndex]);
308 void pme_gpu_realloc_spline_data(PmeGpu* pmeGpu)
310 const int order = pmeGpu->common->pme_order;
311 const int newSplineDataSize = DIM * order * pmeGpu->nAtomsAlloc;
312 GMX_ASSERT(newSplineDataSize > 0, "Bad number of atoms in PME GPU");
313 /* Two arrays of the same size */
314 const bool shouldRealloc = (newSplineDataSize > pmeGpu->archSpecific->splineDataSize);
315 int currentSizeTemp = pmeGpu->archSpecific->splineDataSize;
316 int currentSizeTempAlloc = pmeGpu->archSpecific->splineDataSizeAlloc;
317 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta,
320 ¤tSizeTempAlloc,
321 pmeGpu->archSpecific->deviceContext_);
322 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta,
324 &pmeGpu->archSpecific->splineDataSize,
325 &pmeGpu->archSpecific->splineDataSizeAlloc,
326 pmeGpu->archSpecific->deviceContext_);
327 // the host side reallocation
330 pfree(pmeGpu->staging.h_theta);
331 pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_theta), newSplineDataSize * sizeof(float));
332 pfree(pmeGpu->staging.h_dtheta);
333 pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_dtheta), newSplineDataSize * sizeof(float));
337 void pme_gpu_free_spline_data(const PmeGpu* pmeGpu)
339 /* Two arrays of the same size */
340 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_theta);
341 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_dtheta);
342 pfree(pmeGpu->staging.h_theta);
343 pfree(pmeGpu->staging.h_dtheta);
346 void pme_gpu_realloc_grid_indices(PmeGpu* pmeGpu)
348 const size_t newIndicesSize = DIM * pmeGpu->nAtomsAlloc;
349 GMX_ASSERT(newIndicesSize > 0, "Bad number of atoms in PME GPU");
350 reallocateDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices,
352 &pmeGpu->archSpecific->gridlineIndicesSize,
353 &pmeGpu->archSpecific->gridlineIndicesSizeAlloc,
354 pmeGpu->archSpecific->deviceContext_);
355 pfree(pmeGpu->staging.h_gridlineIndices);
356 pmalloc(reinterpret_cast<void**>(&pmeGpu->staging.h_gridlineIndices), newIndicesSize * sizeof(int));
359 void pme_gpu_free_grid_indices(const PmeGpu* pmeGpu)
361 freeDeviceBuffer(&pmeGpu->kernelParams->atoms.d_gridlineIndices);
362 pfree(pmeGpu->staging.h_gridlineIndices);
365 void pme_gpu_realloc_grids(PmeGpu* pmeGpu)
367 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
369 const int newRealGridSize = kernelParamsPtr->grid.realGridSizePadded[XX]
370 * kernelParamsPtr->grid.realGridSizePadded[YY]
371 * kernelParamsPtr->grid.realGridSizePadded[ZZ];
372 const int newComplexGridSize = kernelParamsPtr->grid.complexGridSizePadded[XX]
373 * kernelParamsPtr->grid.complexGridSizePadded[YY]
374 * kernelParamsPtr->grid.complexGridSizePadded[ZZ] * 2;
375 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
377 // Multiplied by 2 because we count complex grid size for complex numbers, but all allocations/pointers are float
378 if (pmeGpu->archSpecific->performOutOfPlaceFFT)
380 /* 2 separate grids */
381 reallocateDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid[gridIndex],
383 &pmeGpu->archSpecific->complexGridSize[gridIndex],
384 &pmeGpu->archSpecific->complexGridCapacity[gridIndex],
385 pmeGpu->archSpecific->deviceContext_);
386 reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid[gridIndex],
388 &pmeGpu->archSpecific->realGridSize[gridIndex],
389 &pmeGpu->archSpecific->realGridCapacity[gridIndex],
390 pmeGpu->archSpecific->deviceContext_);
394 /* A single buffer so that any grid will fit */
395 const int newGridsSize = std::max(newRealGridSize, newComplexGridSize);
396 reallocateDeviceBuffer(&kernelParamsPtr->grid.d_realGrid[gridIndex],
398 &pmeGpu->archSpecific->realGridSize[gridIndex],
399 &pmeGpu->archSpecific->realGridCapacity[gridIndex],
400 pmeGpu->archSpecific->deviceContext_);
401 kernelParamsPtr->grid.d_fourierGrid[gridIndex] = kernelParamsPtr->grid.d_realGrid[gridIndex];
402 pmeGpu->archSpecific->complexGridSize[gridIndex] =
403 pmeGpu->archSpecific->realGridSize[gridIndex];
404 // the size might get used later for copying the grid
409 void pme_gpu_free_grids(const PmeGpu* pmeGpu)
411 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
413 if (pmeGpu->archSpecific->performOutOfPlaceFFT)
415 freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_fourierGrid[gridIndex]);
417 freeDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid[gridIndex]);
421 void pme_gpu_clear_grids(const PmeGpu* pmeGpu)
423 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
425 clearDeviceBufferAsync(&pmeGpu->kernelParams->grid.d_realGrid[gridIndex],
427 pmeGpu->archSpecific->realGridSize[gridIndex],
428 pmeGpu->archSpecific->pmeStream_);
432 void pme_gpu_realloc_and_copy_fract_shifts(PmeGpu* pmeGpu)
434 pme_gpu_free_fract_shifts(pmeGpu);
436 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
438 const int nx = kernelParamsPtr->grid.realGridSize[XX];
439 const int ny = kernelParamsPtr->grid.realGridSize[YY];
440 const int nz = kernelParamsPtr->grid.realGridSize[ZZ];
441 const int cellCount = c_pmeNeighborUnitcellCount;
442 const int gridDataOffset[DIM] = { 0, cellCount * nx, cellCount * (nx + ny) };
444 memcpy(kernelParamsPtr->grid.tablesOffsets, &gridDataOffset, sizeof(gridDataOffset));
446 const int newFractShiftsSize = cellCount * (nx + ny + nz);
448 initParamLookupTable(&kernelParamsPtr->grid.d_fractShiftsTable,
449 &kernelParamsPtr->fractShiftsTableTexture,
450 pmeGpu->common->fsh.data(),
452 pmeGpu->archSpecific->deviceContext_);
454 initParamLookupTable(&kernelParamsPtr->grid.d_gridlineIndicesTable,
455 &kernelParamsPtr->gridlineIndicesTableTexture,
456 pmeGpu->common->nn.data(),
458 pmeGpu->archSpecific->deviceContext_);
461 void pme_gpu_free_fract_shifts(const PmeGpu* pmeGpu)
463 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
465 destroyParamLookupTable(&kernelParamsPtr->grid.d_fractShiftsTable,
466 &kernelParamsPtr->fractShiftsTableTexture);
467 destroyParamLookupTable(&kernelParamsPtr->grid.d_gridlineIndicesTable,
468 &kernelParamsPtr->gridlineIndicesTableTexture);
469 #elif GMX_GPU_OPENCL || GMX_GPU_SYCL
470 freeDeviceBuffer(&kernelParamsPtr->grid.d_fractShiftsTable);
471 freeDeviceBuffer(&kernelParamsPtr->grid.d_gridlineIndicesTable);
475 bool pme_gpu_stream_query(const PmeGpu* pmeGpu)
477 return haveStreamTasksCompleted(pmeGpu->archSpecific->pmeStream_);
480 void pme_gpu_copy_input_gather_grid(const PmeGpu* pmeGpu, const float* h_grid, const int gridIndex)
482 copyToDeviceBuffer(&pmeGpu->kernelParams->grid.d_realGrid[gridIndex],
485 pmeGpu->archSpecific->realGridSize[gridIndex],
486 pmeGpu->archSpecific->pmeStream_,
487 pmeGpu->settings.transferKind,
491 void pme_gpu_copy_output_spread_grid(const PmeGpu* pmeGpu, float* h_grid, const int gridIndex)
493 copyFromDeviceBuffer(h_grid,
494 &pmeGpu->kernelParams->grid.d_realGrid[gridIndex],
496 pmeGpu->archSpecific->realGridSize[gridIndex],
497 pmeGpu->archSpecific->pmeStream_,
498 pmeGpu->settings.transferKind,
500 pmeGpu->archSpecific->syncSpreadGridD2H.markEvent(pmeGpu->archSpecific->pmeStream_);
503 void pme_gpu_copy_output_spread_atom_data(const PmeGpu* pmeGpu)
505 const size_t splinesCount = DIM * pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order;
506 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
507 copyFromDeviceBuffer(pmeGpu->staging.h_dtheta,
508 &kernelParamsPtr->atoms.d_dtheta,
511 pmeGpu->archSpecific->pmeStream_,
512 pmeGpu->settings.transferKind,
514 copyFromDeviceBuffer(pmeGpu->staging.h_theta,
515 &kernelParamsPtr->atoms.d_theta,
518 pmeGpu->archSpecific->pmeStream_,
519 pmeGpu->settings.transferKind,
521 copyFromDeviceBuffer(pmeGpu->staging.h_gridlineIndices,
522 &kernelParamsPtr->atoms.d_gridlineIndices,
524 kernelParamsPtr->atoms.nAtoms * DIM,
525 pmeGpu->archSpecific->pmeStream_,
526 pmeGpu->settings.transferKind,
530 void pme_gpu_copy_input_gather_atom_data(const PmeGpu* pmeGpu)
532 const size_t splinesCount = DIM * pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order;
533 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
535 // TODO: could clear only the padding and not the whole thing, but this is a test-exclusive code anyway
536 clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_gridlineIndices,
538 pmeGpu->nAtomsAlloc * DIM,
539 pmeGpu->archSpecific->pmeStream_);
540 clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_dtheta,
542 pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM,
543 pmeGpu->archSpecific->pmeStream_);
544 clearDeviceBufferAsync(&kernelParamsPtr->atoms.d_theta,
546 pmeGpu->nAtomsAlloc * pmeGpu->common->pme_order * DIM,
547 pmeGpu->archSpecific->pmeStream_);
549 copyToDeviceBuffer(&kernelParamsPtr->atoms.d_dtheta,
550 pmeGpu->staging.h_dtheta,
553 pmeGpu->archSpecific->pmeStream_,
554 pmeGpu->settings.transferKind,
556 copyToDeviceBuffer(&kernelParamsPtr->atoms.d_theta,
557 pmeGpu->staging.h_theta,
560 pmeGpu->archSpecific->pmeStream_,
561 pmeGpu->settings.transferKind,
563 copyToDeviceBuffer(&kernelParamsPtr->atoms.d_gridlineIndices,
564 pmeGpu->staging.h_gridlineIndices,
566 kernelParamsPtr->atoms.nAtoms * DIM,
567 pmeGpu->archSpecific->pmeStream_,
568 pmeGpu->settings.transferKind,
572 void pme_gpu_sync_spread_grid(const PmeGpu* pmeGpu)
574 pmeGpu->archSpecific->syncSpreadGridD2H.waitForEvent();
577 /*! \brief Internal GPU initialization for PME.
579 * \param[in] pmeGpu GPU PME data.
580 * \param[in] deviceContext GPU context.
581 * \param[in] deviceStream GPU stream.
583 static void pme_gpu_init_internal(PmeGpu* pmeGpu, const DeviceContext& deviceContext, const DeviceStream& deviceStream)
585 /* Allocate the target-specific structures */
586 pmeGpu->archSpecific.reset(new PmeGpuSpecific(deviceContext, deviceStream));
587 pmeGpu->kernelParams.reset(new PmeGpuKernelParams());
589 pmeGpu->archSpecific->performOutOfPlaceFFT = true;
590 /* This should give better performance, according to the cuFFT documentation.
591 * The performance seems to be the same though.
592 * TODO: PME could also try to pick up nice grid sizes (with factors of 2, 3, 5, 7).
596 pmeGpu->maxGridWidthX = deviceContext.deviceInfo().prop.maxGridSize[0];
598 // Use this path for any non-CUDA GPU acceleration
599 // TODO: is there no really global work size limit in OpenCL?
600 pmeGpu->maxGridWidthX = INT32_MAX / 2;
604 void pme_gpu_reinit_3dfft(const PmeGpu* pmeGpu)
606 if (pme_gpu_settings(pmeGpu).performGPUFFT)
608 pmeGpu->archSpecific->fftSetup.resize(0);
609 const bool performOutOfPlaceFFT = pmeGpu->archSpecific->performOutOfPlaceFFT;
610 const bool allocateGrid = false;
611 MPI_Comm comm = MPI_COMM_NULL;
612 std::array<int, 1> gridOffsetsInXForEachRank = { 0 };
613 std::array<int, 1> gridOffsetsInYForEachRank = { 0 };
615 const gmx::FftBackend backend = gmx::FftBackend::Cufft;
617 const gmx::FftBackend backend = gmx::FftBackend::Ocl;
619 # if GMX_SYCL_HIPSYCL
620 const gmx::FftBackend backend = gmx::FftBackend::SyclRocfft;
622 const gmx::FftBackend backend = gmx::FftBackend::Sycl;
625 GMX_RELEASE_ASSERT(false, "Unknown GPU backend");
626 const gmx::FftBackend backend = gmx::FftBackend::Count;
629 PmeGpuGridParams& grid = pme_gpu_get_kernel_params_base_ptr(pmeGpu)->grid;
630 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
632 pmeGpu->archSpecific->fftSetup.push_back(
633 std::make_unique<gmx::Gpu3dFft>(backend,
636 gridOffsetsInXForEachRank,
637 gridOffsetsInYForEachRank,
638 grid.realGridSize[ZZ],
639 performOutOfPlaceFFT,
640 pmeGpu->archSpecific->deviceContext_,
641 pmeGpu->archSpecific->pmeStream_,
643 grid.realGridSizePadded,
644 grid.complexGridSizePadded,
645 &(grid.d_realGrid[gridIndex]),
646 &(grid.d_fourierGrid[gridIndex])));
651 void pme_gpu_destroy_3dfft(const PmeGpu* pmeGpu)
653 pmeGpu->archSpecific->fftSetup.resize(0);
656 void pme_gpu_getEnergyAndVirial(const gmx_pme_t& pme, const float lambda, PmeOutput* output)
658 const PmeGpu* pmeGpu = pme.gpu;
660 GMX_ASSERT(lambda == 1.0 || pmeGpu->common->ngrids == 2,
661 "Invalid combination of lambda and number of grids");
663 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
665 for (int j = 0; j < c_virialAndEnergyCount; j++)
667 GMX_ASSERT(std::isfinite(pmeGpu->staging.h_virialAndEnergy[gridIndex][j]),
668 "PME GPU produces incorrect energy/virial.");
671 for (int dim1 = 0; dim1 < DIM; dim1++)
673 for (int dim2 = 0; dim2 < DIM; dim2++)
675 output->coulombVirial_[dim1][dim2] = 0;
678 output->coulombEnergy_ = 0;
680 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
682 if (pmeGpu->common->ngrids == 2)
684 scale = gridIndex == 0 ? (1.0 - lambda) : lambda;
686 output->coulombVirial_[XX][XX] +=
687 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][0];
688 output->coulombVirial_[YY][YY] +=
689 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][1];
690 output->coulombVirial_[ZZ][ZZ] +=
691 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][2];
692 output->coulombVirial_[XX][YY] +=
693 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][3];
694 output->coulombVirial_[YY][XX] +=
695 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][3];
696 output->coulombVirial_[XX][ZZ] +=
697 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][4];
698 output->coulombVirial_[ZZ][XX] +=
699 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][4];
700 output->coulombVirial_[YY][ZZ] +=
701 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][5];
702 output->coulombVirial_[ZZ][YY] +=
703 scale * 0.25F * pmeGpu->staging.h_virialAndEnergy[gridIndex][5];
704 output->coulombEnergy_ += scale * 0.5F * pmeGpu->staging.h_virialAndEnergy[gridIndex][6];
706 if (pmeGpu->common->ngrids > 1)
708 output->coulombDvdl_ = 0.5F
709 * (pmeGpu->staging.h_virialAndEnergy[FEP_STATE_B][6]
710 - pmeGpu->staging.h_virialAndEnergy[FEP_STATE_A][6]);
714 /*! \brief Sets the force-related members in \p output
716 * \param[in] pmeGpu PME GPU data structure
717 * \param[out] output Pointer to PME output data structure
719 static void pme_gpu_getForceOutput(PmeGpu* pmeGpu, PmeOutput* output)
721 output->haveForceOutput_ = !pmeGpu->settings.useGpuForceReduction;
722 if (output->haveForceOutput_)
724 output->forces_ = pmeGpu->staging.h_forces;
728 PmeOutput pme_gpu_getOutput(const gmx_pme_t& pme, const bool computeEnergyAndVirial, const real lambdaQ)
730 PmeGpu* pmeGpu = pme.gpu;
734 pme_gpu_getForceOutput(pmeGpu, &output);
736 if (computeEnergyAndVirial)
738 if (pme_gpu_settings(pmeGpu).performGPUSolve)
740 pme_gpu_getEnergyAndVirial(pme, lambdaQ, &output);
744 get_pme_ener_vir_q(pme.solve_work, pme.nthread, &output);
750 void pme_gpu_update_input_box(PmeGpu gmx_unused* pmeGpu, const matrix gmx_unused box)
753 GMX_THROW(gmx::NotImplementedError("PME is implemented for single-precision only on GPU"));
756 pmeGpu->common->boxScaler->scaleBox(box, scaledBox);
757 auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
758 kernelParamsPtr->current.boxVolume = scaledBox[XX][XX] * scaledBox[YY][YY] * scaledBox[ZZ][ZZ];
759 GMX_ASSERT(kernelParamsPtr->current.boxVolume != 0.0F, "Zero volume of the unit cell");
761 gmx::invertBoxMatrix(scaledBox, recipBox);
763 /* The GPU recipBox is transposed as compared to the CPU recipBox.
764 * Spread uses matrix columns (while solve and gather use rows).
765 * There is no particular reason for this; it might be further rethought/optimized for better access patterns.
767 const real newRecipBox[DIM][DIM] = { { recipBox[XX][XX], recipBox[YY][XX], recipBox[ZZ][XX] },
768 { 0.0, recipBox[YY][YY], recipBox[ZZ][YY] },
769 { 0.0, 0.0, recipBox[ZZ][ZZ] } };
770 memcpy(kernelParamsPtr->current.recipBox, newRecipBox, sizeof(matrix));
774 /*! \brief \libinternal
775 * (Re-)initializes all the PME GPU data related to the grid size and cut-off.
777 * \param[in] pmeGpu The PME GPU structure.
779 static void pme_gpu_reinit_grids(PmeGpu* pmeGpu)
781 auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
784 pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
785 "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
787 kernelParamsPtr->grid.ewaldFactor =
788 (M_PI * M_PI) / (pmeGpu->common->ewaldcoeff_q * pmeGpu->common->ewaldcoeff_q);
789 /* The grid size variants */
790 for (int i = 0; i < DIM; i++)
792 kernelParamsPtr->grid.realGridSize[i] = pmeGpu->common->nk[i];
793 kernelParamsPtr->grid.realGridSizeFP[i] =
794 static_cast<float>(kernelParamsPtr->grid.realGridSize[i]);
795 kernelParamsPtr->grid.realGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
797 // The complex grid currently uses no padding;
798 // if it starts to do so, then another test should be added for that
799 kernelParamsPtr->grid.complexGridSize[i] = kernelParamsPtr->grid.realGridSize[i];
800 kernelParamsPtr->grid.complexGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
802 /* FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
803 if (!pme_gpu_settings(pmeGpu).performGPUFFT)
805 // This allows for GPU spreading grid and CPU fftgrid to have the same layout, so that we can copy the data directly
806 kernelParamsPtr->grid.realGridSizePadded[ZZ] =
807 (kernelParamsPtr->grid.realGridSize[ZZ] / 2 + 1) * 2;
809 /* GPU FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
810 kernelParamsPtr->grid.complexGridSize[ZZ] /= 2;
811 kernelParamsPtr->grid.complexGridSize[ZZ]++;
812 kernelParamsPtr->grid.complexGridSizePadded[ZZ] = kernelParamsPtr->grid.complexGridSize[ZZ];
814 pme_gpu_realloc_and_copy_fract_shifts(pmeGpu);
815 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
817 pme_gpu_realloc_and_copy_bspline_values(pmeGpu, gridIndex);
820 pme_gpu_realloc_grids(pmeGpu);
821 pme_gpu_reinit_3dfft(pmeGpu);
824 /* Several GPU functions that refer to the CPU PME data live here.
825 * We would like to keep these away from the GPU-framework specific code for clarity,
826 * as well as compilation issues with MPI.
829 /*! \brief \libinternal
830 * Copies everything useful from the PME CPU to the PME GPU structure.
831 * The goal is to minimize interaction with the PME CPU structure in the GPU code.
833 * \param[in] pme The PME structure.
835 static void pme_gpu_copy_common_data_from(const gmx_pme_t* pme)
837 /* TODO: Consider refactoring the CPU PME code to use the same structure,
838 * so that this function becomes 2 lines */
839 PmeGpu* pmeGpu = pme->gpu;
840 pmeGpu->common->ngrids = pme->bFEP_q ? 2 : 1;
841 pmeGpu->common->epsilon_r = pme->epsilon_r;
842 pmeGpu->common->ewaldcoeff_q = pme->ewaldcoeff_q;
843 pmeGpu->common->nk[XX] = pme->nkx;
844 pmeGpu->common->nk[YY] = pme->nky;
845 pmeGpu->common->nk[ZZ] = pme->nkz;
846 pmeGpu->common->pme_order = pme->pme_order;
847 if (pmeGpu->common->pme_order != c_pmeGpuOrder)
849 GMX_THROW(gmx::NotImplementedError("pme_order != 4 is not implemented!"));
851 for (int i = 0; i < DIM; i++)
853 pmeGpu->common->bsp_mod[i].assign(pme->bsp_mod[i], pme->bsp_mod[i] + pmeGpu->common->nk[i]);
855 const int cellCount = c_pmeNeighborUnitcellCount;
856 pmeGpu->common->fsh.resize(0);
857 pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshx, pme->fshx + cellCount * pme->nkx);
858 pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshy, pme->fshy + cellCount * pme->nky);
859 pmeGpu->common->fsh.insert(pmeGpu->common->fsh.end(), pme->fshz, pme->fshz + cellCount * pme->nkz);
860 pmeGpu->common->nn.resize(0);
861 pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnx, pme->nnx + cellCount * pme->nkx);
862 pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nny, pme->nny + cellCount * pme->nky);
863 pmeGpu->common->nn.insert(pmeGpu->common->nn.end(), pme->nnz, pme->nnz + cellCount * pme->nkz);
864 pmeGpu->common->runMode = pme->runMode;
865 pmeGpu->common->isRankPmeOnly = !pme->bPPnode;
866 pmeGpu->common->boxScaler = pme->boxScaler.get();
869 /*! \libinternal \brief
870 * uses heuristics to select the best performing PME gather and scatter kernels
872 * \param[in,out] pmeGpu The PME GPU structure.
874 static void pme_gpu_select_best_performing_pme_spreadgather_kernels(PmeGpu* pmeGpu)
876 if (GMX_GPU_CUDA && pmeGpu->kernelParams->atoms.nAtoms > c_pmeGpuPerformanceAtomLimit)
878 pmeGpu->settings.threadsPerAtom = ThreadsPerAtom::Order;
879 pmeGpu->settings.recalculateSplines = true;
883 pmeGpu->settings.threadsPerAtom = ThreadsPerAtom::OrderSquared;
884 pmeGpu->settings.recalculateSplines = false;
889 /*! \libinternal \brief
890 * Initializes the PME GPU data at the beginning of the run.
891 * TODO: this should become PmeGpu::PmeGpu()
893 * \param[in,out] pme The PME structure.
894 * \param[in] deviceContext The GPU context.
895 * \param[in] deviceStream The GPU stream.
896 * \param[in,out] pmeGpuProgram The handle to the program/kernel data created outside (e.g. in unit tests/runner)
898 static void pme_gpu_init(gmx_pme_t* pme,
899 const DeviceContext& deviceContext,
900 const DeviceStream& deviceStream,
901 const PmeGpuProgram* pmeGpuProgram)
903 pme->gpu = new PmeGpu();
904 PmeGpu* pmeGpu = pme->gpu;
905 changePinningPolicy(&pmeGpu->staging.h_forces, pme_get_pinning_policy());
906 pmeGpu->common = std::make_shared<PmeShared>();
908 /* These settings are set here for the whole run; dynamic ones are set in pme_gpu_reinit() */
909 /* A convenience variable. */
910 pmeGpu->settings.useDecomposition = (pme->nnodes != 1);
911 /* TODO: CPU gather with GPU spread is broken due to different theta/dtheta layout. */
912 pmeGpu->settings.performGPUGather = true;
913 // By default GPU-side reduction is off (explicitly set here for tests, otherwise reset per-step)
914 pmeGpu->settings.useGpuForceReduction = false;
916 pme_gpu_set_testing(pmeGpu, false);
918 GMX_ASSERT(pmeGpuProgram != nullptr, "GPU kernels must be already compiled");
919 pmeGpu->programHandle_ = pmeGpuProgram;
921 pmeGpu->initializedClfftLibrary_ = std::make_unique<gmx::ClfftInitializer>();
923 pme_gpu_init_internal(pmeGpu, deviceContext, deviceStream);
925 pme_gpu_copy_common_data_from(pme);
926 pme_gpu_alloc_energy_virial(pmeGpu);
928 GMX_ASSERT(pmeGpu->common->epsilon_r != 0.0F, "PME GPU: bad electrostatic coefficient");
930 auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
931 kernelParamsPtr->constants.elFactor = gmx::c_one4PiEps0 / pmeGpu->common->epsilon_r;
934 void pme_gpu_get_real_grid_sizes(const PmeGpu* pmeGpu, gmx::IVec* gridSize, gmx::IVec* paddedGridSize)
936 GMX_ASSERT(gridSize != nullptr, "");
937 GMX_ASSERT(paddedGridSize != nullptr, "");
938 GMX_ASSERT(pmeGpu != nullptr, "");
939 auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
940 for (int i = 0; i < DIM; i++)
942 (*gridSize)[i] = kernelParamsPtr->grid.realGridSize[i];
943 (*paddedGridSize)[i] = kernelParamsPtr->grid.realGridSizePadded[i];
947 void pme_gpu_reinit(gmx_pme_t* pme,
948 const DeviceContext* deviceContext,
949 const DeviceStream* deviceStream,
950 const PmeGpuProgram* pmeGpuProgram)
952 GMX_ASSERT(pme != nullptr, "Need valid PME object");
956 GMX_RELEASE_ASSERT(deviceContext != nullptr,
957 "Device context can not be nullptr when setting up PME on GPU.");
958 GMX_RELEASE_ASSERT(deviceStream != nullptr,
959 "Device stream can not be nullptr when setting up PME on GPU.");
960 /* First-time initialization */
961 pme_gpu_init(pme, *deviceContext, *deviceStream, pmeGpuProgram);
965 /* After this call nothing in the GPU code should refer to the gmx_pme_t *pme itself - until the next pme_gpu_reinit */
966 pme_gpu_copy_common_data_from(pme);
968 /* GPU FFT will only get used for a single rank.*/
969 pme->gpu->settings.performGPUFFT =
970 (pme->gpu->common->runMode == PmeRunMode::GPU) && !pme->gpu->settings.useDecomposition;
971 pme->gpu->settings.performGPUSolve = (pme->gpu->common->runMode == PmeRunMode::GPU);
973 /* Reinit active timers */
974 pme_gpu_reinit_timings(pme->gpu);
976 pme_gpu_reinit_grids(pme->gpu);
977 // Note: if timing the reinit launch overhead becomes more relevant
978 // (e.g. with regulat PP-PME re-balancing), we should pass wcycle here.
979 pme_gpu_reinit_computation(pme, nullptr);
980 /* Clear the previous box - doesn't hurt, and forces the PME CPU recipbox
981 * update for mixed mode on grid switch. TODO: use shared recipbox field.
983 std::memset(pme->gpu->common->previousBox, 0, sizeof(pme->gpu->common->previousBox));
986 void pme_gpu_destroy(PmeGpu* pmeGpu)
988 /* Free lots of data */
989 pme_gpu_free_energy_virial(pmeGpu);
990 pme_gpu_free_bspline_values(pmeGpu);
991 pme_gpu_free_forces(pmeGpu);
992 pme_gpu_free_coefficients(pmeGpu);
993 pme_gpu_free_spline_data(pmeGpu);
994 pme_gpu_free_grid_indices(pmeGpu);
995 pme_gpu_free_fract_shifts(pmeGpu);
996 pme_gpu_free_grids(pmeGpu);
998 pme_gpu_destroy_3dfft(pmeGpu);
1003 void pme_gpu_reinit_atoms(PmeGpu* pmeGpu, const int nAtoms, const real* chargesA, const real* chargesB)
1005 auto* kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGpu);
1006 kernelParamsPtr->atoms.nAtoms = nAtoms;
1007 const int block_size = pme_gpu_get_atom_data_block_size();
1008 const int nAtomsNewPadded = ((nAtoms + block_size - 1) / block_size) * block_size;
1009 const bool haveToRealloc = (pmeGpu->nAtomsAlloc < nAtomsNewPadded);
1010 pmeGpu->nAtomsAlloc = nAtomsNewPadded;
1013 GMX_RELEASE_ASSERT(false, "Only single precision supported");
1014 GMX_UNUSED_VALUE(charges);
1017 /* Could also be checked for haveToRealloc, but the copy always needs to be performed */
1018 pme_gpu_realloc_and_copy_input_coefficients(pmeGpu, reinterpret_cast<const float*>(chargesA), gridIndex);
1020 if (chargesB != nullptr)
1022 pme_gpu_realloc_and_copy_input_coefficients(
1023 pmeGpu, reinterpret_cast<const float*>(chargesB), gridIndex);
1027 /* Fill the second set of coefficients with chargesA as well to be able to avoid
1028 * conditionals in the GPU kernels */
1029 /* FIXME: This should be avoided by making a separate templated version of the
1030 * relevant kernel(s) (probably only pme_gather_kernel). That would require a
1031 * reduction of the current number of templated parameters of that kernel. */
1032 pme_gpu_realloc_and_copy_input_coefficients(
1033 pmeGpu, reinterpret_cast<const float*>(chargesA), gridIndex);
1039 pme_gpu_realloc_forces(pmeGpu);
1040 pme_gpu_realloc_spline_data(pmeGpu);
1041 pme_gpu_realloc_grid_indices(pmeGpu);
1043 pme_gpu_select_best_performing_pme_spreadgather_kernels(pmeGpu);
1046 /*! \internal \brief
1047 * Returns raw timing event from the corresponding GpuRegionTimer (if timings are enabled).
1048 * In CUDA result can be nullptr stub, per GpuRegionTimer implementation.
1050 * \param[in] pmeGpu The PME GPU data structure.
1051 * \param[in] pmeStageId The PME GPU stage gtPME_ index from the enum in src/gromacs/timing/gpu_timing.h
1053 static CommandEvent* pme_gpu_fetch_timing_event(const PmeGpu* pmeGpu, PmeStage pmeStageId)
1055 CommandEvent* timingEvent = nullptr;
1056 if (pme_gpu_timings_enabled(pmeGpu))
1058 GMX_ASSERT(pmeStageId < PmeStage::Count, "Wrong PME GPU timing event index");
1059 timingEvent = pmeGpu->archSpecific->timingEvents[pmeStageId].fetchNextEvent();
1064 void pme_gpu_3dfft(const PmeGpu* pmeGpu, gmx_fft_direction dir, const int grid_index)
1066 PmeStage timerId = (dir == GMX_FFT_REAL_TO_COMPLEX) ? PmeStage::FftTransformR2C
1067 : PmeStage::FftTransformC2R;
1069 pme_gpu_start_timing(pmeGpu, timerId);
1070 pmeGpu->archSpecific->fftSetup[grid_index]->perform3dFft(
1071 dir, pme_gpu_fetch_timing_event(pmeGpu, timerId));
1072 pme_gpu_stop_timing(pmeGpu, timerId);
1076 * Given possibly large \p blockCount, returns a compact 1D or 2D grid for kernel scheduling,
1077 * to minimize number of unused blocks.
1079 std::pair<int, int> inline pmeGpuCreateGrid(const PmeGpu* pmeGpu, int blockCount)
1081 // How many maximum widths in X do we need (hopefully just one)
1082 const int minRowCount = (blockCount + pmeGpu->maxGridWidthX - 1) / pmeGpu->maxGridWidthX;
1083 // Trying to make things even
1084 const int colCount = (blockCount + minRowCount - 1) / minRowCount;
1085 GMX_ASSERT((colCount * minRowCount - blockCount) >= 0, "pmeGpuCreateGrid: totally wrong");
1086 GMX_ASSERT((colCount * minRowCount - blockCount) < minRowCount,
1087 "pmeGpuCreateGrid: excessive blocks");
1088 return std::pair<int, int>(colCount, minRowCount);
1092 * Returns a pointer to appropriate spline and spread kernel based on the input bool values
1094 * \param[in] pmeGpu The PME GPU structure.
1095 * \param[in] threadsPerAtom Controls whether we should use order or order*order threads per atom
1096 * \param[in] writeSplinesToGlobal bool controlling if we should write spline data to global memory
1097 * \param[in] numGrids Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1099 * \return Pointer to CUDA kernel
1101 static auto selectSplineAndSpreadKernelPtr(const PmeGpu* pmeGpu,
1102 ThreadsPerAtom threadsPerAtom,
1103 bool writeSplinesToGlobal,
1106 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1107 if (writeSplinesToGlobal)
1109 if (threadsPerAtom == ThreadsPerAtom::Order)
1113 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesThPerAtom4Dual;
1117 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesThPerAtom4Single;
1124 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesDual;
1128 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelWriteSplinesSingle;
1134 if (threadsPerAtom == ThreadsPerAtom::Order)
1138 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Dual;
1142 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Single;
1149 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelDual;
1153 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelSingle;
1162 * Returns a pointer to appropriate spline kernel based on the input bool values
1164 * \param[in] pmeGpu The PME GPU structure.
1165 * \param[in] threadsPerAtom Controls whether we should use order or order*order threads per atom
1166 * \param[in] writeSplinesToGlobal bool controlling if we should write spline data to global memory
1167 * \param[in] numGrids Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1169 * \return Pointer to CUDA kernel
1171 static auto selectSplineKernelPtr(const PmeGpu* pmeGpu,
1172 ThreadsPerAtom threadsPerAtom,
1173 bool gmx_unused writeSplinesToGlobal,
1176 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1178 writeSplinesToGlobal,
1179 "Spline data should always be written to global memory when just calculating splines");
1181 if (threadsPerAtom == ThreadsPerAtom::Order)
1185 kernelPtr = pmeGpu->programHandle_->impl_->splineKernelThPerAtom4Dual;
1189 kernelPtr = pmeGpu->programHandle_->impl_->splineKernelThPerAtom4Single;
1196 kernelPtr = pmeGpu->programHandle_->impl_->splineKernelDual;
1200 kernelPtr = pmeGpu->programHandle_->impl_->splineKernelSingle;
1207 * Returns a pointer to appropriate spread kernel based on the input bool values
1209 * \param[in] pmeGpu The PME GPU structure.
1210 * \param[in] threadsPerAtom Controls whether we should use order or order*order threads per atom
1211 * \param[in] writeSplinesToGlobal bool controlling if we should write spline data to global memory
1212 * \param[in] numGrids Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1214 * \return Pointer to CUDA kernel
1216 static auto selectSpreadKernelPtr(const PmeGpu* pmeGpu,
1217 ThreadsPerAtom threadsPerAtom,
1218 bool writeSplinesToGlobal,
1221 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1222 if (writeSplinesToGlobal)
1224 if (threadsPerAtom == ThreadsPerAtom::Order)
1228 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelThPerAtom4Dual;
1232 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelThPerAtom4Single;
1239 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelDual;
1243 kernelPtr = pmeGpu->programHandle_->impl_->spreadKernelSingle;
1249 /* if we are not saving the spline data we need to recalculate it
1250 using the spline and spread Kernel */
1251 if (threadsPerAtom == ThreadsPerAtom::Order)
1255 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Dual;
1259 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelThPerAtom4Single;
1266 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelDual;
1270 kernelPtr = pmeGpu->programHandle_->impl_->splineAndSpreadKernelSingle;
1277 void pme_gpu_spread(const PmeGpu* pmeGpu,
1278 GpuEventSynchronizer* xReadyOnDevice,
1280 bool computeSplines,
1285 pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
1286 "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
1288 GMX_ASSERT(computeSplines || spreadCharges,
1289 "PME spline/spread kernel has invalid input (nothing to do)");
1290 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1291 GMX_ASSERT(kernelParamsPtr->atoms.nAtoms > 0, "No atom data in PME GPU spread");
1293 const size_t blockSize = pmeGpu->programHandle_->impl_->spreadWorkGroupSize;
1295 const int order = pmeGpu->common->pme_order;
1296 GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1297 const bool writeGlobal = pmeGpu->settings.copyAllOutputs;
1298 const int threadsPerAtom =
1299 (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
1300 const bool recalculateSplines = pmeGpu->settings.recalculateSplines;
1302 GMX_ASSERT(!GMX_GPU_OPENCL || pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::OrderSquared,
1303 "Only 16 threads per atom supported in OpenCL");
1304 GMX_ASSERT(!GMX_GPU_OPENCL || !recalculateSplines,
1305 "Recalculating splines not supported in OpenCL");
1307 const int atomsPerBlock = blockSize / threadsPerAtom;
1309 // TODO: pick smaller block size in runtime if needed
1310 // (e.g. on 660 Ti where 50% occupancy is ~25% faster than 100% occupancy with RNAse (~17.8k atoms))
1311 // If doing so, change atomsPerBlock in the kernels as well.
1312 // TODO: test varying block sizes on modern arch-s as well
1313 // TODO: also consider using cudaFuncSetCacheConfig() for preferring shared memory on older architectures
1314 //(for spline data mostly)
1315 GMX_ASSERT(!(c_pmeAtomDataBlockSize % atomsPerBlock),
1316 "inconsistent atom data padding vs. spreading block size");
1318 // Ensure that coordinates are ready on the device before launching spread;
1319 // only needed with CUDA on PP+PME ranks, not on separate PME ranks, in unit tests
1320 // nor in OpenCL as these cases use a single stream (hence xReadyOnDevice == nullptr).
1321 GMX_ASSERT(!GMX_GPU_CUDA || xReadyOnDevice != nullptr || pmeGpu->common->isRankPmeOnly
1322 || pme_gpu_settings(pmeGpu).copyAllOutputs,
1323 "Need a valid coordinate synchronizer on PP+PME ranks with CUDA.");
1327 xReadyOnDevice->enqueueWaitEvent(pmeGpu->archSpecific->pmeStream_);
1330 const int blockCount = pmeGpu->nAtomsAlloc / atomsPerBlock;
1331 auto dimGrid = pmeGpuCreateGrid(pmeGpu, blockCount);
1333 if (pmeGpu->common->ngrids == 1)
1335 kernelParamsPtr->current.scale = 1.0;
1339 kernelParamsPtr->current.scale = 1.0 - lambda;
1342 KernelLaunchConfig config;
1343 config.blockSize[0] = order;
1344 config.blockSize[1] = (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? 1 : order);
1345 config.blockSize[2] = atomsPerBlock;
1346 config.gridSize[0] = dimGrid.first;
1347 config.gridSize[1] = dimGrid.second;
1350 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1355 timingId = PmeStage::SplineAndSpread;
1356 kernelPtr = selectSplineAndSpreadKernelPtr(pmeGpu,
1357 pmeGpu->settings.threadsPerAtom,
1358 writeGlobal || (!recalculateSplines),
1359 pmeGpu->common->ngrids);
1363 timingId = PmeStage::Spline;
1364 kernelPtr = selectSplineKernelPtr(pmeGpu,
1365 pmeGpu->settings.threadsPerAtom,
1366 writeGlobal || (!recalculateSplines),
1367 pmeGpu->common->ngrids);
1372 timingId = PmeStage::Spread;
1373 kernelPtr = selectSpreadKernelPtr(pmeGpu,
1374 pmeGpu->settings.threadsPerAtom,
1375 writeGlobal || (!recalculateSplines),
1376 pmeGpu->common->ngrids);
1380 pme_gpu_start_timing(pmeGpu, timingId);
1381 auto* timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1382 #if c_canEmbedBuffers
1383 const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1385 const auto kernelArgs =
1386 prepareGpuKernelArguments(kernelPtr,
1389 &kernelParamsPtr->atoms.d_theta,
1390 &kernelParamsPtr->atoms.d_dtheta,
1391 &kernelParamsPtr->atoms.d_gridlineIndices,
1392 &kernelParamsPtr->grid.d_realGrid[FEP_STATE_A],
1393 &kernelParamsPtr->grid.d_realGrid[FEP_STATE_B],
1394 &kernelParamsPtr->grid.d_fractShiftsTable,
1395 &kernelParamsPtr->grid.d_gridlineIndicesTable,
1396 &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_A],
1397 &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_B],
1398 &kernelParamsPtr->atoms.d_coordinates);
1402 kernelPtr, config, pmeGpu->archSpecific->pmeStream_, timingEvent, "PME spline/spread", kernelArgs);
1403 pme_gpu_stop_timing(pmeGpu, timingId);
1405 const auto& settings = pmeGpu->settings;
1406 const bool copyBackGrid = spreadCharges && (!settings.performGPUFFT || settings.copyAllOutputs);
1409 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
1411 float* h_grid = h_grids[gridIndex];
1412 pme_gpu_copy_output_spread_grid(pmeGpu, h_grid, gridIndex);
1415 const bool copyBackAtomData =
1416 computeSplines && (!settings.performGPUGather || settings.copyAllOutputs);
1417 if (copyBackAtomData)
1419 pme_gpu_copy_output_spread_atom_data(pmeGpu);
1423 void pme_gpu_solve(const PmeGpu* pmeGpu,
1424 const int gridIndex,
1426 GridOrdering gridOrdering,
1427 bool computeEnergyAndVirial)
1430 pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
1431 "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
1432 GMX_ASSERT(gridIndex < pmeGpu->common->ngrids,
1433 "Invalid combination of gridIndex and number of grids");
1435 const auto& settings = pmeGpu->settings;
1436 const bool copyInputAndOutputGrid = !settings.performGPUFFT || settings.copyAllOutputs;
1438 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1440 float* h_gridFloat = reinterpret_cast<float*>(h_grid);
1441 if (copyInputAndOutputGrid)
1443 copyToDeviceBuffer(&kernelParamsPtr->grid.d_fourierGrid[gridIndex],
1446 pmeGpu->archSpecific->complexGridSize[gridIndex],
1447 pmeGpu->archSpecific->pmeStream_,
1448 pmeGpu->settings.transferKind,
1452 int majorDim = -1, middleDim = -1, minorDim = -1;
1453 switch (gridOrdering)
1455 case GridOrdering::YZX:
1461 case GridOrdering::XYZ:
1467 default: GMX_ASSERT(false, "Implement grid ordering here and below for the kernel launch");
1470 const int maxBlockSize = pmeGpu->programHandle_->impl_->solveMaxWorkGroupSize;
1472 const int gridLineSize = pmeGpu->kernelParams->grid.complexGridSize[minorDim];
1473 const int gridLinesPerBlock = std::max(maxBlockSize / gridLineSize, 1);
1474 const int blocksPerGridLine = (gridLineSize + maxBlockSize - 1) / maxBlockSize;
1476 if (blocksPerGridLine == 1)
1478 cellsPerBlock = gridLineSize * gridLinesPerBlock;
1482 cellsPerBlock = (gridLineSize + blocksPerGridLine - 1) / blocksPerGridLine;
1484 const int warpSize = pmeGpu->programHandle_->warpSize();
1485 const int blockSize = (cellsPerBlock + warpSize - 1) / warpSize * warpSize;
1487 static_assert(!GMX_GPU_CUDA || c_solveMaxWarpsPerBlock / 2 >= 4,
1488 "The CUDA solve energy kernels needs at least 4 warps. "
1489 "Here we launch at least half of the max warps.");
1491 KernelLaunchConfig config;
1492 config.blockSize[0] = blockSize;
1493 config.gridSize[0] = blocksPerGridLine;
1494 // rounding up to full warps so that shuffle operations produce defined results
1495 config.gridSize[1] = (pmeGpu->kernelParams->grid.complexGridSize[middleDim] + gridLinesPerBlock - 1)
1496 / gridLinesPerBlock;
1497 config.gridSize[2] = pmeGpu->kernelParams->grid.complexGridSize[majorDim];
1499 PmeStage timingId = PmeStage::Solve;
1500 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1501 if (gridOrdering == GridOrdering::YZX)
1505 kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveYZXEnergyKernelA
1506 : pmeGpu->programHandle_->impl_->solveYZXKernelA;
1510 kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveYZXEnergyKernelB
1511 : pmeGpu->programHandle_->impl_->solveYZXKernelB;
1514 else if (gridOrdering == GridOrdering::XYZ)
1518 kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveXYZEnergyKernelA
1519 : pmeGpu->programHandle_->impl_->solveXYZKernelA;
1523 kernelPtr = computeEnergyAndVirial ? pmeGpu->programHandle_->impl_->solveXYZEnergyKernelB
1524 : pmeGpu->programHandle_->impl_->solveXYZKernelB;
1528 pme_gpu_start_timing(pmeGpu, timingId);
1529 auto* timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1530 #if c_canEmbedBuffers
1531 const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1533 const auto kernelArgs =
1534 prepareGpuKernelArguments(kernelPtr,
1537 &kernelParamsPtr->grid.d_splineModuli[gridIndex],
1538 &kernelParamsPtr->constants.d_virialAndEnergy[gridIndex],
1539 &kernelParamsPtr->grid.d_fourierGrid[gridIndex]);
1541 launchGpuKernel(kernelPtr, config, pmeGpu->archSpecific->pmeStream_, timingEvent, "PME solve", kernelArgs);
1542 pme_gpu_stop_timing(pmeGpu, timingId);
1544 if (computeEnergyAndVirial)
1546 copyFromDeviceBuffer(pmeGpu->staging.h_virialAndEnergy[gridIndex],
1547 &kernelParamsPtr->constants.d_virialAndEnergy[gridIndex],
1549 c_virialAndEnergyCount,
1550 pmeGpu->archSpecific->pmeStream_,
1551 pmeGpu->settings.transferKind,
1555 if (copyInputAndOutputGrid)
1557 copyFromDeviceBuffer(h_gridFloat,
1558 &kernelParamsPtr->grid.d_fourierGrid[gridIndex],
1560 pmeGpu->archSpecific->complexGridSize[gridIndex],
1561 pmeGpu->archSpecific->pmeStream_,
1562 pmeGpu->settings.transferKind,
1568 * Returns a pointer to appropriate gather kernel based on the inputvalues
1570 * \param[in] pmeGpu The PME GPU structure.
1571 * \param[in] threadsPerAtom Controls whether we should use order or order*order threads per atom
1572 * \param[in] readSplinesFromGlobal bool controlling if we should write spline data to global memory
1573 * \param[in] numGrids Number of grids to use. numGrids == 2 if Coulomb is perturbed.
1575 * \return Pointer to CUDA kernel
1577 inline auto selectGatherKernelPtr(const PmeGpu* pmeGpu,
1578 ThreadsPerAtom threadsPerAtom,
1579 bool readSplinesFromGlobal,
1583 PmeGpuProgramImpl::PmeKernelHandle kernelPtr = nullptr;
1585 if (readSplinesFromGlobal)
1587 if (threadsPerAtom == ThreadsPerAtom::Order)
1591 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesThPerAtom4Dual;
1595 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesThPerAtom4Single;
1602 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesDual;
1606 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelReadSplinesSingle;
1612 if (threadsPerAtom == ThreadsPerAtom::Order)
1616 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelThPerAtom4Dual;
1620 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelThPerAtom4Single;
1627 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelDual;
1631 kernelPtr = pmeGpu->programHandle_->impl_->gatherKernelSingle;
1638 void pme_gpu_gather(PmeGpu* pmeGpu, real** h_grids, const float lambda)
1641 pmeGpu->common->ngrids == 1 || pmeGpu->common->ngrids == 2,
1642 "Only one (normal Coulomb PME) or two (FEP coulomb PME) PME grids can be used on GPU");
1644 const auto& settings = pmeGpu->settings;
1646 if (!settings.performGPUFFT || settings.copyAllOutputs)
1648 for (int gridIndex = 0; gridIndex < pmeGpu->common->ngrids; gridIndex++)
1650 float* h_grid = const_cast<float*>(h_grids[gridIndex]);
1651 pme_gpu_copy_input_gather_grid(pmeGpu, h_grid, gridIndex);
1655 if (settings.copyAllOutputs)
1657 pme_gpu_copy_input_gather_atom_data(pmeGpu);
1660 /* Set if we have unit tests */
1661 const bool readGlobal = pmeGpu->settings.copyAllOutputs;
1662 const size_t blockSize = pmeGpu->programHandle_->impl_->gatherWorkGroupSize;
1663 const int order = pmeGpu->common->pme_order;
1664 GMX_ASSERT(order == c_pmeGpuOrder, "Only PME order 4 is implemented");
1665 const int threadsPerAtom =
1666 (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
1667 const bool recalculateSplines = pmeGpu->settings.recalculateSplines;
1669 GMX_ASSERT(!GMX_GPU_OPENCL || pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::OrderSquared,
1670 "Only 16 threads per atom supported in OpenCL");
1671 GMX_ASSERT(!GMX_GPU_OPENCL || !recalculateSplines,
1672 "Recalculating splines not supported in OpenCL");
1674 const int atomsPerBlock = blockSize / threadsPerAtom;
1676 GMX_ASSERT(!(c_pmeAtomDataBlockSize % atomsPerBlock),
1677 "inconsistent atom data padding vs. gathering block size");
1679 const int blockCount = pmeGpu->nAtomsAlloc / atomsPerBlock;
1680 auto dimGrid = pmeGpuCreateGrid(pmeGpu, blockCount);
1682 KernelLaunchConfig config;
1683 config.blockSize[0] = order;
1684 config.blockSize[1] = (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? 1 : order);
1685 config.blockSize[2] = atomsPerBlock;
1686 config.gridSize[0] = dimGrid.first;
1687 config.gridSize[1] = dimGrid.second;
1689 // TODO test different cache configs
1691 PmeStage timingId = PmeStage::Gather;
1692 PmeGpuProgramImpl::PmeKernelHandle kernelPtr =
1693 selectGatherKernelPtr(pmeGpu,
1694 pmeGpu->settings.threadsPerAtom,
1695 readGlobal || (!recalculateSplines),
1696 pmeGpu->common->ngrids);
1697 // TODO design kernel selection getters and make PmeGpu a friend of PmeGpuProgramImpl
1699 pme_gpu_start_timing(pmeGpu, timingId);
1700 auto* timingEvent = pme_gpu_fetch_timing_event(pmeGpu, timingId);
1701 auto* kernelParamsPtr = pmeGpu->kernelParams.get();
1702 if (pmeGpu->common->ngrids == 1)
1704 kernelParamsPtr->current.scale = 1.0;
1708 kernelParamsPtr->current.scale = 1.0 - lambda;
1711 #if c_canEmbedBuffers
1712 const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, kernelParamsPtr);
1714 const auto kernelArgs =
1715 prepareGpuKernelArguments(kernelPtr,
1718 &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_A],
1719 &kernelParamsPtr->atoms.d_coefficients[FEP_STATE_B],
1720 &kernelParamsPtr->grid.d_realGrid[FEP_STATE_A],
1721 &kernelParamsPtr->grid.d_realGrid[FEP_STATE_B],
1722 &kernelParamsPtr->atoms.d_theta,
1723 &kernelParamsPtr->atoms.d_dtheta,
1724 &kernelParamsPtr->atoms.d_gridlineIndices,
1725 &kernelParamsPtr->atoms.d_forces);
1727 launchGpuKernel(kernelPtr, config, pmeGpu->archSpecific->pmeStream_, timingEvent, "PME gather", kernelArgs);
1728 pme_gpu_stop_timing(pmeGpu, timingId);
1730 if (pmeGpu->settings.useGpuForceReduction)
1732 pmeGpu->archSpecific->pmeForcesReady.markEvent(pmeGpu->archSpecific->pmeStream_);
1736 pme_gpu_copy_output_forces(pmeGpu);
1740 DeviceBuffer<gmx::RVec> pme_gpu_get_kernelparam_forces(const PmeGpu* pmeGpu)
1742 if (pmeGpu && pmeGpu->kernelParams)
1744 return pmeGpu->kernelParams->atoms.d_forces;
1748 return DeviceBuffer<gmx::RVec>{};
1752 void pme_gpu_set_kernelparam_coordinates(const PmeGpu* pmeGpu, DeviceBuffer<gmx::RVec> d_x)
1754 GMX_ASSERT(pmeGpu && pmeGpu->kernelParams,
1755 "PME GPU device buffer can not be set in non-GPU builds or before the GPU PME was "
1758 GMX_ASSERT(checkDeviceBuffer(d_x, pmeGpu->kernelParams->atoms.nAtoms),
1759 "The device-side buffer can not be set.");
1761 pmeGpu->kernelParams->atoms.d_coordinates = d_x;
1764 GpuEventSynchronizer* pme_gpu_get_forces_ready_synchronizer(const PmeGpu* pmeGpu)
1766 if (pmeGpu && pmeGpu->kernelParams)
1768 return &pmeGpu->archSpecific->pmeForcesReady;