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.
38 * Implements common routines for PME tests.
40 * \author Aleksei Iupinov <a.yupinov@gmail.com>
41 * \ingroup module_ewald
45 #include "pmetestcommon.h"
51 #include "gromacs/domdec/domdec.h"
52 #include "gromacs/ewald/pme_gather.h"
53 #include "gromacs/ewald/pme_gpu_calculate_splines.h"
54 #include "gromacs/ewald/pme_gpu_constants.h"
55 #include "gromacs/ewald/pme_gpu_internal.h"
56 #include "gromacs/ewald/pme_gpu_staging.h"
57 #include "gromacs/ewald/pme_grid.h"
58 #include "gromacs/ewald/pme_internal.h"
59 #include "gromacs/ewald/pme_redistribute.h"
60 #include "gromacs/ewald/pme_solve.h"
61 #include "gromacs/ewald/pme_spread.h"
62 #include "gromacs/fft/parallel_3dfft.h"
63 #include "gromacs/gpu_utils/gpu_utils.h"
64 #include "gromacs/hardware/device_management.h"
65 #include "gromacs/math/invertmatrix.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/pbcutil/pbc.h"
68 #include "gromacs/topology/topology.h"
69 #include "gromacs/utility/exceptions.h"
70 #include "gromacs/utility/gmxassert.h"
71 #include "gromacs/utility/logger.h"
72 #include "gromacs/utility/stringutil.h"
74 #include "testutils/test_hardware_environment.h"
75 #include "testutils/testasserts.h"
84 bool pmeSupportsInputForMode(const gmx_hw_info_t& hwinfo, const t_inputrec* inputRec, CodePath mode)
89 case CodePath::CPU: implemented = true; break;
92 implemented = (pme_gpu_supports_build(nullptr) && pme_gpu_supports_hardware(hwinfo, nullptr)
93 && pme_gpu_supports_input(*inputRec, nullptr));
96 default: GMX_THROW(InternalError("Test not implemented for this mode"));
101 uint64_t getSplineModuliDoublePrecisionUlps(int splineOrder)
103 /* Arbitrary ulp tolerance for sine/cosine implementation. It's
104 * hard to know what to pick without testing lots of
105 * implementations. */
106 const uint64_t sineUlps = 10;
107 return 4 * (splineOrder - 2) + 2 * sineUlps * splineOrder;
110 //! PME initialization
111 PmeSafePointer pmeInitWrapper(const t_inputrec* inputRec,
113 const DeviceContext* deviceContext,
114 const DeviceStream* deviceStream,
115 const PmeGpuProgram* pmeGpuProgram,
116 const Matrix3x3& box,
117 const real ewaldCoeff_q,
118 const real ewaldCoeff_lj)
120 const MDLogger dummyLogger;
121 const auto runMode = (mode == CodePath::CPU) ? PmeRunMode::CPU : PmeRunMode::Mixed;
122 t_commrec dummyCommrec = { 0 };
123 NumPmeDomains numPmeDomains = { 1, 1 };
124 gmx_pme_t* pmeDataRaw = gmx_pme_init(&dummyCommrec,
139 PmeSafePointer pme(pmeDataRaw); // taking ownership
141 // TODO get rid of this with proper matrix type
143 for (int i = 0; i < DIM; i++)
145 for (int j = 0; j < DIM; j++)
147 boxTemp[i][j] = box[i * DIM + j];
150 const char* boxError = check_box(PbcType::Unset, boxTemp);
151 GMX_RELEASE_ASSERT(boxError == nullptr, boxError);
155 case CodePath::CPU: invertBoxMatrix(boxTemp, pme->recipbox); break;
158 pme_gpu_set_testing(pme->gpu, true);
159 pme_gpu_update_input_box(pme->gpu, boxTemp);
162 default: GMX_THROW(InternalError("Test not implemented for this mode"));
168 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec)
170 const Matrix3x3 defaultBox = { { 1.0F, 0.0F, 0.0F, 0.0F, 1.0F, 0.0F, 0.0F, 0.0F, 1.0F } };
171 return pmeInitWrapper(inputRec, CodePath::CPU, nullptr, nullptr, nullptr, defaultBox, 0.0F, 0.0F);
174 //! Make a GPU state-propagator manager
175 std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme_t& pme,
176 const DeviceContext* deviceContext,
177 const DeviceStream* deviceStream)
179 // TODO: Pin the host buffer and use async memory copies
180 // TODO: Special constructor for PME-only rank / PME-tests is used here. There should be a mechanism to
181 // restrict one from using other constructor here.
182 return std::make_unique<StatePropagatorDataGpu>(
183 deviceStream, *deviceContext, GpuApiCallBehavior::Sync, pme_gpu_get_block_size(&pme), nullptr);
186 //! PME initialization with atom data
187 void pmeInitAtoms(gmx_pme_t* pme,
188 StatePropagatorDataGpu* stateGpu,
190 const CoordinatesVector& coordinates,
191 const ChargesVector& charges)
193 const index atomCount = coordinates.size();
194 GMX_RELEASE_ASSERT(atomCount == charges.ssize(), "Mismatch in atom data");
195 PmeAtomComm* atc = nullptr;
200 atc = &(pme->atc[0]);
201 atc->x = coordinates;
202 atc->coefficient = charges;
203 gmx_pme_reinit_atoms(pme, atomCount, charges, {});
204 /* With decomposition there would be more boilerplate atc code here, e.g. do_redist_pos_coeffs */
208 // TODO: Avoid use of atc in the GPU code path
209 atc = &(pme->atc[0]);
210 // We need to set atc->n for passing the size in the tests
211 atc->setNumAtoms(atomCount);
212 gmx_pme_reinit_atoms(pme, atomCount, charges, {});
214 stateGpu->reinit(atomCount, atomCount);
215 stateGpu->copyCoordinatesToGpu(arrayRefFromArray(coordinates.data(), coordinates.size()),
216 gmx::AtomLocality::Local);
217 pme_gpu_set_kernelparam_coordinates(pme->gpu, stateGpu->getCoordinates());
221 default: GMX_THROW(InternalError("Test not implemented for this mode"));
225 //! Getting local PME real grid pointer for test I/O
226 static real* pmeGetRealGridInternal(const gmx_pme_t* pme)
228 const size_t gridIndex = 0;
229 return pme->fftgrid[gridIndex];
232 //! Getting local PME real grid dimensions
233 static void pmeGetRealGridSizesInternal(const gmx_pme_t* pme,
235 IVec& gridSize, //NOLINT(google-runtime-references)
236 IVec& paddedGridSize) //NOLINT(google-runtime-references)
238 const size_t gridIndex = 0;
239 IVec gridOffsetUnused;
243 gmx_parallel_3dfft_real_limits(
244 pme->pfft_setup[gridIndex], gridSize, gridOffsetUnused, paddedGridSize);
248 pme_gpu_get_real_grid_sizes(pme->gpu, &gridSize, &paddedGridSize);
251 default: GMX_THROW(InternalError("Test not implemented for this mode"));
255 //! Getting local PME complex grid pointer for test I/O
256 static t_complex* pmeGetComplexGridInternal(const gmx_pme_t* pme)
258 const size_t gridIndex = 0;
259 return pme->cfftgrid[gridIndex];
262 //! Getting local PME complex grid dimensions
263 static void pmeGetComplexGridSizesInternal(const gmx_pme_t* pme,
264 IVec& gridSize, //NOLINT(google-runtime-references)
265 IVec& paddedGridSize) //NOLINT(google-runtime-references)
267 const size_t gridIndex = 0;
268 IVec gridOffsetUnused, complexOrderUnused;
269 gmx_parallel_3dfft_complex_limits(
270 pme->pfft_setup[gridIndex], complexOrderUnused, gridSize, gridOffsetUnused, paddedGridSize); // TODO: what about YZX ordering?
273 //! Getting the PME grid memory buffer and its sizes - template definition
274 template<typename ValueType>
275 static void pmeGetGridAndSizesInternal(const gmx_pme_t* /*unused*/,
277 ValueType*& /*unused*/, //NOLINT(google-runtime-references)
278 IVec& /*unused*/, //NOLINT(google-runtime-references)
279 IVec& /*unused*/) //NOLINT(google-runtime-references)
281 GMX_THROW(InternalError("Deleted function call"));
282 // explicitly deleting general template does not compile in clang, see https://llvm.org/bugs/show_bug.cgi?id=17537
285 //! Getting the PME real grid memory buffer and its sizes
287 void pmeGetGridAndSizesInternal<real>(const gmx_pme_t* pme, CodePath mode, real*& grid, IVec& gridSize, IVec& paddedGridSize)
289 grid = pmeGetRealGridInternal(pme);
290 pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSize);
293 //! Getting the PME complex grid memory buffer and its sizes
295 void pmeGetGridAndSizesInternal<t_complex>(const gmx_pme_t* pme,
299 IVec& paddedGridSize)
301 grid = pmeGetComplexGridInternal(pme);
302 pmeGetComplexGridSizesInternal(pme, gridSize, paddedGridSize);
305 //! PME spline calculation and charge spreading
306 void pmePerformSplineAndSpread(gmx_pme_t* pme,
307 CodePath mode, // TODO const qualifiers elsewhere
311 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
312 PmeAtomComm* atc = &(pme->atc[0]);
313 const size_t gridIndex = 0;
314 const bool computeSplinesForZeroCharges = true;
315 real** fftgrid = spreadCharges ? pme->fftgrid : nullptr;
316 real* pmegrid = pme->pmegrid[gridIndex].grid.grid;
323 &pme->pmegrid[gridIndex],
326 fftgrid != nullptr ? fftgrid[gridIndex] : nullptr,
327 computeSplinesForZeroCharges,
329 if (spreadCharges && !pme->bUseThreads)
331 wrap_periodic_pmegrid(pme, pmegrid);
332 copy_pmegrid_to_fftgrid(
333 pme, pmegrid, fftgrid != nullptr ? fftgrid[gridIndex] : nullptr, gridIndex);
337 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
338 * double precision. GPUs are not used with double precision anyhow. */
342 const real lambdaQ = 1.0;
343 // no synchronization needed as x is transferred in the PME stream
344 GpuEventSynchronizer* xReadyOnDevice = nullptr;
345 pme_gpu_spread(pme->gpu, xReadyOnDevice, fftgrid, computeSplines, spreadCharges, lambdaQ);
350 default: GMX_THROW(InternalError("Test not implemented for this mode"));
354 //! Getting the internal spline data buffer pointer
355 static real* pmeGetSplineDataInternal(const gmx_pme_t* pme, PmeSplineDataType type, int dimIndex)
357 GMX_ASSERT((0 <= dimIndex) && (dimIndex < DIM), "Invalid dimension index");
358 const PmeAtomComm* atc = &(pme->atc[0]);
359 const size_t threadIndex = 0;
360 real* splineBuffer = nullptr;
363 case PmeSplineDataType::Values:
364 splineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
367 case PmeSplineDataType::Derivatives:
368 splineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
371 default: GMX_THROW(InternalError("Unknown spline data type"));
377 void pmePerformSolve(const gmx_pme_t* pme,
379 PmeSolveAlgorithm method,
381 GridOrdering gridOrdering,
382 bool computeEnergyAndVirial)
384 t_complex* h_grid = pmeGetComplexGridInternal(pme);
385 const bool useLorentzBerthelot = false;
386 const size_t threadIndex = 0;
387 const size_t gridIndex = 0;
391 if (gridOrdering != GridOrdering::YZX)
393 GMX_THROW(InternalError("Test not implemented for this mode"));
397 case PmeSolveAlgorithm::Coulomb:
398 solve_pme_yzx(pme, h_grid, cellVolume, computeEnergyAndVirial, pme->nthread, threadIndex);
401 case PmeSolveAlgorithm::LennardJones:
402 solve_pme_lj_yzx(pme,
406 computeEnergyAndVirial,
411 default: GMX_THROW(InternalError("Test not implemented for this mode"));
418 case PmeSolveAlgorithm::Coulomb:
419 pme_gpu_solve(pme->gpu, gridIndex, h_grid, gridOrdering, computeEnergyAndVirial);
422 default: GMX_THROW(InternalError("Test not implemented for this mode"));
426 default: GMX_THROW(InternalError("Test not implemented for this mode"));
430 //! PME force gathering
431 void pmePerformGather(gmx_pme_t* pme, CodePath mode, ForcesVector& forces)
433 PmeAtomComm* atc = &(pme->atc[0]);
434 const index atomCount = atc->numAtoms();
435 GMX_RELEASE_ASSERT(forces.ssize() == atomCount, "Invalid force buffer size");
436 const real scale = 1.0;
437 const size_t threadIndex = 0;
438 const size_t gridIndex = 0;
439 real* pmegrid = pme->pmegrid[gridIndex].grid.grid;
440 real** fftgrid = pme->fftgrid;
446 if (atc->nthread == 1)
448 // something which is normally done in serial spline computation (make_thread_local_ind())
449 atc->spline[threadIndex].n = atomCount;
451 copy_fftgrid_to_pmegrid(pme, fftgrid[gridIndex], pmegrid, gridIndex, pme->nthread, threadIndex);
452 unwrap_periodic_pmegrid(pme, pmegrid);
453 gather_f_bsplines(pme, pmegrid, true, atc, &atc->spline[threadIndex], scale);
456 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
457 * double precision. GPUs are not used with double precision anyhow. */
461 // Variable initialization needs a non-switch scope
462 const bool computeEnergyAndVirial = false;
463 const real lambdaQ = 1.0;
464 PmeOutput output = pme_gpu_getOutput(*pme, computeEnergyAndVirial, lambdaQ);
465 GMX_ASSERT(forces.size() == output.forces_.size(),
466 "Size of force buffers did not match");
467 pme_gpu_gather(pme->gpu, fftgrid, lambdaQ);
468 std::copy(std::begin(output.forces_), std::end(output.forces_), std::begin(forces));
473 default: GMX_THROW(InternalError("Test not implemented for this mode"));
477 //! PME test finalization before fetching the outputs
478 void pmeFinalizeTest(const gmx_pme_t* pme, CodePath mode)
482 case CodePath::CPU: break;
484 case CodePath::GPU: pme_gpu_synchronize(pme->gpu); break;
486 default: GMX_THROW(InternalError("Test not implemented for this mode"));
490 //! A binary enum for spline data layout transformation
491 enum class PmeLayoutTransform
497 /*! \brief Gets a unique index to an element in a spline parameter buffer.
499 * These theta/dtheta buffers are laid out for GPU spread/gather
500 * kernels. The index is wrt the execution block, in range(0,
501 * atomsPerBlock * order * DIM).
503 * This is a wrapper, only used in unit tests.
504 * \param[in] order PME order
505 * \param[in] splineIndex Spline contribution index (from 0 to \p order - 1)
506 * \param[in] dimIndex Dimension index (from 0 to 2)
507 * \param[in] atomIndex Atom index wrt the block.
508 * \param[in] atomsPerWarp Number of atoms processed by a warp.
510 * \returns Index into theta or dtheta array using GPU layout.
512 static int getSplineParamFullIndex(int order, int splineIndex, int dimIndex, int atomIndex, int atomsPerWarp)
514 if (order != c_pmeGpuOrder)
518 constexpr int fixedOrder = c_pmeGpuOrder;
519 GMX_UNUSED_VALUE(fixedOrder);
521 const int atomWarpIndex = atomIndex % atomsPerWarp;
522 const int warpIndex = atomIndex / atomsPerWarp;
523 int indexBase, result;
524 switch (atomsPerWarp)
527 indexBase = getSplineParamIndexBase<fixedOrder, 1>(warpIndex, atomWarpIndex);
528 result = getSplineParamIndex<fixedOrder, 1>(indexBase, dimIndex, splineIndex);
532 indexBase = getSplineParamIndexBase<fixedOrder, 2>(warpIndex, atomWarpIndex);
533 result = getSplineParamIndex<fixedOrder, 2>(indexBase, dimIndex, splineIndex);
537 indexBase = getSplineParamIndexBase<fixedOrder, 4>(warpIndex, atomWarpIndex);
538 result = getSplineParamIndex<fixedOrder, 4>(indexBase, dimIndex, splineIndex);
542 indexBase = getSplineParamIndexBase<fixedOrder, 8>(warpIndex, atomWarpIndex);
543 result = getSplineParamIndex<fixedOrder, 8>(indexBase, dimIndex, splineIndex);
547 GMX_THROW(NotImplementedError(
548 formatString("Test function call not unrolled for atomsPerWarp = %d in "
549 "getSplineParamFullIndex",
555 /*!\brief Return the number of atoms per warp */
556 static int pme_gpu_get_atoms_per_warp(const PmeGpu* pmeGpu)
558 const int order = pmeGpu->common->pme_order;
559 const int threadsPerAtom =
560 (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
561 return pmeGpu->programHandle_->warpSize() / threadsPerAtom;
564 /*! \brief Rearranges the atom spline data between the GPU and host layouts.
565 * Only used for test purposes so far, likely to be horribly slow.
567 * \param[in] pmeGpu The PME GPU structure.
568 * \param[out] atc The PME CPU atom data structure (with a single-threaded layout).
569 * \param[in] type The spline data type (values or derivatives).
570 * \param[in] dimIndex Dimension index.
571 * \param[in] transform Layout transform type
573 static void pme_gpu_transform_spline_atom_data(const PmeGpu* pmeGpu,
574 const PmeAtomComm* atc,
575 PmeSplineDataType type,
577 PmeLayoutTransform transform)
579 // The GPU atom spline data is laid out in a different way currently than the CPU one.
580 // This function converts the data from GPU to CPU layout (in the host memory).
581 // It is only intended for testing purposes so far.
582 // Ideally we should use similar layouts on CPU and GPU if we care about mixed modes and their
583 // performance (e.g. spreading on GPU, gathering on CPU).
584 GMX_RELEASE_ASSERT(atc->nthread == 1, "Only the serial PME data layout is supported");
585 const uintmax_t threadIndex = 0;
586 const auto atomCount = atc->numAtoms();
587 const auto atomsPerWarp = pme_gpu_get_atoms_per_warp(pmeGpu);
588 const auto pmeOrder = pmeGpu->common->pme_order;
589 GMX_ASSERT(pmeOrder == c_pmeGpuOrder, "Only PME order 4 is implemented");
591 real* cpuSplineBuffer;
592 float* h_splineBuffer;
595 case PmeSplineDataType::Values:
596 cpuSplineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
597 h_splineBuffer = pmeGpu->staging.h_theta;
600 case PmeSplineDataType::Derivatives:
601 cpuSplineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
602 h_splineBuffer = pmeGpu->staging.h_dtheta;
605 default: GMX_THROW(InternalError("Unknown spline data type"));
608 for (auto atomIndex = 0; atomIndex < atomCount; atomIndex++)
610 for (auto orderIndex = 0; orderIndex < pmeOrder; orderIndex++)
612 const auto gpuValueIndex =
613 getSplineParamFullIndex(pmeOrder, orderIndex, dimIndex, atomIndex, atomsPerWarp);
614 const auto cpuValueIndex = atomIndex * pmeOrder + orderIndex;
615 GMX_ASSERT(cpuValueIndex < atomCount * pmeOrder,
616 "Atom spline data index out of bounds (while transforming GPU data layout "
620 case PmeLayoutTransform::GpuToHost:
621 cpuSplineBuffer[cpuValueIndex] = h_splineBuffer[gpuValueIndex];
624 case PmeLayoutTransform::HostToGpu:
625 h_splineBuffer[gpuValueIndex] = cpuSplineBuffer[cpuValueIndex];
628 default: GMX_THROW(InternalError("Unknown layout transform"));
634 //! Setting atom spline values/derivatives to be used in spread/gather
635 void pmeSetSplineData(const gmx_pme_t* pme,
637 const SplineParamsDimVector& splineValues,
638 PmeSplineDataType type,
641 const PmeAtomComm* atc = &(pme->atc[0]);
642 const index atomCount = atc->numAtoms();
643 const index pmeOrder = pme->pme_order;
644 const index dimSize = pmeOrder * atomCount;
645 GMX_RELEASE_ASSERT(dimSize == splineValues.ssize(), "Mismatch in spline data");
646 real* splineBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
651 std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
655 std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
656 pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::HostToGpu);
659 default: GMX_THROW(InternalError("Test not implemented for this mode"));
663 //! Setting gridline indices to be used in spread/gather
664 void pmeSetGridLineIndices(gmx_pme_t* pme, CodePath mode, const GridLineIndicesVector& gridLineIndices)
666 PmeAtomComm* atc = &(pme->atc[0]);
667 const index atomCount = atc->numAtoms();
668 GMX_RELEASE_ASSERT(atomCount == gridLineIndices.ssize(), "Mismatch in gridline indices size");
670 IVec paddedGridSizeUnused, gridSize(0, 0, 0);
671 pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSizeUnused);
673 for (const auto& index : gridLineIndices)
675 for (int i = 0; i < DIM; i++)
677 GMX_RELEASE_ASSERT((0 <= index[i]) && (index[i] < gridSize[i]),
678 "Invalid gridline index");
685 memcpy(pme_gpu_staging(pme->gpu).h_gridlineIndices,
686 gridLineIndices.data(),
687 atomCount * sizeof(gridLineIndices[0]));
691 atc->idx.resize(gridLineIndices.size());
692 std::copy(gridLineIndices.begin(), gridLineIndices.end(), atc->idx.begin());
694 default: GMX_THROW(InternalError("Test not implemented for this mode"));
698 //! Getting plain index into the complex 3d grid
699 inline size_t pmeGetGridPlainIndexInternal(const IVec& index, const IVec& paddedGridSize, GridOrdering gridOrdering)
702 switch (gridOrdering)
704 case GridOrdering::YZX:
705 result = (index[YY] * paddedGridSize[ZZ] + index[ZZ]) * paddedGridSize[XX] + index[XX];
708 case GridOrdering::XYZ:
709 result = (index[XX] * paddedGridSize[YY] + index[YY]) * paddedGridSize[ZZ] + index[ZZ];
712 default: GMX_THROW(InternalError("Test not implemented for this mode"));
717 //! Setting real or complex grid
718 template<typename ValueType>
719 static void pmeSetGridInternal(const gmx_pme_t* pme,
721 GridOrdering gridOrdering,
722 const SparseGridValuesInput<ValueType>& gridValues)
724 IVec gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
726 pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
730 case CodePath::GPU: // intentional absence of break, the grid will be copied from the host buffer in testing mode
732 std::memset(grid, 0, paddedGridSize[XX] * paddedGridSize[YY] * paddedGridSize[ZZ] * sizeof(ValueType));
733 for (const auto& gridValue : gridValues)
735 for (int i = 0; i < DIM; i++)
737 GMX_RELEASE_ASSERT((0 <= gridValue.first[i]) && (gridValue.first[i] < gridSize[i]),
738 "Invalid grid value index");
740 const size_t gridValueIndex =
741 pmeGetGridPlainIndexInternal(gridValue.first, paddedGridSize, gridOrdering);
742 grid[gridValueIndex] = gridValue.second;
746 default: GMX_THROW(InternalError("Test not implemented for this mode"));
750 //! Setting real grid to be used in gather
751 void pmeSetRealGrid(const gmx_pme_t* pme, CodePath mode, const SparseRealGridValuesInput& gridValues)
753 pmeSetGridInternal<real>(pme, mode, GridOrdering::XYZ, gridValues);
756 //! Setting complex grid to be used in solve
757 void pmeSetComplexGrid(const gmx_pme_t* pme,
759 GridOrdering gridOrdering,
760 const SparseComplexGridValuesInput& gridValues)
762 pmeSetGridInternal<t_complex>(pme, mode, gridOrdering, gridValues);
765 //! Getting the single dimension's spline values or derivatives
766 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t* pme, CodePath mode, PmeSplineDataType type, int dimIndex)
768 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
769 const PmeAtomComm* atc = &(pme->atc[0]);
770 const size_t atomCount = atc->numAtoms();
771 const size_t pmeOrder = pme->pme_order;
772 const size_t dimSize = pmeOrder * atomCount;
774 real* sourceBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
775 SplineParamsDimVector result;
779 pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::GpuToHost);
780 result = arrayRefFromArray(sourceBuffer, dimSize);
783 case CodePath::CPU: result = arrayRefFromArray(sourceBuffer, dimSize); break;
785 default: GMX_THROW(InternalError("Test not implemented for this mode"));
790 //! Getting the gridline indices
791 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t* pme, CodePath mode)
793 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
794 const PmeAtomComm* atc = &(pme->atc[0]);
795 const size_t atomCount = atc->numAtoms();
797 GridLineIndicesVector gridLineIndices;
801 gridLineIndices = arrayRefFromArray(
802 reinterpret_cast<IVec*>(pme_gpu_staging(pme->gpu).h_gridlineIndices), atomCount);
805 case CodePath::CPU: gridLineIndices = atc->idx; break;
807 default: GMX_THROW(InternalError("Test not implemented for this mode"));
809 return gridLineIndices;
812 //! Getting real or complex grid - only non zero values
813 template<typename ValueType>
814 static SparseGridValuesOutput<ValueType> pmeGetGridInternal(const gmx_pme_t* pme,
816 GridOrdering gridOrdering)
818 IVec gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
820 pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
821 SparseGridValuesOutput<ValueType> gridValues;
824 case CodePath::GPU: // intentional absence of break
827 for (int ix = 0; ix < gridSize[XX]; ix++)
829 for (int iy = 0; iy < gridSize[YY]; iy++)
831 for (int iz = 0; iz < gridSize[ZZ]; iz++)
833 IVec temp(ix, iy, iz);
834 const size_t gridValueIndex =
835 pmeGetGridPlainIndexInternal(temp, paddedGridSize, gridOrdering);
836 const ValueType value = grid[gridValueIndex];
837 if (value != ValueType{})
839 auto key = formatString("Cell %d %d %d", ix, iy, iz);
840 gridValues[key] = value;
847 default: GMX_THROW(InternalError("Test not implemented for this mode"));
852 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
853 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t* pme, CodePath mode)
855 return pmeGetGridInternal<real>(pme, mode, GridOrdering::XYZ);
858 //! Getting the complex grid output of pmePerformSolve()
859 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t* pme, CodePath mode, GridOrdering gridOrdering)
861 return pmeGetGridInternal<t_complex>(pme, mode, gridOrdering);
864 //! Getting the reciprocal energy and virial
865 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t* pme, CodePath mode, PmeSolveAlgorithm method)
868 const real lambdaQ = 1.0;
874 case PmeSolveAlgorithm::Coulomb:
875 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output);
878 case PmeSolveAlgorithm::LennardJones:
879 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output);
882 default: GMX_THROW(InternalError("Test not implemented for this mode"));
888 case PmeSolveAlgorithm::Coulomb:
889 pme_gpu_getEnergyAndVirial(*pme, lambdaQ, &output);
892 default: GMX_THROW(InternalError("Test not implemented for this mode"));
896 default: GMX_THROW(InternalError("Test not implemented for this mode"));
901 const char* codePathToString(CodePath codePath)
905 case CodePath::CPU: return "CPU";
906 case CodePath::GPU: return "GPU";
907 default: GMX_THROW(NotImplementedError("This CodePath should support codePathToString"));
911 PmeTestHardwareContext::PmeTestHardwareContext() : codePath_(CodePath::CPU) {}
913 PmeTestHardwareContext::PmeTestHardwareContext(TestDevice* testDevice) :
914 codePath_(CodePath::GPU), testDevice_(testDevice)
916 setActiveDevice(testDevice_->deviceInfo());
917 pmeGpuProgram_ = buildPmeGpuProgram(testDevice_->deviceContext());
920 //! Returns a human-readable context description line
921 std::string PmeTestHardwareContext::description() const
925 case CodePath::CPU: return "CPU";
926 case CodePath::GPU: return "GPU (" + testDevice_->description() + ")";
927 default: return "Unknown code path.";
931 void PmeTestHardwareContext::activate() const
933 if (codePath_ == CodePath::GPU)
935 setActiveDevice(testDevice_->deviceInfo());
939 std::vector<std::unique_ptr<PmeTestHardwareContext>> createPmeTestHardwareContextList()
941 std::vector<std::unique_ptr<PmeTestHardwareContext>> pmeTestHardwareContextList;
943 pmeTestHardwareContextList.emplace_back(std::make_unique<PmeTestHardwareContext>());
945 const auto& testDeviceList = getTestHardwareEnvironment()->getTestDeviceList();
946 for (const auto& testDevice : testDeviceList)
948 pmeTestHardwareContextList.emplace_back(std::make_unique<PmeTestHardwareContext>(testDevice.get()));
950 return pmeTestHardwareContextList;