2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019,2020, 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/device_context.h"
63 #include "gromacs/gpu_utils/gpu_utils.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"
81 bool pmeSupportsInputForMode(const gmx_hw_info_t& hwinfo, const t_inputrec* inputRec, CodePath mode)
86 case CodePath::CPU: implemented = true; break;
89 implemented = (pme_gpu_supports_build(nullptr) && pme_gpu_supports_hardware(hwinfo, nullptr)
90 && pme_gpu_supports_input(*inputRec, nullptr));
93 default: GMX_THROW(InternalError("Test not implemented for this mode"));
98 uint64_t getSplineModuliDoublePrecisionUlps(int splineOrder)
100 /* Arbitrary ulp tolerance for sine/cosine implementation. It's
101 * hard to know what to pick without testing lots of
102 * implementations. */
103 const uint64_t sineUlps = 10;
104 return 4 * (splineOrder - 2) + 2 * sineUlps * splineOrder;
107 //! PME initialization
108 PmeSafePointer pmeInitWrapper(const t_inputrec* inputRec,
110 const DeviceContext* deviceContext,
111 const DeviceStream* deviceStream,
112 const PmeGpuProgram* pmeGpuProgram,
113 const Matrix3x3& box,
114 const real ewaldCoeff_q,
115 const real ewaldCoeff_lj)
117 const MDLogger dummyLogger;
118 const auto runMode = (mode == CodePath::CPU) ? PmeRunMode::CPU : PmeRunMode::Mixed;
119 t_commrec dummyCommrec = { 0 };
120 NumPmeDomains numPmeDomains = { 1, 1 };
121 gmx_pme_t* pmeDataRaw = gmx_pme_init(&dummyCommrec, numPmeDomains, inputRec, false, false, true,
122 ewaldCoeff_q, ewaldCoeff_lj, 1, runMode, nullptr,
123 deviceContext, deviceStream, pmeGpuProgram, dummyLogger);
124 PmeSafePointer pme(pmeDataRaw); // taking ownership
126 // TODO get rid of this with proper matrix type
128 for (int i = 0; i < DIM; i++)
130 for (int j = 0; j < DIM; j++)
132 boxTemp[i][j] = box[i * DIM + j];
135 const char* boxError = check_box(PbcType::Unset, boxTemp);
136 GMX_RELEASE_ASSERT(boxError == nullptr, boxError);
140 case CodePath::CPU: invertBoxMatrix(boxTemp, pme->recipbox); break;
143 pme_gpu_set_testing(pme->gpu, true);
144 pme_gpu_update_input_box(pme->gpu, boxTemp);
147 default: GMX_THROW(InternalError("Test not implemented for this mode"));
153 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec)
155 const Matrix3x3 defaultBox = { { 1.0F, 0.0F, 0.0F, 0.0F, 1.0F, 0.0F, 0.0F, 0.0F, 1.0F } };
156 return pmeInitWrapper(inputRec, CodePath::CPU, nullptr, nullptr, nullptr, defaultBox, 0.0F, 0.0F);
159 //! Make a GPU state-propagator manager
160 std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme_t& pme,
161 const DeviceContext* deviceContext,
162 const DeviceStream* deviceStream)
164 // TODO: Pin the host buffer and use async memory copies
165 // TODO: Special constructor for PME-only rank / PME-tests is used here. There should be a mechanism to
166 // restrict one from using other constructor here.
167 return std::make_unique<StatePropagatorDataGpu>(deviceStream, *deviceContext, GpuApiCallBehavior::Sync,
168 pme_gpu_get_block_size(&pme), nullptr);
171 //! PME initialization with atom data
172 void pmeInitAtoms(gmx_pme_t* pme,
173 StatePropagatorDataGpu* stateGpu,
175 const CoordinatesVector& coordinates,
176 const ChargesVector& charges)
178 const index atomCount = coordinates.size();
179 GMX_RELEASE_ASSERT(atomCount == charges.ssize(), "Mismatch in atom data");
180 PmeAtomComm* atc = nullptr;
185 atc = &(pme->atc[0]);
186 atc->x = coordinates;
187 atc->coefficient = charges;
188 gmx_pme_reinit_atoms(pme, atomCount, charges.data(), nullptr);
189 /* With decomposition there would be more boilerplate atc code here, e.g. do_redist_pos_coeffs */
193 // TODO: Avoid use of atc in the GPU code path
194 atc = &(pme->atc[0]);
195 // We need to set atc->n for passing the size in the tests
196 atc->setNumAtoms(atomCount);
197 gmx_pme_reinit_atoms(pme, atomCount, charges.data(), nullptr);
199 stateGpu->reinit(atomCount, atomCount);
200 stateGpu->copyCoordinatesToGpu(arrayRefFromArray(coordinates.data(), coordinates.size()),
201 gmx::AtomLocality::All);
202 pme_gpu_set_kernelparam_coordinates(pme->gpu, stateGpu->getCoordinates());
206 default: GMX_THROW(InternalError("Test not implemented for this mode"));
210 //! Getting local PME real grid pointer for test I/O
211 static real* pmeGetRealGridInternal(const gmx_pme_t* pme)
213 const size_t gridIndex = 0;
214 return pme->fftgrid[gridIndex];
217 //! Getting local PME real grid dimensions
218 static void pmeGetRealGridSizesInternal(const gmx_pme_t* pme,
220 IVec& gridSize, //NOLINT(google-runtime-references)
221 IVec& paddedGridSize) //NOLINT(google-runtime-references)
223 const size_t gridIndex = 0;
224 IVec gridOffsetUnused;
228 gmx_parallel_3dfft_real_limits(pme->pfft_setup[gridIndex], gridSize, gridOffsetUnused,
233 pme_gpu_get_real_grid_sizes(pme->gpu, &gridSize, &paddedGridSize);
236 default: GMX_THROW(InternalError("Test not implemented for this mode"));
240 //! Getting local PME complex grid pointer for test I/O
241 static t_complex* pmeGetComplexGridInternal(const gmx_pme_t* pme)
243 const size_t gridIndex = 0;
244 return pme->cfftgrid[gridIndex];
247 //! Getting local PME complex grid dimensions
248 static void pmeGetComplexGridSizesInternal(const gmx_pme_t* pme,
249 IVec& gridSize, //NOLINT(google-runtime-references)
250 IVec& paddedGridSize) //NOLINT(google-runtime-references)
252 const size_t gridIndex = 0;
253 IVec gridOffsetUnused, complexOrderUnused;
254 gmx_parallel_3dfft_complex_limits(pme->pfft_setup[gridIndex], complexOrderUnused, gridSize,
255 gridOffsetUnused, paddedGridSize); // TODO: what about YZX ordering?
258 //! Getting the PME grid memory buffer and its sizes - template definition
259 template<typename ValueType>
260 static void pmeGetGridAndSizesInternal(const gmx_pme_t* /*unused*/,
262 ValueType*& /*unused*/, //NOLINT(google-runtime-references)
263 IVec& /*unused*/, //NOLINT(google-runtime-references)
264 IVec& /*unused*/) //NOLINT(google-runtime-references)
266 GMX_THROW(InternalError("Deleted function call"));
267 // explicitly deleting general template does not compile in clang/icc, see https://llvm.org/bugs/show_bug.cgi?id=17537
270 //! Getting the PME real grid memory buffer and its sizes
272 void pmeGetGridAndSizesInternal<real>(const gmx_pme_t* pme, CodePath mode, real*& grid, IVec& gridSize, IVec& paddedGridSize)
274 grid = pmeGetRealGridInternal(pme);
275 pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSize);
278 //! Getting the PME complex grid memory buffer and its sizes
280 void pmeGetGridAndSizesInternal<t_complex>(const gmx_pme_t* pme,
284 IVec& paddedGridSize)
286 grid = pmeGetComplexGridInternal(pme);
287 pmeGetComplexGridSizesInternal(pme, gridSize, paddedGridSize);
290 //! PME spline calculation and charge spreading
291 void pmePerformSplineAndSpread(gmx_pme_t* pme,
292 CodePath mode, // TODO const qualifiers elsewhere
296 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
297 PmeAtomComm* atc = &(pme->atc[0]);
298 const size_t gridIndex = 0;
299 const bool computeSplinesForZeroCharges = true;
300 real** fftgrid = spreadCharges ? pme->fftgrid : nullptr;
301 real* pmegrid = pme->pmegrid[gridIndex].grid.grid;
306 spread_on_grid(pme, atc, &pme->pmegrid[gridIndex], computeSplines, spreadCharges,
307 fftgrid != nullptr ? fftgrid[gridIndex] : nullptr,
308 computeSplinesForZeroCharges, gridIndex);
309 if (spreadCharges && !pme->bUseThreads)
311 wrap_periodic_pmegrid(pme, pmegrid);
312 copy_pmegrid_to_fftgrid(
313 pme, pmegrid, fftgrid != nullptr ? fftgrid[gridIndex] : nullptr, gridIndex);
317 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
318 * double precision. GPUs are not used with double precision anyhow. */
322 const real lambdaQ = 1.0;
323 // no synchronization needed as x is transferred in the PME stream
324 GpuEventSynchronizer* xReadyOnDevice = nullptr;
325 pme_gpu_spread(pme->gpu, xReadyOnDevice, fftgrid, computeSplines, spreadCharges, lambdaQ);
330 default: GMX_THROW(InternalError("Test not implemented for this mode"));
334 //! Getting the internal spline data buffer pointer
335 static real* pmeGetSplineDataInternal(const gmx_pme_t* pme, PmeSplineDataType type, int dimIndex)
337 GMX_ASSERT((0 <= dimIndex) && (dimIndex < DIM), "Invalid dimension index");
338 const PmeAtomComm* atc = &(pme->atc[0]);
339 const size_t threadIndex = 0;
340 real* splineBuffer = nullptr;
343 case PmeSplineDataType::Values:
344 splineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
347 case PmeSplineDataType::Derivatives:
348 splineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
351 default: GMX_THROW(InternalError("Unknown spline data type"));
357 void pmePerformSolve(const gmx_pme_t* pme,
359 PmeSolveAlgorithm method,
361 GridOrdering gridOrdering,
362 bool computeEnergyAndVirial)
364 t_complex* h_grid = pmeGetComplexGridInternal(pme);
365 const bool useLorentzBerthelot = false;
366 const size_t threadIndex = 0;
367 const size_t gridIndex = 0;
371 if (gridOrdering != GridOrdering::YZX)
373 GMX_THROW(InternalError("Test not implemented for this mode"));
377 case PmeSolveAlgorithm::Coulomb:
378 solve_pme_yzx(pme, h_grid, cellVolume, computeEnergyAndVirial, pme->nthread, threadIndex);
381 case PmeSolveAlgorithm::LennardJones:
382 solve_pme_lj_yzx(pme, &h_grid, useLorentzBerthelot, cellVolume,
383 computeEnergyAndVirial, pme->nthread, threadIndex);
386 default: GMX_THROW(InternalError("Test not implemented for this mode"));
393 case PmeSolveAlgorithm::Coulomb:
394 pme_gpu_solve(pme->gpu, gridIndex, h_grid, gridOrdering, computeEnergyAndVirial);
397 default: GMX_THROW(InternalError("Test not implemented for this mode"));
401 default: GMX_THROW(InternalError("Test not implemented for this mode"));
405 //! PME force gathering
406 void pmePerformGather(gmx_pme_t* pme, CodePath mode, ForcesVector& forces)
408 PmeAtomComm* atc = &(pme->atc[0]);
409 const index atomCount = atc->numAtoms();
410 GMX_RELEASE_ASSERT(forces.ssize() == atomCount, "Invalid force buffer size");
411 const real scale = 1.0;
412 const size_t threadIndex = 0;
413 const size_t gridIndex = 0;
414 real* pmegrid = pme->pmegrid[gridIndex].grid.grid;
415 real** fftgrid = pme->fftgrid;
421 if (atc->nthread == 1)
423 // something which is normally done in serial spline computation (make_thread_local_ind())
424 atc->spline[threadIndex].n = atomCount;
426 copy_fftgrid_to_pmegrid(pme, fftgrid[gridIndex], pmegrid, gridIndex, pme->nthread, threadIndex);
427 unwrap_periodic_pmegrid(pme, pmegrid);
428 gather_f_bsplines(pme, pmegrid, true, atc, &atc->spline[threadIndex], scale);
431 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
432 * double precision. GPUs are not used with double precision anyhow. */
436 // Variable initialization needs a non-switch scope
437 const bool computeEnergyAndVirial = false;
438 const real lambdaQ = 1.0;
439 PmeOutput output = pme_gpu_getOutput(*pme, computeEnergyAndVirial, lambdaQ);
440 GMX_ASSERT(forces.size() == output.forces_.size(),
441 "Size of force buffers did not match");
442 pme_gpu_gather(pme->gpu, fftgrid, lambdaQ);
443 std::copy(std::begin(output.forces_), std::end(output.forces_), std::begin(forces));
448 default: GMX_THROW(InternalError("Test not implemented for this mode"));
452 //! PME test finalization before fetching the outputs
453 void pmeFinalizeTest(const gmx_pme_t* pme, CodePath mode)
457 case CodePath::CPU: break;
459 case CodePath::GPU: pme_gpu_synchronize(pme->gpu); break;
461 default: GMX_THROW(InternalError("Test not implemented for this mode"));
465 //! A binary enum for spline data layout transformation
466 enum class PmeLayoutTransform
472 /*! \brief Gets a unique index to an element in a spline parameter buffer.
474 * These theta/dtheta buffers are laid out for GPU spread/gather
475 * kernels. The index is wrt the execution block, in range(0,
476 * atomsPerBlock * order * DIM).
478 * This is a wrapper, only used in unit tests.
479 * \param[in] order PME order
480 * \param[in] splineIndex Spline contribution index (from 0 to \p order - 1)
481 * \param[in] dimIndex Dimension index (from 0 to 2)
482 * \param[in] atomIndex Atom index wrt the block.
483 * \param[in] atomsPerWarp Number of atoms processed by a warp.
485 * \returns Index into theta or dtheta array using GPU layout.
487 static int getSplineParamFullIndex(int order, int splineIndex, int dimIndex, int atomIndex, int atomsPerWarp)
489 if (order != c_pmeGpuOrder)
493 constexpr int fixedOrder = c_pmeGpuOrder;
494 GMX_UNUSED_VALUE(fixedOrder);
496 const int atomWarpIndex = atomIndex % atomsPerWarp;
497 const int warpIndex = atomIndex / atomsPerWarp;
498 int indexBase, result;
499 switch (atomsPerWarp)
502 indexBase = getSplineParamIndexBase<fixedOrder, 1>(warpIndex, atomWarpIndex);
503 result = getSplineParamIndex<fixedOrder, 1>(indexBase, dimIndex, splineIndex);
507 indexBase = getSplineParamIndexBase<fixedOrder, 2>(warpIndex, atomWarpIndex);
508 result = getSplineParamIndex<fixedOrder, 2>(indexBase, dimIndex, splineIndex);
512 indexBase = getSplineParamIndexBase<fixedOrder, 4>(warpIndex, atomWarpIndex);
513 result = getSplineParamIndex<fixedOrder, 4>(indexBase, dimIndex, splineIndex);
517 indexBase = getSplineParamIndexBase<fixedOrder, 8>(warpIndex, atomWarpIndex);
518 result = getSplineParamIndex<fixedOrder, 8>(indexBase, dimIndex, splineIndex);
522 GMX_THROW(NotImplementedError(
523 formatString("Test function call not unrolled for atomsPerWarp = %d in "
524 "getSplineParamFullIndex",
530 /*!\brief Return the number of atoms per warp */
531 static int pme_gpu_get_atoms_per_warp(const PmeGpu* pmeGpu)
533 const int order = pmeGpu->common->pme_order;
534 const int threadsPerAtom =
535 (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
536 return pmeGpu->programHandle_->warpSize() / threadsPerAtom;
539 /*! \brief Rearranges the atom spline data between the GPU and host layouts.
540 * Only used for test purposes so far, likely to be horribly slow.
542 * \param[in] pmeGpu The PME GPU structure.
543 * \param[out] atc The PME CPU atom data structure (with a single-threaded layout).
544 * \param[in] type The spline data type (values or derivatives).
545 * \param[in] dimIndex Dimension index.
546 * \param[in] transform Layout transform type
548 static void pme_gpu_transform_spline_atom_data(const PmeGpu* pmeGpu,
549 const PmeAtomComm* atc,
550 PmeSplineDataType type,
552 PmeLayoutTransform transform)
554 // The GPU atom spline data is laid out in a different way currently than the CPU one.
555 // This function converts the data from GPU to CPU layout (in the host memory).
556 // It is only intended for testing purposes so far.
557 // Ideally we should use similar layouts on CPU and GPU if we care about mixed modes and their
558 // performance (e.g. spreading on GPU, gathering on CPU).
559 GMX_RELEASE_ASSERT(atc->nthread == 1, "Only the serial PME data layout is supported");
560 const uintmax_t threadIndex = 0;
561 const auto atomCount = atc->numAtoms();
562 const auto atomsPerWarp = pme_gpu_get_atoms_per_warp(pmeGpu);
563 const auto pmeOrder = pmeGpu->common->pme_order;
564 GMX_ASSERT(pmeOrder == c_pmeGpuOrder, "Only PME order 4 is implemented");
566 real* cpuSplineBuffer;
567 float* h_splineBuffer;
570 case PmeSplineDataType::Values:
571 cpuSplineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
572 h_splineBuffer = pmeGpu->staging.h_theta;
575 case PmeSplineDataType::Derivatives:
576 cpuSplineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
577 h_splineBuffer = pmeGpu->staging.h_dtheta;
580 default: GMX_THROW(InternalError("Unknown spline data type"));
583 for (auto atomIndex = 0; atomIndex < atomCount; atomIndex++)
585 for (auto orderIndex = 0; orderIndex < pmeOrder; orderIndex++)
587 const auto gpuValueIndex =
588 getSplineParamFullIndex(pmeOrder, orderIndex, dimIndex, atomIndex, atomsPerWarp);
589 const auto cpuValueIndex = atomIndex * pmeOrder + orderIndex;
590 GMX_ASSERT(cpuValueIndex < atomCount * pmeOrder,
591 "Atom spline data index out of bounds (while transforming GPU data layout "
595 case PmeLayoutTransform::GpuToHost:
596 cpuSplineBuffer[cpuValueIndex] = h_splineBuffer[gpuValueIndex];
599 case PmeLayoutTransform::HostToGpu:
600 h_splineBuffer[gpuValueIndex] = cpuSplineBuffer[cpuValueIndex];
603 default: GMX_THROW(InternalError("Unknown layout transform"));
609 //! Setting atom spline values/derivatives to be used in spread/gather
610 void pmeSetSplineData(const gmx_pme_t* pme,
612 const SplineParamsDimVector& splineValues,
613 PmeSplineDataType type,
616 const PmeAtomComm* atc = &(pme->atc[0]);
617 const index atomCount = atc->numAtoms();
618 const index pmeOrder = pme->pme_order;
619 const index dimSize = pmeOrder * atomCount;
620 GMX_RELEASE_ASSERT(dimSize == splineValues.ssize(), "Mismatch in spline data");
621 real* splineBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
626 std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
630 std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
631 pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::HostToGpu);
634 default: GMX_THROW(InternalError("Test not implemented for this mode"));
638 //! Setting gridline indices to be used in spread/gather
639 void pmeSetGridLineIndices(gmx_pme_t* pme, CodePath mode, const GridLineIndicesVector& gridLineIndices)
641 PmeAtomComm* atc = &(pme->atc[0]);
642 const index atomCount = atc->numAtoms();
643 GMX_RELEASE_ASSERT(atomCount == gridLineIndices.ssize(), "Mismatch in gridline indices size");
645 IVec paddedGridSizeUnused, gridSize(0, 0, 0);
646 pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSizeUnused);
648 for (const auto& index : gridLineIndices)
650 for (int i = 0; i < DIM; i++)
652 GMX_RELEASE_ASSERT((0 <= index[i]) && (index[i] < gridSize[i]),
653 "Invalid gridline index");
660 memcpy(pme_gpu_staging(pme->gpu).h_gridlineIndices, gridLineIndices.data(),
661 atomCount * sizeof(gridLineIndices[0]));
665 atc->idx.resize(gridLineIndices.size());
666 std::copy(gridLineIndices.begin(), gridLineIndices.end(), atc->idx.begin());
668 default: GMX_THROW(InternalError("Test not implemented for this mode"));
672 //! Getting plain index into the complex 3d grid
673 inline size_t pmeGetGridPlainIndexInternal(const IVec& index, const IVec& paddedGridSize, GridOrdering gridOrdering)
676 switch (gridOrdering)
678 case GridOrdering::YZX:
679 result = (index[YY] * paddedGridSize[ZZ] + index[ZZ]) * paddedGridSize[XX] + index[XX];
682 case GridOrdering::XYZ:
683 result = (index[XX] * paddedGridSize[YY] + index[YY]) * paddedGridSize[ZZ] + index[ZZ];
686 default: GMX_THROW(InternalError("Test not implemented for this mode"));
691 //! Setting real or complex grid
692 template<typename ValueType>
693 static void pmeSetGridInternal(const gmx_pme_t* pme,
695 GridOrdering gridOrdering,
696 const SparseGridValuesInput<ValueType>& gridValues)
698 IVec gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
700 pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
704 case CodePath::GPU: // intentional absence of break, the grid will be copied from the host buffer in testing mode
707 paddedGridSize[XX] * paddedGridSize[YY] * paddedGridSize[ZZ] * sizeof(ValueType));
708 for (const auto& gridValue : gridValues)
710 for (int i = 0; i < DIM; i++)
712 GMX_RELEASE_ASSERT((0 <= gridValue.first[i]) && (gridValue.first[i] < gridSize[i]),
713 "Invalid grid value index");
715 const size_t gridValueIndex =
716 pmeGetGridPlainIndexInternal(gridValue.first, paddedGridSize, gridOrdering);
717 grid[gridValueIndex] = gridValue.second;
721 default: GMX_THROW(InternalError("Test not implemented for this mode"));
725 //! Setting real grid to be used in gather
726 void pmeSetRealGrid(const gmx_pme_t* pme, CodePath mode, const SparseRealGridValuesInput& gridValues)
728 pmeSetGridInternal<real>(pme, mode, GridOrdering::XYZ, gridValues);
731 //! Setting complex grid to be used in solve
732 void pmeSetComplexGrid(const gmx_pme_t* pme,
734 GridOrdering gridOrdering,
735 const SparseComplexGridValuesInput& gridValues)
737 pmeSetGridInternal<t_complex>(pme, mode, gridOrdering, gridValues);
740 //! Getting the single dimension's spline values or derivatives
741 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t* pme, CodePath mode, PmeSplineDataType type, int dimIndex)
743 GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
744 const PmeAtomComm* atc = &(pme->atc[0]);
745 const size_t atomCount = atc->numAtoms();
746 const size_t pmeOrder = pme->pme_order;
747 const size_t dimSize = pmeOrder * atomCount;
749 real* sourceBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
750 SplineParamsDimVector result;
754 pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::GpuToHost);
755 result = arrayRefFromArray(sourceBuffer, dimSize);
758 case CodePath::CPU: result = arrayRefFromArray(sourceBuffer, dimSize); break;
760 default: GMX_THROW(InternalError("Test not implemented for this mode"));
765 //! Getting the gridline indices
766 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t* pme, CodePath mode)
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();
772 GridLineIndicesVector gridLineIndices;
776 gridLineIndices = arrayRefFromArray(
777 reinterpret_cast<IVec*>(pme_gpu_staging(pme->gpu).h_gridlineIndices), atomCount);
780 case CodePath::CPU: gridLineIndices = atc->idx; break;
782 default: GMX_THROW(InternalError("Test not implemented for this mode"));
784 return gridLineIndices;
787 //! Getting real or complex grid - only non zero values
788 template<typename ValueType>
789 static SparseGridValuesOutput<ValueType> pmeGetGridInternal(const gmx_pme_t* pme,
791 GridOrdering gridOrdering)
793 IVec gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
795 pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
796 SparseGridValuesOutput<ValueType> gridValues;
799 case CodePath::GPU: // intentional absence of break
802 for (int ix = 0; ix < gridSize[XX]; ix++)
804 for (int iy = 0; iy < gridSize[YY]; iy++)
806 for (int iz = 0; iz < gridSize[ZZ]; iz++)
808 IVec temp(ix, iy, iz);
809 const size_t gridValueIndex =
810 pmeGetGridPlainIndexInternal(temp, paddedGridSize, gridOrdering);
811 const ValueType value = grid[gridValueIndex];
812 if (value != ValueType{})
814 auto key = formatString("Cell %d %d %d", ix, iy, iz);
815 gridValues[key] = value;
822 default: GMX_THROW(InternalError("Test not implemented for this mode"));
827 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
828 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t* pme, CodePath mode)
830 return pmeGetGridInternal<real>(pme, mode, GridOrdering::XYZ);
833 //! Getting the complex grid output of pmePerformSolve()
834 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t* pme, CodePath mode, GridOrdering gridOrdering)
836 return pmeGetGridInternal<t_complex>(pme, mode, gridOrdering);
839 //! Getting the reciprocal energy and virial
840 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t* pme, CodePath mode, PmeSolveAlgorithm method)
843 const real lambdaQ = 1.0;
849 case PmeSolveAlgorithm::Coulomb:
850 get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output);
853 case PmeSolveAlgorithm::LennardJones:
854 get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output);
857 default: GMX_THROW(InternalError("Test not implemented for this mode"));
863 case PmeSolveAlgorithm::Coulomb:
864 pme_gpu_getEnergyAndVirial(*pme, lambdaQ, &output);
867 default: GMX_THROW(InternalError("Test not implemented for this mode"));
871 default: GMX_THROW(InternalError("Test not implemented for this mode"));
876 const char* codePathToString(CodePath codePath)
880 case CodePath::CPU: return "CPU";
881 case CodePath::GPU: return "GPU";
882 default: GMX_THROW(NotImplementedError("This CodePath should support codePathToString"));
886 PmeTestHardwareContext::PmeTestHardwareContext() : codePath_(CodePath::CPU) {}
888 PmeTestHardwareContext::PmeTestHardwareContext(TestDevice* testDevice) :
889 codePath_(CodePath::CPU),
890 testDevice_(testDevice)
892 setActiveDevice(testDevice_->deviceInfo());
893 pmeGpuProgram_ = buildPmeGpuProgram(testDevice_->deviceContext());
896 //! Returns a human-readable context description line
897 std::string PmeTestHardwareContext::description() const
901 case CodePath::CPU: return "CPU";
902 case CodePath::GPU: return "GPU (" + testDevice_->description() + ")";
903 default: return "Unknown code path.";
907 void PmeTestHardwareContext::activate() const
909 if (codePath_ == CodePath::GPU)
911 setActiveDevice(testDevice_->deviceInfo());
915 std::vector<std::unique_ptr<PmeTestHardwareContext>> createPmeTestHardwareContextList()
917 std::vector<std::unique_ptr<PmeTestHardwareContext>> pmeTestHardwareContextList;
919 pmeTestHardwareContextList.emplace_back(std::make_unique<PmeTestHardwareContext>());
921 const auto& testDeviceList = getTestHardwareEnvironment()->getTestDeviceList();
922 for (const auto& testDevice : testDeviceList)
924 pmeTestHardwareContextList.emplace_back(std::make_unique<PmeTestHardwareContext>(testDevice.get()));
926 return pmeTestHardwareContextList;