2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019,2020,2021, 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.
37 * Implements common routines for PME tests.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
40 * \ingroup module_ewald
44 #include "pmetestcommon.h"
50 #include "gromacs/domdec/domdec.h"
51 #include "gromacs/ewald/pme_gather.h"
52 #include "gromacs/ewald/pme_gpu_calculate_splines.h"
53 #include "gromacs/ewald/pme_gpu_constants.h"
54 #include "gromacs/ewald/pme_gpu_internal.h"
55 #include "gromacs/ewald/pme_gpu_staging.h"
56 #include "gromacs/ewald/pme_grid.h"
57 #include "gromacs/ewald/pme_internal.h"
58 #include "gromacs/ewald/pme_redistribute.h"
59 #include "gromacs/ewald/pme_solve.h"
60 #include "gromacs/ewald/pme_spread.h"
61 #include "gromacs/fft/parallel_3dfft.h"
62 #include "gromacs/gpu_utils/gpu_utils.h"
63 #include "gromacs/hardware/device_management.h"
64 #include "gromacs/math/invertmatrix.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/topology/topology.h"
68 #include "gromacs/utility/exceptions.h"
69 #include "gromacs/utility/gmxassert.h"
70 #include "gromacs/utility/logger.h"
71 #include "gromacs/utility/stringutil.h"
73 #include "testutils/test_hardware_environment.h"
74 #include "testutils/testasserts.h"
83 bool pmeSupportsInputForMode(const gmx_hw_info_t& hwinfo, const t_inputrec* inputRec, CodePath mode)
88 case CodePath::CPU: implemented = true; break;
91 implemented = (pme_gpu_supports_build(nullptr) && pme_gpu_supports_hardware(hwinfo, nullptr)
92 && pme_gpu_supports_input(*inputRec, nullptr));
95 default: GMX_THROW(InternalError("Test not implemented for this mode"));
100 uint64_t getSplineModuliDoublePrecisionUlps(int splineOrder)
102 /* Arbitrary ulp tolerance for sine/cosine implementation. It's
103 * hard to know what to pick without testing lots of
104 * implementations. */
105 const uint64_t sineUlps = 10;
106 return 4 * (splineOrder - 2) + 2 * sineUlps * splineOrder;
109 //! PME initialization
110 PmeSafePointer pmeInitWrapper(const t_inputrec* inputRec,
112 const DeviceContext* deviceContext,
113 const DeviceStream* deviceStream,
114 const PmeGpuProgram* pmeGpuProgram,
115 const Matrix3x3& box,
116 const real ewaldCoeff_q,
117 const real ewaldCoeff_lj)
119 const MDLogger dummyLogger;
120 const auto runMode = (mode == CodePath::CPU) ? PmeRunMode::CPU : PmeRunMode::Mixed;
121 t_commrec dummyCommrec = { 0 };
122 NumPmeDomains numPmeDomains = { 1, 1 };
123 gmx_pme_t* pmeDataRaw = gmx_pme_init(&dummyCommrec,
138 PmeSafePointer pme(pmeDataRaw); // taking ownership
140 // TODO get rid of this with proper matrix type
142 for (int i = 0; i < DIM; i++)
144 for (int j = 0; j < DIM; j++)
146 boxTemp[i][j] = box[i * DIM + j];
149 const char* boxError = check_box(PbcType::Unset, boxTemp);
150 GMX_RELEASE_ASSERT(boxError == nullptr, boxError);
154 case CodePath::CPU: invertBoxMatrix(boxTemp, pme->recipbox); break;
157 pme_gpu_set_testing(pme->gpu, true);
158 pme_gpu_update_input_box(pme->gpu, boxTemp);
161 default: GMX_THROW(InternalError("Test not implemented for this mode"));
167 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec)
169 const Matrix3x3 defaultBox = { { 1.0F, 0.0F, 0.0F, 0.0F, 1.0F, 0.0F, 0.0F, 0.0F, 1.0F } };
170 return pmeInitWrapper(inputRec, CodePath::CPU, nullptr, nullptr, nullptr, defaultBox, 0.0F, 0.0F);
173 //! Make a GPU state-propagator manager
174 std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme_t& pme,
175 const DeviceContext* deviceContext,
176 const DeviceStream* deviceStream)
178 // TODO: Pin the host buffer and use async memory copies
179 // TODO: Special constructor for PME-only rank / PME-tests is used here. There should be a mechanism to
180 // restrict one from using other constructor here.
181 return std::make_unique<StatePropagatorDataGpu>(
182 deviceStream, *deviceContext, GpuApiCallBehavior::Sync, pme_gpu_get_block_size(&pme), nullptr);
185 //! PME initialization with atom data
186 void pmeInitAtoms(gmx_pme_t* pme,
187 StatePropagatorDataGpu* stateGpu,
189 const CoordinatesVector& coordinates,
190 const ChargesVector& charges)
192 const index atomCount = coordinates.size();
193 GMX_RELEASE_ASSERT(atomCount == charges.ssize(), "Mismatch in atom data");
194 PmeAtomComm* atc = nullptr;
199 atc = &(pme->atc[0]);
200 atc->x = coordinates;
201 atc->coefficient = charges;
202 gmx_pme_reinit_atoms(pme, atomCount, charges.data(), nullptr);
203 /* With decomposition there would be more boilerplate atc code here, e.g. do_redist_pos_coeffs */
207 // TODO: Avoid use of atc in the GPU code path
208 atc = &(pme->atc[0]);
209 // We need to set atc->n for passing the size in the tests
210 atc->setNumAtoms(atomCount);
211 gmx_pme_reinit_atoms(pme, atomCount, charges.data(), nullptr);
213 stateGpu->reinit(atomCount, atomCount);
214 stateGpu->copyCoordinatesToGpu(arrayRefFromArray(coordinates.data(), coordinates.size()),
215 gmx::AtomLocality::All);
216 pme_gpu_set_kernelparam_coordinates(pme->gpu, stateGpu->getCoordinates());
220 default: GMX_THROW(InternalError("Test not implemented for this mode"));
224 //! Getting local PME real grid pointer for test I/O
225 static real* pmeGetRealGridInternal(const gmx_pme_t* pme)
227 const size_t gridIndex = 0;
228 return pme->fftgrid[gridIndex];
231 //! Getting local PME real grid dimensions
232 static void pmeGetRealGridSizesInternal(const gmx_pme_t* pme,
234 IVec& gridSize, //NOLINT(google-runtime-references)
235 IVec& paddedGridSize) //NOLINT(google-runtime-references)
237 const size_t gridIndex = 0;
238 IVec gridOffsetUnused;
242 gmx_parallel_3dfft_real_limits(
243 pme->pfft_setup[gridIndex], gridSize, gridOffsetUnused, paddedGridSize);
247 pme_gpu_get_real_grid_sizes(pme->gpu, &gridSize, &paddedGridSize);
250 default: GMX_THROW(InternalError("Test not implemented for this mode"));
254 //! Getting local PME complex grid pointer for test I/O
255 static t_complex* pmeGetComplexGridInternal(const gmx_pme_t* pme)
257 const size_t gridIndex = 0;
258 return pme->cfftgrid[gridIndex];
261 //! Getting local PME complex grid dimensions
262 static void pmeGetComplexGridSizesInternal(const gmx_pme_t* pme,
263 IVec& gridSize, //NOLINT(google-runtime-references)
264 IVec& paddedGridSize) //NOLINT(google-runtime-references)
266 const size_t gridIndex = 0;
267 IVec gridOffsetUnused, complexOrderUnused;
268 gmx_parallel_3dfft_complex_limits(
269 pme->pfft_setup[gridIndex], complexOrderUnused, gridSize, gridOffsetUnused, paddedGridSize); // TODO: what about YZX ordering?
272 //! Getting the PME grid memory buffer and its sizes - template definition
273 template<typename ValueType>
274 static void pmeGetGridAndSizesInternal(const gmx_pme_t* /*unused*/,
276 ValueType*& /*unused*/, //NOLINT(google-runtime-references)
277 IVec& /*unused*/, //NOLINT(google-runtime-references)
278 IVec& /*unused*/) //NOLINT(google-runtime-references)
280 GMX_THROW(InternalError("Deleted function call"));
281 // explicitly deleting general template does not compile in clang, see https://llvm.org/bugs/show_bug.cgi?id=17537
284 //! Getting the PME real grid memory buffer and its sizes
286 void pmeGetGridAndSizesInternal<real>(const gmx_pme_t* pme, CodePath mode, real*& grid, IVec& gridSize, IVec& paddedGridSize)
288 grid = pmeGetRealGridInternal(pme);
289 pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSize);
292 //! Getting the PME complex grid memory buffer and its sizes
294 void pmeGetGridAndSizesInternal<t_complex>(const gmx_pme_t* pme,
298 IVec& paddedGridSize)
300 grid = pmeGetComplexGridInternal(pme);
301 pmeGetComplexGridSizesInternal(pme, gridSize, paddedGridSize);
304 //! PME spline calculation and charge spreading
305 void pmePerformSplineAndSpread(gmx_pme_t* pme,
306 CodePath mode, // TODO const qualifiers elsewhere
310 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
311 PmeAtomComm* atc = &(pme->atc[0]);
312 const size_t gridIndex = 0;
313 const bool computeSplinesForZeroCharges = true;
314 real** fftgrid = spreadCharges ? pme->fftgrid : nullptr;
315 real* pmegrid = pme->pmegrid[gridIndex].grid.grid;
322 &pme->pmegrid[gridIndex],
325 fftgrid != nullptr ? fftgrid[gridIndex] : nullptr,
326 computeSplinesForZeroCharges,
328 if (spreadCharges && !pme->bUseThreads)
330 wrap_periodic_pmegrid(pme, pmegrid);
331 copy_pmegrid_to_fftgrid(
332 pme, pmegrid, fftgrid != nullptr ? fftgrid[gridIndex] : nullptr, gridIndex);
336 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
337 * double precision. GPUs are not used with double precision anyhow. */
341 const real lambdaQ = 1.0;
342 // no synchronization needed as x is transferred in the PME stream
343 GpuEventSynchronizer* xReadyOnDevice = nullptr;
344 pme_gpu_spread(pme->gpu, xReadyOnDevice, fftgrid, computeSplines, spreadCharges, lambdaQ);
349 default: GMX_THROW(InternalError("Test not implemented for this mode"));
353 //! Getting the internal spline data buffer pointer
354 static real* pmeGetSplineDataInternal(const gmx_pme_t* pme, PmeSplineDataType type, int dimIndex)
356 GMX_ASSERT((0 <= dimIndex) && (dimIndex < DIM), "Invalid dimension index");
357 const PmeAtomComm* atc = &(pme->atc[0]);
358 const size_t threadIndex = 0;
359 real* splineBuffer = nullptr;
362 case PmeSplineDataType::Values:
363 splineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
366 case PmeSplineDataType::Derivatives:
367 splineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
370 default: GMX_THROW(InternalError("Unknown spline data type"));
376 void pmePerformSolve(const gmx_pme_t* pme,
378 PmeSolveAlgorithm method,
380 GridOrdering gridOrdering,
381 bool computeEnergyAndVirial)
383 t_complex* h_grid = pmeGetComplexGridInternal(pme);
384 const bool useLorentzBerthelot = false;
385 const size_t threadIndex = 0;
386 const size_t gridIndex = 0;
390 if (gridOrdering != GridOrdering::YZX)
392 GMX_THROW(InternalError("Test not implemented for this mode"));
396 case PmeSolveAlgorithm::Coulomb:
397 solve_pme_yzx(pme, h_grid, cellVolume, computeEnergyAndVirial, pme->nthread, threadIndex);
400 case PmeSolveAlgorithm::LennardJones:
401 solve_pme_lj_yzx(pme,
405 computeEnergyAndVirial,
410 default: GMX_THROW(InternalError("Test not implemented for this mode"));
417 case PmeSolveAlgorithm::Coulomb:
418 pme_gpu_solve(pme->gpu, gridIndex, h_grid, gridOrdering, computeEnergyAndVirial);
421 default: GMX_THROW(InternalError("Test not implemented for this mode"));
425 default: GMX_THROW(InternalError("Test not implemented for this mode"));
429 //! PME force gathering
430 void pmePerformGather(gmx_pme_t* pme, CodePath mode, ForcesVector& forces)
432 PmeAtomComm* atc = &(pme->atc[0]);
433 const index atomCount = atc->numAtoms();
434 GMX_RELEASE_ASSERT(forces.ssize() == atomCount, "Invalid force buffer size");
435 const real scale = 1.0;
436 const size_t threadIndex = 0;
437 const size_t gridIndex = 0;
438 real* pmegrid = pme->pmegrid[gridIndex].grid.grid;
439 real** fftgrid = pme->fftgrid;
445 if (atc->nthread == 1)
447 // something which is normally done in serial spline computation (make_thread_local_ind())
448 atc->spline[threadIndex].n = atomCount;
450 copy_fftgrid_to_pmegrid(pme, fftgrid[gridIndex], pmegrid, gridIndex, pme->nthread, threadIndex);
451 unwrap_periodic_pmegrid(pme, pmegrid);
452 gather_f_bsplines(pme, pmegrid, true, atc, &atc->spline[threadIndex], scale);
455 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
456 * double precision. GPUs are not used with double precision anyhow. */
460 // Variable initialization needs a non-switch scope
461 const bool computeEnergyAndVirial = false;
462 const real lambdaQ = 1.0;
463 PmeOutput output = pme_gpu_getOutput(*pme, computeEnergyAndVirial, lambdaQ);
464 GMX_ASSERT(forces.size() == output.forces_.size(),
465 "Size of force buffers did not match");
466 pme_gpu_gather(pme->gpu, fftgrid, lambdaQ);
467 std::copy(std::begin(output.forces_), std::end(output.forces_), std::begin(forces));
472 default: GMX_THROW(InternalError("Test not implemented for this mode"));
476 //! PME test finalization before fetching the outputs
477 void pmeFinalizeTest(const gmx_pme_t* pme, CodePath mode)
481 case CodePath::CPU: break;
483 case CodePath::GPU: pme_gpu_synchronize(pme->gpu); break;
485 default: GMX_THROW(InternalError("Test not implemented for this mode"));
489 //! A binary enum for spline data layout transformation
490 enum class PmeLayoutTransform
496 /*! \brief Gets a unique index to an element in a spline parameter buffer.
498 * These theta/dtheta buffers are laid out for GPU spread/gather
499 * kernels. The index is wrt the execution block, in range(0,
500 * atomsPerBlock * order * DIM).
502 * This is a wrapper, only used in unit tests.
503 * \param[in] order PME order
504 * \param[in] splineIndex Spline contribution index (from 0 to \p order - 1)
505 * \param[in] dimIndex Dimension index (from 0 to 2)
506 * \param[in] atomIndex Atom index wrt the block.
507 * \param[in] atomsPerWarp Number of atoms processed by a warp.
509 * \returns Index into theta or dtheta array using GPU layout.
511 static int getSplineParamFullIndex(int order, int splineIndex, int dimIndex, int atomIndex, int atomsPerWarp)
513 if (order != c_pmeGpuOrder)
517 constexpr int fixedOrder = c_pmeGpuOrder;
518 GMX_UNUSED_VALUE(fixedOrder);
520 const int atomWarpIndex = atomIndex % atomsPerWarp;
521 const int warpIndex = atomIndex / atomsPerWarp;
522 int indexBase, result;
523 switch (atomsPerWarp)
526 indexBase = getSplineParamIndexBase<fixedOrder, 1>(warpIndex, atomWarpIndex);
527 result = getSplineParamIndex<fixedOrder, 1>(indexBase, dimIndex, splineIndex);
531 indexBase = getSplineParamIndexBase<fixedOrder, 2>(warpIndex, atomWarpIndex);
532 result = getSplineParamIndex<fixedOrder, 2>(indexBase, dimIndex, splineIndex);
536 indexBase = getSplineParamIndexBase<fixedOrder, 4>(warpIndex, atomWarpIndex);
537 result = getSplineParamIndex<fixedOrder, 4>(indexBase, dimIndex, splineIndex);
541 indexBase = getSplineParamIndexBase<fixedOrder, 8>(warpIndex, atomWarpIndex);
542 result = getSplineParamIndex<fixedOrder, 8>(indexBase, dimIndex, splineIndex);
546 GMX_THROW(NotImplementedError(
547 formatString("Test function call not unrolled for atomsPerWarp = %d in "
548 "getSplineParamFullIndex",
554 /*!\brief Return the number of atoms per warp */
555 static int pme_gpu_get_atoms_per_warp(const PmeGpu* pmeGpu)
557 const int order = pmeGpu->common->pme_order;
558 const int threadsPerAtom =
559 (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
560 return pmeGpu->programHandle_->warpSize() / threadsPerAtom;
563 /*! \brief Rearranges the atom spline data between the GPU and host layouts.
564 * Only used for test purposes so far, likely to be horribly slow.
566 * \param[in] pmeGpu The PME GPU structure.
567 * \param[out] atc The PME CPU atom data structure (with a single-threaded layout).
568 * \param[in] type The spline data type (values or derivatives).
569 * \param[in] dimIndex Dimension index.
570 * \param[in] transform Layout transform type
572 static void pme_gpu_transform_spline_atom_data(const PmeGpu* pmeGpu,
573 const PmeAtomComm* atc,
574 PmeSplineDataType type,
576 PmeLayoutTransform transform)
578 // The GPU atom spline data is laid out in a different way currently than the CPU one.
579 // This function converts the data from GPU to CPU layout (in the host memory).
580 // It is only intended for testing purposes so far.
581 // Ideally we should use similar layouts on CPU and GPU if we care about mixed modes and their
582 // performance (e.g. spreading on GPU, gathering on CPU).
583 GMX_RELEASE_ASSERT(atc->nthread == 1, "Only the serial PME data layout is supported");
584 const uintmax_t threadIndex = 0;
585 const auto atomCount = atc->numAtoms();
586 const auto atomsPerWarp = pme_gpu_get_atoms_per_warp(pmeGpu);
587 const auto pmeOrder = pmeGpu->common->pme_order;
588 GMX_ASSERT(pmeOrder == c_pmeGpuOrder, "Only PME order 4 is implemented");
590 real* cpuSplineBuffer;
591 float* h_splineBuffer;
594 case PmeSplineDataType::Values:
595 cpuSplineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
596 h_splineBuffer = pmeGpu->staging.h_theta;
599 case PmeSplineDataType::Derivatives:
600 cpuSplineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
601 h_splineBuffer = pmeGpu->staging.h_dtheta;
604 default: GMX_THROW(InternalError("Unknown spline data type"));
607 for (auto atomIndex = 0; atomIndex < atomCount; atomIndex++)
609 for (auto orderIndex = 0; orderIndex < pmeOrder; orderIndex++)
611 const auto gpuValueIndex =
612 getSplineParamFullIndex(pmeOrder, orderIndex, dimIndex, atomIndex, atomsPerWarp);
613 const auto cpuValueIndex = atomIndex * pmeOrder + orderIndex;
614 GMX_ASSERT(cpuValueIndex < atomCount * pmeOrder,
615 "Atom spline data index out of bounds (while transforming GPU data layout "
619 case PmeLayoutTransform::GpuToHost:
620 cpuSplineBuffer[cpuValueIndex] = h_splineBuffer[gpuValueIndex];
623 case PmeLayoutTransform::HostToGpu:
624 h_splineBuffer[gpuValueIndex] = cpuSplineBuffer[cpuValueIndex];
627 default: GMX_THROW(InternalError("Unknown layout transform"));
633 //! Setting atom spline values/derivatives to be used in spread/gather
634 void pmeSetSplineData(const gmx_pme_t* pme,
636 const SplineParamsDimVector& splineValues,
637 PmeSplineDataType type,
640 const PmeAtomComm* atc = &(pme->atc[0]);
641 const index atomCount = atc->numAtoms();
642 const index pmeOrder = pme->pme_order;
643 const index dimSize = pmeOrder * atomCount;
644 GMX_RELEASE_ASSERT(dimSize == splineValues.ssize(), "Mismatch in spline data");
645 real* splineBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
650 std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
654 std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
655 pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::HostToGpu);
658 default: GMX_THROW(InternalError("Test not implemented for this mode"));
662 //! Setting gridline indices to be used in spread/gather
663 void pmeSetGridLineIndices(gmx_pme_t* pme, CodePath mode, const GridLineIndicesVector& gridLineIndices)
665 PmeAtomComm* atc = &(pme->atc[0]);
666 const index atomCount = atc->numAtoms();
667 GMX_RELEASE_ASSERT(atomCount == gridLineIndices.ssize(), "Mismatch in gridline indices size");
669 IVec paddedGridSizeUnused, gridSize(0, 0, 0);
670 pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSizeUnused);
672 for (const auto& index : gridLineIndices)
674 for (int i = 0; i < DIM; i++)
676 GMX_RELEASE_ASSERT((0 <= index[i]) && (index[i] < gridSize[i]),
677 "Invalid gridline index");
684 memcpy(pme_gpu_staging(pme->gpu).h_gridlineIndices,
685 gridLineIndices.data(),
686 atomCount * sizeof(gridLineIndices[0]));
690 atc->idx.resize(gridLineIndices.size());
691 std::copy(gridLineIndices.begin(), gridLineIndices.end(), atc->idx.begin());
693 default: GMX_THROW(InternalError("Test not implemented for this mode"));
697 //! Getting plain index into the complex 3d grid
698 inline size_t pmeGetGridPlainIndexInternal(const IVec& index, const IVec& paddedGridSize, GridOrdering gridOrdering)
701 switch (gridOrdering)
703 case GridOrdering::YZX:
704 result = (index[YY] * paddedGridSize[ZZ] + index[ZZ]) * paddedGridSize[XX] + index[XX];
707 case GridOrdering::XYZ:
708 result = (index[XX] * paddedGridSize[YY] + index[YY]) * paddedGridSize[ZZ] + index[ZZ];
711 default: GMX_THROW(InternalError("Test not implemented for this mode"));
716 //! Setting real or complex grid
717 template<typename ValueType>
718 static void pmeSetGridInternal(const gmx_pme_t* pme,
720 GridOrdering gridOrdering,
721 const SparseGridValuesInput<ValueType>& gridValues)
723 IVec gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
725 pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
729 case CodePath::GPU: // intentional absence of break, the grid will be copied from the host buffer in testing mode
731 std::memset(grid, 0, paddedGridSize[XX] * paddedGridSize[YY] * paddedGridSize[ZZ] * sizeof(ValueType));
732 for (const auto& gridValue : gridValues)
734 for (int i = 0; i < DIM; i++)
736 GMX_RELEASE_ASSERT((0 <= gridValue.first[i]) && (gridValue.first[i] < gridSize[i]),
737 "Invalid grid value index");
739 const size_t gridValueIndex =
740 pmeGetGridPlainIndexInternal(gridValue.first, paddedGridSize, gridOrdering);
741 grid[gridValueIndex] = gridValue.second;
745 default: GMX_THROW(InternalError("Test not implemented for this mode"));
749 //! Setting real grid to be used in gather
750 void pmeSetRealGrid(const gmx_pme_t* pme, CodePath mode, const SparseRealGridValuesInput& gridValues)
752 pmeSetGridInternal<real>(pme, mode, GridOrdering::XYZ, gridValues);
755 //! Setting complex grid to be used in solve
756 void pmeSetComplexGrid(const gmx_pme_t* pme,
758 GridOrdering gridOrdering,
759 const SparseComplexGridValuesInput& gridValues)
761 pmeSetGridInternal<t_complex>(pme, mode, gridOrdering, gridValues);
764 //! Getting the single dimension's spline values or derivatives
765 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t* pme, CodePath mode, PmeSplineDataType type, int dimIndex)
767 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
768 const PmeAtomComm* atc = &(pme->atc[0]);
769 const size_t atomCount = atc->numAtoms();
770 const size_t pmeOrder = pme->pme_order;
771 const size_t dimSize = pmeOrder * atomCount;
773 real* sourceBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
774 SplineParamsDimVector result;
778 pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::GpuToHost);
779 result = arrayRefFromArray(sourceBuffer, dimSize);
782 case CodePath::CPU: result = arrayRefFromArray(sourceBuffer, dimSize); break;
784 default: GMX_THROW(InternalError("Test not implemented for this mode"));
789 //! Getting the gridline indices
790 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t* pme, CodePath mode)
792 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
793 const PmeAtomComm* atc = &(pme->atc[0]);
794 const size_t atomCount = atc->numAtoms();
796 GridLineIndicesVector gridLineIndices;
800 gridLineIndices = arrayRefFromArray(
801 reinterpret_cast<IVec*>(pme_gpu_staging(pme->gpu).h_gridlineIndices), atomCount);
804 case CodePath::CPU: gridLineIndices = atc->idx; break;
806 default: GMX_THROW(InternalError("Test not implemented for this mode"));
808 return gridLineIndices;
811 //! Getting real or complex grid - only non zero values
812 template<typename ValueType>
813 static SparseGridValuesOutput<ValueType> pmeGetGridInternal(const gmx_pme_t* pme,
815 GridOrdering gridOrdering)
817 IVec gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
819 pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
820 SparseGridValuesOutput<ValueType> gridValues;
823 case CodePath::GPU: // intentional absence of break
826 for (int ix = 0; ix < gridSize[XX]; ix++)
828 for (int iy = 0; iy < gridSize[YY]; iy++)
830 for (int iz = 0; iz < gridSize[ZZ]; iz++)
832 IVec temp(ix, iy, iz);
833 const size_t gridValueIndex =
834 pmeGetGridPlainIndexInternal(temp, paddedGridSize, gridOrdering);
835 const ValueType value = grid[gridValueIndex];
836 if (value != ValueType{})
838 auto key = formatString("Cell %d %d %d", ix, iy, iz);
839 gridValues[key] = value;
846 default: GMX_THROW(InternalError("Test not implemented for this mode"));
851 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
852 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t* pme, CodePath mode)
854 return pmeGetGridInternal<real>(pme, mode, GridOrdering::XYZ);
857 //! Getting the complex grid output of pmePerformSolve()
858 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t* pme, CodePath mode, GridOrdering gridOrdering)
860 return pmeGetGridInternal<t_complex>(pme, mode, gridOrdering);
863 //! Getting the reciprocal energy and virial
864 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t* pme, CodePath mode, PmeSolveAlgorithm method)
867 const real lambdaQ = 1.0;
873 case PmeSolveAlgorithm::Coulomb:
874 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output);
877 case PmeSolveAlgorithm::LennardJones:
878 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output);
881 default: GMX_THROW(InternalError("Test not implemented for this mode"));
887 case PmeSolveAlgorithm::Coulomb:
888 pme_gpu_getEnergyAndVirial(*pme, lambdaQ, &output);
891 default: GMX_THROW(InternalError("Test not implemented for this mode"));
895 default: GMX_THROW(InternalError("Test not implemented for this mode"));
900 const char* codePathToString(CodePath codePath)
904 case CodePath::CPU: return "CPU";
905 case CodePath::GPU: return "GPU";
906 default: GMX_THROW(NotImplementedError("This CodePath should support codePathToString"));
910 PmeTestHardwareContext::PmeTestHardwareContext() : codePath_(CodePath::CPU) {}
912 PmeTestHardwareContext::PmeTestHardwareContext(TestDevice* testDevice) :
913 codePath_(CodePath::CPU),
914 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;