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