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