eb1e170411bbf090582852ed283fe65824919143
[alexxy/gromacs.git] / src / gromacs / ewald / tests / pmetestcommon.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements common routines for PME tests.
38  *
39  * \author Aleksei Iupinov <a.yupinov@gmail.com>
40  * \ingroup module_ewald
41  */
42 #include "gmxpre.h"
43
44 #include "pmetestcommon.h"
45
46 #include <cstring>
47
48 #include <algorithm>
49
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"
72
73 #include "testutils/test_hardware_environment.h"
74 #include "testutils/testasserts.h"
75
76 class DeviceContext;
77
78 namespace gmx
79 {
80 namespace test
81 {
82
83 bool pmeSupportsInputForMode(const gmx_hw_info_t& hwinfo, const t_inputrec* inputRec, CodePath mode)
84 {
85     bool implemented;
86     switch (mode)
87     {
88         case CodePath::CPU: implemented = true; break;
89
90         case CodePath::GPU:
91             implemented = (pme_gpu_supports_build(nullptr) && pme_gpu_supports_hardware(hwinfo, nullptr)
92                            && pme_gpu_supports_input(*inputRec, nullptr));
93             break;
94
95         default: GMX_THROW(InternalError("Test not implemented for this mode"));
96     }
97     return implemented;
98 }
99
100 uint64_t getSplineModuliDoublePrecisionUlps(int splineOrder)
101 {
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;
107 }
108
109 //! PME initialization
110 PmeSafePointer pmeInitWrapper(const t_inputrec*    inputRec,
111                               const CodePath       mode,
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)
118 {
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, numPmeDomains, inputRec, false, false, true,
124                                          ewaldCoeff_q, ewaldCoeff_lj, 1, runMode, nullptr,
125                                          deviceContext, deviceStream, pmeGpuProgram, dummyLogger);
126     PmeSafePointer pme(pmeDataRaw); // taking ownership
127
128     // TODO get rid of this with proper matrix type
129     matrix boxTemp;
130     for (int i = 0; i < DIM; i++)
131     {
132         for (int j = 0; j < DIM; j++)
133         {
134             boxTemp[i][j] = box[i * DIM + j];
135         }
136     }
137     const char* boxError = check_box(PbcType::Unset, boxTemp);
138     GMX_RELEASE_ASSERT(boxError == nullptr, boxError);
139
140     switch (mode)
141     {
142         case CodePath::CPU: invertBoxMatrix(boxTemp, pme->recipbox); break;
143
144         case CodePath::GPU:
145             pme_gpu_set_testing(pme->gpu, true);
146             pme_gpu_update_input_box(pme->gpu, boxTemp);
147             break;
148
149         default: GMX_THROW(InternalError("Test not implemented for this mode"));
150     }
151
152     return pme;
153 }
154
155 PmeSafePointer pmeInitEmpty(const t_inputrec* inputRec)
156 {
157     const Matrix3x3 defaultBox = { { 1.0F, 0.0F, 0.0F, 0.0F, 1.0F, 0.0F, 0.0F, 0.0F, 1.0F } };
158     return pmeInitWrapper(inputRec, CodePath::CPU, nullptr, nullptr, nullptr, defaultBox, 0.0F, 0.0F);
159 }
160
161 //! Make a GPU state-propagator manager
162 std::unique_ptr<StatePropagatorDataGpu> makeStatePropagatorDataGpu(const gmx_pme_t&     pme,
163                                                                    const DeviceContext* deviceContext,
164                                                                    const DeviceStream* deviceStream)
165 {
166     // TODO: Pin the host buffer and use async memory copies
167     // TODO: Special constructor for PME-only rank / PME-tests is used here. There should be a mechanism to
168     //       restrict one from using other constructor here.
169     return std::make_unique<StatePropagatorDataGpu>(deviceStream, *deviceContext, GpuApiCallBehavior::Sync,
170                                                     pme_gpu_get_block_size(&pme), nullptr);
171 }
172
173 //! PME initialization with atom data
174 void pmeInitAtoms(gmx_pme_t*               pme,
175                   StatePropagatorDataGpu*  stateGpu,
176                   const CodePath           mode,
177                   const CoordinatesVector& coordinates,
178                   const ChargesVector&     charges)
179 {
180     const index atomCount = coordinates.size();
181     GMX_RELEASE_ASSERT(atomCount == charges.ssize(), "Mismatch in atom data");
182     PmeAtomComm* atc = nullptr;
183
184     switch (mode)
185     {
186         case CodePath::CPU:
187             atc              = &(pme->atc[0]);
188             atc->x           = coordinates;
189             atc->coefficient = charges;
190             gmx_pme_reinit_atoms(pme, atomCount, charges.data(), nullptr);
191             /* With decomposition there would be more boilerplate atc code here, e.g. do_redist_pos_coeffs */
192             break;
193
194         case CodePath::GPU:
195             // TODO: Avoid use of atc in the GPU code path
196             atc = &(pme->atc[0]);
197             // We need to set atc->n for passing the size in the tests
198             atc->setNumAtoms(atomCount);
199             gmx_pme_reinit_atoms(pme, atomCount, charges.data(), nullptr);
200
201             stateGpu->reinit(atomCount, atomCount);
202             stateGpu->copyCoordinatesToGpu(arrayRefFromArray(coordinates.data(), coordinates.size()),
203                                            gmx::AtomLocality::All);
204             pme_gpu_set_kernelparam_coordinates(pme->gpu, stateGpu->getCoordinates());
205
206             break;
207
208         default: GMX_THROW(InternalError("Test not implemented for this mode"));
209     }
210 }
211
212 //! Getting local PME real grid pointer for test I/O
213 static real* pmeGetRealGridInternal(const gmx_pme_t* pme)
214 {
215     const size_t gridIndex = 0;
216     return pme->fftgrid[gridIndex];
217 }
218
219 //! Getting local PME real grid dimensions
220 static void pmeGetRealGridSizesInternal(const gmx_pme_t* pme,
221                                         CodePath         mode,
222                                         IVec& gridSize,       //NOLINT(google-runtime-references)
223                                         IVec& paddedGridSize) //NOLINT(google-runtime-references)
224 {
225     const size_t gridIndex = 0;
226     IVec         gridOffsetUnused;
227     switch (mode)
228     {
229         case CodePath::CPU:
230             gmx_parallel_3dfft_real_limits(pme->pfft_setup[gridIndex], gridSize, gridOffsetUnused,
231                                            paddedGridSize);
232             break;
233
234         case CodePath::GPU:
235             pme_gpu_get_real_grid_sizes(pme->gpu, &gridSize, &paddedGridSize);
236             break;
237
238         default: GMX_THROW(InternalError("Test not implemented for this mode"));
239     }
240 }
241
242 //! Getting local PME complex grid pointer for test I/O
243 static t_complex* pmeGetComplexGridInternal(const gmx_pme_t* pme)
244 {
245     const size_t gridIndex = 0;
246     return pme->cfftgrid[gridIndex];
247 }
248
249 //! Getting local PME complex grid dimensions
250 static void pmeGetComplexGridSizesInternal(const gmx_pme_t* pme,
251                                            IVec& gridSize,       //NOLINT(google-runtime-references)
252                                            IVec& paddedGridSize) //NOLINT(google-runtime-references)
253 {
254     const size_t gridIndex = 0;
255     IVec         gridOffsetUnused, complexOrderUnused;
256     gmx_parallel_3dfft_complex_limits(pme->pfft_setup[gridIndex], complexOrderUnused, gridSize,
257                                       gridOffsetUnused, paddedGridSize); // TODO: what about YZX ordering?
258 }
259
260 //! Getting the PME grid memory buffer and its sizes - template definition
261 template<typename ValueType>
262 static void pmeGetGridAndSizesInternal(const gmx_pme_t* /*unused*/,
263                                        CodePath /*unused*/,
264                                        ValueType*& /*unused*/, //NOLINT(google-runtime-references)
265                                        IVec& /*unused*/,       //NOLINT(google-runtime-references)
266                                        IVec& /*unused*/)       //NOLINT(google-runtime-references)
267 {
268     GMX_THROW(InternalError("Deleted function call"));
269     // explicitly deleting general template does not compile in clang/icc, see https://llvm.org/bugs/show_bug.cgi?id=17537
270 }
271
272 //! Getting the PME real grid memory buffer and its sizes
273 template<>
274 void pmeGetGridAndSizesInternal<real>(const gmx_pme_t* pme, CodePath mode, real*& grid, IVec& gridSize, IVec& paddedGridSize)
275 {
276     grid = pmeGetRealGridInternal(pme);
277     pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSize);
278 }
279
280 //! Getting the PME complex grid memory buffer and its sizes
281 template<>
282 void pmeGetGridAndSizesInternal<t_complex>(const gmx_pme_t* pme,
283                                            CodePath /*unused*/,
284                                            t_complex*& grid,
285                                            IVec&       gridSize,
286                                            IVec&       paddedGridSize)
287 {
288     grid = pmeGetComplexGridInternal(pme);
289     pmeGetComplexGridSizesInternal(pme, gridSize, paddedGridSize);
290 }
291
292 //! PME spline calculation and charge spreading
293 void pmePerformSplineAndSpread(gmx_pme_t* pme,
294                                CodePath   mode, // TODO const qualifiers elsewhere
295                                bool       computeSplines,
296                                bool       spreadCharges)
297 {
298     GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
299     PmeAtomComm* atc                          = &(pme->atc[0]);
300     const size_t gridIndex                    = 0;
301     const bool   computeSplinesForZeroCharges = true;
302     real**       fftgrid                      = spreadCharges ? pme->fftgrid : nullptr;
303     real*        pmegrid                      = pme->pmegrid[gridIndex].grid.grid;
304
305     switch (mode)
306     {
307         case CodePath::CPU:
308             spread_on_grid(pme, atc, &pme->pmegrid[gridIndex], computeSplines, spreadCharges,
309                            fftgrid != nullptr ? fftgrid[gridIndex] : nullptr,
310                            computeSplinesForZeroCharges, gridIndex);
311             if (spreadCharges && !pme->bUseThreads)
312             {
313                 wrap_periodic_pmegrid(pme, pmegrid);
314                 copy_pmegrid_to_fftgrid(
315                         pme, pmegrid, fftgrid != nullptr ? fftgrid[gridIndex] : nullptr, gridIndex);
316             }
317             break;
318
319 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
320  * double precision. GPUs are not used with double precision anyhow. */
321 #if !GMX_DOUBLE
322         case CodePath::GPU:
323         {
324             const real lambdaQ = 1.0;
325             // no synchronization needed as x is transferred in the PME stream
326             GpuEventSynchronizer* xReadyOnDevice = nullptr;
327             pme_gpu_spread(pme->gpu, xReadyOnDevice, fftgrid, computeSplines, spreadCharges, lambdaQ);
328         }
329         break;
330 #endif
331
332         default: GMX_THROW(InternalError("Test not implemented for this mode"));
333     }
334 }
335
336 //! Getting the internal spline data buffer pointer
337 static real* pmeGetSplineDataInternal(const gmx_pme_t* pme, PmeSplineDataType type, int dimIndex)
338 {
339     GMX_ASSERT((0 <= dimIndex) && (dimIndex < DIM), "Invalid dimension index");
340     const PmeAtomComm* atc          = &(pme->atc[0]);
341     const size_t       threadIndex  = 0;
342     real*              splineBuffer = nullptr;
343     switch (type)
344     {
345         case PmeSplineDataType::Values:
346             splineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
347             break;
348
349         case PmeSplineDataType::Derivatives:
350             splineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
351             break;
352
353         default: GMX_THROW(InternalError("Unknown spline data type"));
354     }
355     return splineBuffer;
356 }
357
358 //! PME solving
359 void pmePerformSolve(const gmx_pme_t*  pme,
360                      CodePath          mode,
361                      PmeSolveAlgorithm method,
362                      real              cellVolume,
363                      GridOrdering      gridOrdering,
364                      bool              computeEnergyAndVirial)
365 {
366     t_complex*   h_grid              = pmeGetComplexGridInternal(pme);
367     const bool   useLorentzBerthelot = false;
368     const size_t threadIndex         = 0;
369     const size_t gridIndex           = 0;
370     switch (mode)
371     {
372         case CodePath::CPU:
373             if (gridOrdering != GridOrdering::YZX)
374             {
375                 GMX_THROW(InternalError("Test not implemented for this mode"));
376             }
377             switch (method)
378             {
379                 case PmeSolveAlgorithm::Coulomb:
380                     solve_pme_yzx(pme, h_grid, cellVolume, computeEnergyAndVirial, pme->nthread, threadIndex);
381                     break;
382
383                 case PmeSolveAlgorithm::LennardJones:
384                     solve_pme_lj_yzx(pme, &h_grid, useLorentzBerthelot, cellVolume,
385                                      computeEnergyAndVirial, pme->nthread, threadIndex);
386                     break;
387
388                 default: GMX_THROW(InternalError("Test not implemented for this mode"));
389             }
390             break;
391
392         case CodePath::GPU:
393             switch (method)
394             {
395                 case PmeSolveAlgorithm::Coulomb:
396                     pme_gpu_solve(pme->gpu, gridIndex, h_grid, gridOrdering, computeEnergyAndVirial);
397                     break;
398
399                 default: GMX_THROW(InternalError("Test not implemented for this mode"));
400             }
401             break;
402
403         default: GMX_THROW(InternalError("Test not implemented for this mode"));
404     }
405 }
406
407 //! PME force gathering
408 void pmePerformGather(gmx_pme_t* pme, CodePath mode, ForcesVector& forces)
409 {
410     PmeAtomComm* atc       = &(pme->atc[0]);
411     const index  atomCount = atc->numAtoms();
412     GMX_RELEASE_ASSERT(forces.ssize() == atomCount, "Invalid force buffer size");
413     const real   scale       = 1.0;
414     const size_t threadIndex = 0;
415     const size_t gridIndex   = 0;
416     real*        pmegrid     = pme->pmegrid[gridIndex].grid.grid;
417     real**       fftgrid     = pme->fftgrid;
418
419     switch (mode)
420     {
421         case CodePath::CPU:
422             atc->f = forces;
423             if (atc->nthread == 1)
424             {
425                 // something which is normally done in serial spline computation (make_thread_local_ind())
426                 atc->spline[threadIndex].n = atomCount;
427             }
428             copy_fftgrid_to_pmegrid(pme, fftgrid[gridIndex], pmegrid, gridIndex, pme->nthread, threadIndex);
429             unwrap_periodic_pmegrid(pme, pmegrid);
430             gather_f_bsplines(pme, pmegrid, true, atc, &atc->spline[threadIndex], scale);
431             break;
432
433 /* The compiler will complain about passing fftgrid (converting double ** to float **) if using
434  * double precision. GPUs are not used with double precision anyhow. */
435 #if !GMX_DOUBLE
436         case CodePath::GPU:
437         {
438             // Variable initialization needs a non-switch scope
439             const bool computeEnergyAndVirial = false;
440             const real lambdaQ                = 1.0;
441             PmeOutput  output = pme_gpu_getOutput(*pme, computeEnergyAndVirial, lambdaQ);
442             GMX_ASSERT(forces.size() == output.forces_.size(),
443                        "Size of force buffers did not match");
444             pme_gpu_gather(pme->gpu, fftgrid, lambdaQ);
445             std::copy(std::begin(output.forces_), std::end(output.forces_), std::begin(forces));
446         }
447         break;
448 #endif
449
450         default: GMX_THROW(InternalError("Test not implemented for this mode"));
451     }
452 }
453
454 //! PME test finalization before fetching the outputs
455 void pmeFinalizeTest(const gmx_pme_t* pme, CodePath mode)
456 {
457     switch (mode)
458     {
459         case CodePath::CPU: break;
460
461         case CodePath::GPU: pme_gpu_synchronize(pme->gpu); break;
462
463         default: GMX_THROW(InternalError("Test not implemented for this mode"));
464     }
465 }
466
467 //! A binary enum for spline data layout transformation
468 enum class PmeLayoutTransform
469 {
470     GpuToHost,
471     HostToGpu
472 };
473
474 /*! \brief Gets a unique index to an element in a spline parameter buffer.
475  *
476  * These theta/dtheta buffers are laid out for GPU spread/gather
477  * kernels. The index is wrt the execution block, in range(0,
478  * atomsPerBlock * order * DIM).
479  *
480  * This is a wrapper, only used in unit tests.
481  * \param[in] order            PME order
482  * \param[in] splineIndex      Spline contribution index (from 0 to \p order - 1)
483  * \param[in] dimIndex         Dimension index (from 0 to 2)
484  * \param[in] atomIndex        Atom index wrt the block.
485  * \param[in] atomsPerWarp     Number of atoms processed by a warp.
486  *
487  * \returns Index into theta or dtheta array using GPU layout.
488  */
489 static int getSplineParamFullIndex(int order, int splineIndex, int dimIndex, int atomIndex, int atomsPerWarp)
490 {
491     if (order != c_pmeGpuOrder)
492     {
493         throw order;
494     }
495     constexpr int fixedOrder = c_pmeGpuOrder;
496     GMX_UNUSED_VALUE(fixedOrder);
497
498     const int atomWarpIndex = atomIndex % atomsPerWarp;
499     const int warpIndex     = atomIndex / atomsPerWarp;
500     int       indexBase, result;
501     switch (atomsPerWarp)
502     {
503         case 1:
504             indexBase = getSplineParamIndexBase<fixedOrder, 1>(warpIndex, atomWarpIndex);
505             result    = getSplineParamIndex<fixedOrder, 1>(indexBase, dimIndex, splineIndex);
506             break;
507
508         case 2:
509             indexBase = getSplineParamIndexBase<fixedOrder, 2>(warpIndex, atomWarpIndex);
510             result    = getSplineParamIndex<fixedOrder, 2>(indexBase, dimIndex, splineIndex);
511             break;
512
513         case 4:
514             indexBase = getSplineParamIndexBase<fixedOrder, 4>(warpIndex, atomWarpIndex);
515             result    = getSplineParamIndex<fixedOrder, 4>(indexBase, dimIndex, splineIndex);
516             break;
517
518         case 8:
519             indexBase = getSplineParamIndexBase<fixedOrder, 8>(warpIndex, atomWarpIndex);
520             result    = getSplineParamIndex<fixedOrder, 8>(indexBase, dimIndex, splineIndex);
521             break;
522
523         default:
524             GMX_THROW(NotImplementedError(
525                     formatString("Test function call not unrolled for atomsPerWarp = %d in "
526                                  "getSplineParamFullIndex",
527                                  atomsPerWarp)));
528     }
529     return result;
530 }
531
532 /*!\brief Return the number of atoms per warp */
533 static int pme_gpu_get_atoms_per_warp(const PmeGpu* pmeGpu)
534 {
535     const int order = pmeGpu->common->pme_order;
536     const int threadsPerAtom =
537             (pmeGpu->settings.threadsPerAtom == ThreadsPerAtom::Order ? order : order * order);
538     return pmeGpu->programHandle_->warpSize() / threadsPerAtom;
539 }
540
541 /*! \brief Rearranges the atom spline data between the GPU and host layouts.
542  * Only used for test purposes so far, likely to be horribly slow.
543  *
544  * \param[in]  pmeGpu     The PME GPU structure.
545  * \param[out] atc        The PME CPU atom data structure (with a single-threaded layout).
546  * \param[in]  type       The spline data type (values or derivatives).
547  * \param[in]  dimIndex   Dimension index.
548  * \param[in]  transform  Layout transform type
549  */
550 static void pme_gpu_transform_spline_atom_data(const PmeGpu*      pmeGpu,
551                                                const PmeAtomComm* atc,
552                                                PmeSplineDataType  type,
553                                                int                dimIndex,
554                                                PmeLayoutTransform transform)
555 {
556     // The GPU atom spline data is laid out in a different way currently than the CPU one.
557     // This function converts the data from GPU to CPU layout (in the host memory).
558     // It is only intended for testing purposes so far.
559     // Ideally we should use similar layouts on CPU and GPU if we care about mixed modes and their
560     // performance (e.g. spreading on GPU, gathering on CPU).
561     GMX_RELEASE_ASSERT(atc->nthread == 1, "Only the serial PME data layout is supported");
562     const uintmax_t threadIndex  = 0;
563     const auto      atomCount    = atc->numAtoms();
564     const auto      atomsPerWarp = pme_gpu_get_atoms_per_warp(pmeGpu);
565     const auto      pmeOrder     = pmeGpu->common->pme_order;
566     GMX_ASSERT(pmeOrder == c_pmeGpuOrder, "Only PME order 4 is implemented");
567
568     real*  cpuSplineBuffer;
569     float* h_splineBuffer;
570     switch (type)
571     {
572         case PmeSplineDataType::Values:
573             cpuSplineBuffer = atc->spline[threadIndex].theta.coefficients[dimIndex];
574             h_splineBuffer  = pmeGpu->staging.h_theta;
575             break;
576
577         case PmeSplineDataType::Derivatives:
578             cpuSplineBuffer = atc->spline[threadIndex].dtheta.coefficients[dimIndex];
579             h_splineBuffer  = pmeGpu->staging.h_dtheta;
580             break;
581
582         default: GMX_THROW(InternalError("Unknown spline data type"));
583     }
584
585     for (auto atomIndex = 0; atomIndex < atomCount; atomIndex++)
586     {
587         for (auto orderIndex = 0; orderIndex < pmeOrder; orderIndex++)
588         {
589             const auto gpuValueIndex =
590                     getSplineParamFullIndex(pmeOrder, orderIndex, dimIndex, atomIndex, atomsPerWarp);
591             const auto cpuValueIndex = atomIndex * pmeOrder + orderIndex;
592             GMX_ASSERT(cpuValueIndex < atomCount * pmeOrder,
593                        "Atom spline data index out of bounds (while transforming GPU data layout "
594                        "for host)");
595             switch (transform)
596             {
597                 case PmeLayoutTransform::GpuToHost:
598                     cpuSplineBuffer[cpuValueIndex] = h_splineBuffer[gpuValueIndex];
599                     break;
600
601                 case PmeLayoutTransform::HostToGpu:
602                     h_splineBuffer[gpuValueIndex] = cpuSplineBuffer[cpuValueIndex];
603                     break;
604
605                 default: GMX_THROW(InternalError("Unknown layout transform"));
606             }
607         }
608     }
609 }
610
611 //! Setting atom spline values/derivatives to be used in spread/gather
612 void pmeSetSplineData(const gmx_pme_t*             pme,
613                       CodePath                     mode,
614                       const SplineParamsDimVector& splineValues,
615                       PmeSplineDataType            type,
616                       int                          dimIndex)
617 {
618     const PmeAtomComm* atc       = &(pme->atc[0]);
619     const index        atomCount = atc->numAtoms();
620     const index        pmeOrder  = pme->pme_order;
621     const index        dimSize   = pmeOrder * atomCount;
622     GMX_RELEASE_ASSERT(dimSize == splineValues.ssize(), "Mismatch in spline data");
623     real* splineBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
624
625     switch (mode)
626     {
627         case CodePath::CPU:
628             std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
629             break;
630
631         case CodePath::GPU:
632             std::copy(splineValues.begin(), splineValues.end(), splineBuffer);
633             pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::HostToGpu);
634             break;
635
636         default: GMX_THROW(InternalError("Test not implemented for this mode"));
637     }
638 }
639
640 //! Setting gridline indices to be used in spread/gather
641 void pmeSetGridLineIndices(gmx_pme_t* pme, CodePath mode, const GridLineIndicesVector& gridLineIndices)
642 {
643     PmeAtomComm* atc       = &(pme->atc[0]);
644     const index  atomCount = atc->numAtoms();
645     GMX_RELEASE_ASSERT(atomCount == gridLineIndices.ssize(), "Mismatch in gridline indices size");
646
647     IVec paddedGridSizeUnused, gridSize(0, 0, 0);
648     pmeGetRealGridSizesInternal(pme, mode, gridSize, paddedGridSizeUnused);
649
650     for (const auto& index : gridLineIndices)
651     {
652         for (int i = 0; i < DIM; i++)
653         {
654             GMX_RELEASE_ASSERT((0 <= index[i]) && (index[i] < gridSize[i]),
655                                "Invalid gridline index");
656         }
657     }
658
659     switch (mode)
660     {
661         case CodePath::GPU:
662             memcpy(pme_gpu_staging(pme->gpu).h_gridlineIndices, gridLineIndices.data(),
663                    atomCount * sizeof(gridLineIndices[0]));
664             break;
665
666         case CodePath::CPU:
667             atc->idx.resize(gridLineIndices.size());
668             std::copy(gridLineIndices.begin(), gridLineIndices.end(), atc->idx.begin());
669             break;
670         default: GMX_THROW(InternalError("Test not implemented for this mode"));
671     }
672 }
673
674 //! Getting plain index into the complex 3d grid
675 inline size_t pmeGetGridPlainIndexInternal(const IVec& index, const IVec& paddedGridSize, GridOrdering gridOrdering)
676 {
677     size_t result;
678     switch (gridOrdering)
679     {
680         case GridOrdering::YZX:
681             result = (index[YY] * paddedGridSize[ZZ] + index[ZZ]) * paddedGridSize[XX] + index[XX];
682             break;
683
684         case GridOrdering::XYZ:
685             result = (index[XX] * paddedGridSize[YY] + index[YY]) * paddedGridSize[ZZ] + index[ZZ];
686             break;
687
688         default: GMX_THROW(InternalError("Test not implemented for this mode"));
689     }
690     return result;
691 }
692
693 //! Setting real or complex grid
694 template<typename ValueType>
695 static void pmeSetGridInternal(const gmx_pme_t*                        pme,
696                                CodePath                                mode,
697                                GridOrdering                            gridOrdering,
698                                const SparseGridValuesInput<ValueType>& gridValues)
699 {
700     IVec       gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
701     ValueType* grid;
702     pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
703
704     switch (mode)
705     {
706         case CodePath::GPU: // intentional absence of break, the grid will be copied from the host buffer in testing mode
707         case CodePath::CPU:
708             std::memset(grid, 0,
709                         paddedGridSize[XX] * paddedGridSize[YY] * paddedGridSize[ZZ] * sizeof(ValueType));
710             for (const auto& gridValue : gridValues)
711             {
712                 for (int i = 0; i < DIM; i++)
713                 {
714                     GMX_RELEASE_ASSERT((0 <= gridValue.first[i]) && (gridValue.first[i] < gridSize[i]),
715                                        "Invalid grid value index");
716                 }
717                 const size_t gridValueIndex =
718                         pmeGetGridPlainIndexInternal(gridValue.first, paddedGridSize, gridOrdering);
719                 grid[gridValueIndex] = gridValue.second;
720             }
721             break;
722
723         default: GMX_THROW(InternalError("Test not implemented for this mode"));
724     }
725 }
726
727 //! Setting real grid to be used in gather
728 void pmeSetRealGrid(const gmx_pme_t* pme, CodePath mode, const SparseRealGridValuesInput& gridValues)
729 {
730     pmeSetGridInternal<real>(pme, mode, GridOrdering::XYZ, gridValues);
731 }
732
733 //! Setting complex grid to be used in solve
734 void pmeSetComplexGrid(const gmx_pme_t*                    pme,
735                        CodePath                            mode,
736                        GridOrdering                        gridOrdering,
737                        const SparseComplexGridValuesInput& gridValues)
738 {
739     pmeSetGridInternal<t_complex>(pme, mode, gridOrdering, gridValues);
740 }
741
742 //! Getting the single dimension's spline values or derivatives
743 SplineParamsDimVector pmeGetSplineData(const gmx_pme_t* pme, CodePath mode, PmeSplineDataType type, int dimIndex)
744 {
745     GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
746     const PmeAtomComm* atc       = &(pme->atc[0]);
747     const size_t       atomCount = atc->numAtoms();
748     const size_t       pmeOrder  = pme->pme_order;
749     const size_t       dimSize   = pmeOrder * atomCount;
750
751     real*                 sourceBuffer = pmeGetSplineDataInternal(pme, type, dimIndex);
752     SplineParamsDimVector result;
753     switch (mode)
754     {
755         case CodePath::GPU:
756             pme_gpu_transform_spline_atom_data(pme->gpu, atc, type, dimIndex, PmeLayoutTransform::GpuToHost);
757             result = arrayRefFromArray(sourceBuffer, dimSize);
758             break;
759
760         case CodePath::CPU: result = arrayRefFromArray(sourceBuffer, dimSize); break;
761
762         default: GMX_THROW(InternalError("Test not implemented for this mode"));
763     }
764     return result;
765 }
766
767 //! Getting the gridline indices
768 GridLineIndicesVector pmeGetGridlineIndices(const gmx_pme_t* pme, CodePath mode)
769 {
770     GMX_RELEASE_ASSERT(pme != nullptr, "PME data is not initialized");
771     const PmeAtomComm* atc       = &(pme->atc[0]);
772     const size_t       atomCount = atc->numAtoms();
773
774     GridLineIndicesVector gridLineIndices;
775     switch (mode)
776     {
777         case CodePath::GPU:
778             gridLineIndices = arrayRefFromArray(
779                     reinterpret_cast<IVec*>(pme_gpu_staging(pme->gpu).h_gridlineIndices), atomCount);
780             break;
781
782         case CodePath::CPU: gridLineIndices = atc->idx; break;
783
784         default: GMX_THROW(InternalError("Test not implemented for this mode"));
785     }
786     return gridLineIndices;
787 }
788
789 //! Getting real or complex grid - only non zero values
790 template<typename ValueType>
791 static SparseGridValuesOutput<ValueType> pmeGetGridInternal(const gmx_pme_t* pme,
792                                                             CodePath         mode,
793                                                             GridOrdering     gridOrdering)
794 {
795     IVec       gridSize(0, 0, 0), paddedGridSize(0, 0, 0);
796     ValueType* grid;
797     pmeGetGridAndSizesInternal<ValueType>(pme, mode, grid, gridSize, paddedGridSize);
798     SparseGridValuesOutput<ValueType> gridValues;
799     switch (mode)
800     {
801         case CodePath::GPU: // intentional absence of break
802         case CodePath::CPU:
803             gridValues.clear();
804             for (int ix = 0; ix < gridSize[XX]; ix++)
805             {
806                 for (int iy = 0; iy < gridSize[YY]; iy++)
807                 {
808                     for (int iz = 0; iz < gridSize[ZZ]; iz++)
809                     {
810                         IVec         temp(ix, iy, iz);
811                         const size_t gridValueIndex =
812                                 pmeGetGridPlainIndexInternal(temp, paddedGridSize, gridOrdering);
813                         const ValueType value = grid[gridValueIndex];
814                         if (value != ValueType{})
815                         {
816                             auto key        = formatString("Cell %d %d %d", ix, iy, iz);
817                             gridValues[key] = value;
818                         }
819                     }
820                 }
821             }
822             break;
823
824         default: GMX_THROW(InternalError("Test not implemented for this mode"));
825     }
826     return gridValues;
827 }
828
829 //! Getting the real grid (spreading output of pmePerformSplineAndSpread())
830 SparseRealGridValuesOutput pmeGetRealGrid(const gmx_pme_t* pme, CodePath mode)
831 {
832     return pmeGetGridInternal<real>(pme, mode, GridOrdering::XYZ);
833 }
834
835 //! Getting the complex grid output of pmePerformSolve()
836 SparseComplexGridValuesOutput pmeGetComplexGrid(const gmx_pme_t* pme, CodePath mode, GridOrdering gridOrdering)
837 {
838     return pmeGetGridInternal<t_complex>(pme, mode, gridOrdering);
839 }
840
841 //! Getting the reciprocal energy and virial
842 PmeOutput pmeGetReciprocalEnergyAndVirial(const gmx_pme_t* pme, CodePath mode, PmeSolveAlgorithm method)
843 {
844     PmeOutput  output;
845     const real lambdaQ = 1.0;
846     switch (mode)
847     {
848         case CodePath::CPU:
849             switch (method)
850             {
851                 case PmeSolveAlgorithm::Coulomb:
852                     get_pme_ener_vir_q(pme->solve_work, pme->nthread, &output);
853                     break;
854
855                 case PmeSolveAlgorithm::LennardJones:
856                     get_pme_ener_vir_lj(pme->solve_work, pme->nthread, &output);
857                     break;
858
859                 default: GMX_THROW(InternalError("Test not implemented for this mode"));
860             }
861             break;
862         case CodePath::GPU:
863             switch (method)
864             {
865                 case PmeSolveAlgorithm::Coulomb:
866                     pme_gpu_getEnergyAndVirial(*pme, lambdaQ, &output);
867                     break;
868
869                 default: GMX_THROW(InternalError("Test not implemented for this mode"));
870             }
871             break;
872
873         default: GMX_THROW(InternalError("Test not implemented for this mode"));
874     }
875     return output;
876 }
877
878 const char* codePathToString(CodePath codePath)
879 {
880     switch (codePath)
881     {
882         case CodePath::CPU: return "CPU";
883         case CodePath::GPU: return "GPU";
884         default: GMX_THROW(NotImplementedError("This CodePath should support codePathToString"));
885     }
886 }
887
888 PmeTestHardwareContext::PmeTestHardwareContext() : codePath_(CodePath::CPU) {}
889
890 PmeTestHardwareContext::PmeTestHardwareContext(TestDevice* testDevice) :
891     codePath_(CodePath::CPU),
892     testDevice_(testDevice)
893 {
894     setActiveDevice(testDevice_->deviceInfo());
895     pmeGpuProgram_ = buildPmeGpuProgram(testDevice_->deviceContext());
896 }
897
898 //! Returns a human-readable context description line
899 std::string PmeTestHardwareContext::description() const
900 {
901     switch (codePath_)
902     {
903         case CodePath::CPU: return "CPU";
904         case CodePath::GPU: return "GPU (" + testDevice_->description() + ")";
905         default: return "Unknown code path.";
906     }
907 }
908
909 void PmeTestHardwareContext::activate() const
910 {
911     if (codePath_ == CodePath::GPU)
912     {
913         setActiveDevice(testDevice_->deviceInfo());
914     }
915 }
916
917 std::vector<std::unique_ptr<PmeTestHardwareContext>> createPmeTestHardwareContextList()
918 {
919     std::vector<std::unique_ptr<PmeTestHardwareContext>> pmeTestHardwareContextList;
920     // Add CPU
921     pmeTestHardwareContextList.emplace_back(std::make_unique<PmeTestHardwareContext>());
922     // Add GPU devices
923     const auto& testDeviceList = getTestHardwareEnvironment()->getTestDeviceList();
924     for (const auto& testDevice : testDeviceList)
925     {
926         pmeTestHardwareContextList.emplace_back(std::make_unique<PmeTestHardwareContext>(testDevice.get()));
927     }
928     return pmeTestHardwareContextList;
929 }
930
931 } // namespace test
932 } // namespace gmx