Refactor PME tests for better usability
[alexxy/gromacs.git] / src / gromacs / ewald / tests / pmegathertest.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020,2021, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
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 PME force gathering tests.
38  *
39  * \author Aleksei Iupinov <a.yupinov@gmail.com>
40  * \author Mark Abraham <mark.j.abraham@gmail.com>
41  * \ingroup module_ewald
42  */
43
44 #include "gmxpre.h"
45
46 #include <string>
47
48 #include <gmock/gmock.h>
49
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/utility/stringutil.h"
52
53 #include "testutils/refdata.h"
54 #include "testutils/test_hardware_environment.h"
55 #include "testutils/testasserts.h"
56
57 #include "pmetestcommon.h"
58
59 namespace gmx
60 {
61 namespace test
62 {
63 namespace
64 {
65
66 //! A couple of valid inputs for grid sizes
67 std::vector<IVec> const c_inputGridSizes{ IVec{ 16, 12, 14 }, IVec{ 13, 15, 11 } };
68
69 //! A structure for all the spline data which depends in size both on the PME order and atom count
70 struct SplineData
71 {
72     //! Spline values
73     SplineParamsVector splineValues;
74     //! Spline derivatives
75     SplineParamsVector splineDerivatives;
76 };
77
78 //! Return synthetic spline data to gather
79 SplineData getSplineData(const int pmeOrder, const int atomCount)
80 {
81     // Spline values/derivatives below are also generated randomly, so
82     // they are bogus, but that should not affect the reproducibility,
83     // which is what we're after.
84
85     //! A lot of random input spline values - should have at least (max PME order = 5) * (DIM = 3) * (total unique atom number in all test cases = 13) values
86     static const std::vector<real> s_sampleSplineValuesFull{
87         0.12F, 0.81F, 0.29F, 0.22F, 0.13F, 0.19F, 0.12F, 0.8F,  0.44F, 0.38F, 0.32F, 0.36F, 0.27F,
88         0.11F, 0.17F, 0.94F, 0.07F, 0.9F,  0.98F, 0.96F, 0.07F, 0.94F, 0.77F, 0.24F, 0.84F, 0.16F,
89         0.77F, 0.57F, 0.52F, 0.27F, 0.39F, 0.45F, 0.6F,  0.59F, 0.44F, 0.91F, 0.97F, 0.43F, 0.24F,
90         0.52F, 0.73F, 0.55F, 0.99F, 0.39F, 0.97F, 0.35F, 0.1F,  0.68F, 0.19F, 0.1F,  0.77F, 0.2F,
91         0.43F, 0.69F, 0.76F, 0.32F, 0.31F, 0.94F, 0.53F, 0.6F,  0.93F, 0.57F, 0.94F, 0.88F, 0.75F,
92         0.77F, 0.91F, 0.72F, 0.07F, 0.78F, 0.09F, 0.02F, 0.48F, 0.97F, 0.89F, 0.39F, 0.48F, 0.19F,
93         0.02F, 0.92F, 0.8F,  0.41F, 0.53F, 0.32F, 0.38F, 0.58F, 0.36F, 0.46F, 0.92F, 0.91F, 0.01F,
94         0.86F, 0.54F, 0.86F, 0.94F, 0.37F, 0.35F, 0.81F, 0.89F, 0.48F, 0.34F, 0.18F, 0.11F, 0.02F,
95         0.87F, 0.95F, 0.66F, 0.67F, 0.38F, 0.45F, 0.04F, 0.94F, 0.54F, 0.76F, 0.58F, 0.83F, 0.31F,
96         0.73F, 0.71F, 0.06F, 0.35F, 0.32F, 0.35F, 0.61F, 0.27F, 0.98F, 0.83F, 0.11F, 0.3F,  0.42F,
97         0.95F, 0.69F, 0.58F, 0.29F, 0.1F,  0.68F, 0.94F, 0.62F, 0.51F, 0.47F, 0.04F, 0.47F, 0.34F,
98         0.71F, 0.52F, 0.19F, 0.69F, 0.5F,  0.59F, 0.05F, 0.74F, 0.11F, 0.4F,  0.81F, 0.24F, 0.53F,
99         0.71F, 0.07F, 0.17F, 0.41F, 0.23F, 0.78F, 0.27F, 0.1F,  0.71F, 0.36F, 0.67F, 0.6F,  0.94F,
100         0.69F, 0.19F, 0.58F, 0.68F, 0.5F,  0.62F, 0.38F, 0.29F, 0.44F, 0.04F, 0.89F, 0.0F,  0.76F,
101         0.22F, 0.16F, 0.08F, 0.62F, 0.51F, 0.62F, 0.83F, 0.72F, 0.96F, 0.99F, 0.4F,  0.79F, 0.83F,
102         0.21F, 0.43F, 0.32F, 0.44F, 0.72F, 0.21F, 0.4F,  0.93F, 0.07F, 0.11F, 0.41F, 0.24F, 0.04F,
103         0.36F, 0.15F, 0.92F, 0.08F, 0.99F, 0.35F, 0.42F, 0.7F,  0.17F, 0.39F, 0.69F, 0.0F,  0.86F,
104         0.89F, 0.59F, 0.81F, 0.77F, 0.15F, 0.89F, 0.17F, 0.76F, 0.67F, 0.58F, 0.78F, 0.26F, 0.19F,
105         0.69F, 0.18F, 0.46F, 0.6F,  0.69F, 0.23F, 0.34F, 0.3F,  0.64F, 0.34F, 0.6F,  0.99F, 0.69F,
106         0.57F, 0.75F, 0.07F, 0.36F, 0.75F, 0.81F, 0.8F,  0.42F, 0.09F, 0.94F, 0.66F, 0.35F, 0.67F,
107         0.34F, 0.66F, 0.02F, 0.47F, 0.78F, 0.21F, 0.02F, 0.18F, 0.42F, 0.2F,  0.46F, 0.34F, 0.4F,
108         0.46F, 0.96F, 0.86F, 0.25F, 0.25F, 0.22F, 0.37F, 0.59F, 0.19F, 0.45F, 0.61F, 0.04F, 0.71F,
109         0.77F, 0.51F, 0.77F, 0.15F, 0.78F, 0.36F, 0.62F, 0.24F, 0.86F, 0.2F,  0.77F, 0.08F, 0.09F,
110         0.3F,  0.0F,  0.6F,  0.99F, 0.69F,
111     };
112
113     //! A lot of random input spline derivatives - should have at least (max PME order = 5) * (DIM = 3) * (total unique atom number in all test cases = 13) values
114     static const std::vector<real> s_sampleSplineDerivativesFull{
115         0.82F, 0.88F, 0.83F, 0.11F, 0.93F, 0.32F, 0.71F, 0.37F, 0.69F, 0.88F, 0.11F, 0.38F, 0.25F,
116         0.5F,  0.36F, 0.81F, 0.78F, 0.31F, 0.66F, 0.32F, 0.27F, 0.35F, 0.53F, 0.83F, 0.08F, 0.08F,
117         0.94F, 0.71F, 0.65F, 0.24F, 0.13F, 0.01F, 0.33F, 0.65F, 0.24F, 0.53F, 0.45F, 0.84F, 0.33F,
118         0.97F, 0.31F, 0.7F,  0.03F, 0.31F, 0.41F, 0.76F, 0.12F, 0.3F,  0.57F, 0.65F, 0.87F, 0.99F,
119         0.42F, 0.97F, 0.32F, 0.39F, 0.73F, 0.23F, 0.03F, 0.67F, 0.97F, 0.57F, 0.42F, 0.38F, 0.54F,
120         0.17F, 0.53F, 0.54F, 0.18F, 0.8F,  0.76F, 0.13F, 0.29F, 0.83F, 0.77F, 0.56F, 0.4F,  0.87F,
121         0.36F, 0.18F, 0.59F, 0.04F, 0.05F, 0.61F, 0.26F, 0.91F, 0.62F, 0.16F, 0.89F, 0.23F, 0.26F,
122         0.59F, 0.33F, 0.2F,  0.49F, 0.41F, 0.25F, 0.4F,  0.16F, 0.83F, 0.44F, 0.82F, 0.21F, 0.95F,
123         0.14F, 0.8F,  0.37F, 0.31F, 0.41F, 0.53F, 0.15F, 0.85F, 0.78F, 0.17F, 0.92F, 0.03F, 0.13F,
124         0.2F,  0.03F, 0.33F, 0.87F, 0.38F, 0,     0.08F, 0.79F, 0.36F, 0.53F, 0.05F, 0.07F, 0.94F,
125         0.23F, 0.85F, 0.13F, 0.27F, 0.23F, 0.22F, 0.26F, 0.38F, 0.15F, 0.48F, 0.18F, 0.33F, 0.23F,
126         0.62F, 0.1F,  0.36F, 0.99F, 0.07F, 0.02F, 0.04F, 0.09F, 0.29F, 0.52F, 0.29F, 0.83F, 0.97F,
127         0.61F, 0.81F, 0.49F, 0.56F, 0.08F, 0.09F, 0.03F, 0.65F, 0.46F, 0.1F,  0.06F, 0.06F, 0.39F,
128         0.29F, 0.04F, 0.03F, 0.1F,  0.83F, 0.94F, 0.59F, 0.97F, 0.82F, 0.2F,  0.66F, 0.23F, 0.11F,
129         0.03F, 0.16F, 0.27F, 0.53F, 0.94F, 0.46F, 0.43F, 0.29F, 0.97F, 0.64F, 0.46F, 0.37F, 0.43F,
130         0.48F, 0.37F, 0.93F, 0.5F,  0.2F,  0.92F, 0.09F, 0.74F, 0.55F, 0.44F, 0.05F, 0.13F, 0.17F,
131         0.79F, 0.44F, 0.11F, 0.6F,  0.64F, 0.05F, 0.96F, 0.3F,  0.45F, 0.47F, 0.42F, 0.74F, 0.91F,
132         0.06F, 0.89F, 0.24F, 0.26F, 0.68F, 0.4F,  0.88F, 0.5F,  0.65F, 0.48F, 0.15F, 0.0F,  0.41F,
133         0.67F, 0.4F,  0.31F, 0.73F, 0.77F, 0.36F, 0.26F, 0.74F, 0.46F, 0.56F, 0.78F, 0.92F, 0.32F,
134         0.9F,  0.06F, 0.55F, 0.6F,  0.13F, 0.38F, 0.93F, 0.5F,  0.92F, 0.96F, 0.82F, 0.0F,  0.04F,
135         0.9F,  0.55F, 0.97F, 1.0F,  0.23F, 0.46F, 0.52F, 0.49F, 0.0F,  0.32F, 0.16F, 0.4F,  0.62F,
136         0.36F, 0.03F, 0.63F, 0.16F, 0.58F, 0.97F, 0.03F, 0.44F, 0.07F, 0.22F, 0.75F, 0.32F, 0.61F,
137         0.94F, 0.33F, 0.7F,  0.57F, 0.5F,  0.84F, 0.7F,  0.47F, 0.18F, 0.09F, 0.25F, 0.77F, 0.94F,
138         0.85F, 0.09F, 0.83F, 0.02F, 0.91F,
139     };
140
141     SplineData splineData;
142     const int  dimSize = atomCount * pmeOrder;
143     for (int dimIndex = 0; dimIndex < DIM; dimIndex++)
144     {
145         splineData.splineValues[dimIndex] =
146                 SplineParamsDimVector(s_sampleSplineValuesFull).subArray(dimIndex * dimSize, dimSize);
147         splineData.splineDerivatives[dimIndex] =
148                 SplineParamsDimVector(s_sampleSplineDerivativesFull).subArray(dimIndex * dimSize, dimSize);
149     }
150     return splineData;
151 }
152
153 //! Two input grids - only non-zero values have to be listed
154 const std::map<std::string, SparseRealGridValuesInput> c_inputGrids = {
155     { "first",
156       SparseRealGridValuesInput{
157               { IVec{ 0, 0, 0 }, 3.5F },
158               { IVec{ 7, 0, 0 }, -2.5F },
159               { IVec{ 3, 5, 7 }, -0.006F },
160               { IVec{ 1, 6, 7 }, -2.76F },
161               { IVec{ 3, 1, 2 }, 0.6F },
162               { IVec{ 6, 2, 4 }, 7.1F },
163               { IVec{ 4, 5, 6 }, 4.1F },
164               { IVec{ 4, 4, 6 }, -3.7F },
165       } },
166     { "second",
167       SparseRealGridValuesInput{
168               { IVec{ 0, 4, 0 }, 6.F },
169               { IVec{ 4, 2, 7 }, 13.76F },
170               { IVec{ 0, 6, 7 }, 3.6F },
171               { IVec{ 3, 2, 8 }, 0.61F },
172               { IVec{ 5, 4, 3 }, 2.1F },
173               { IVec{ 2, 5, 10 }, 3.6F },
174               { IVec{ 5, 3, 6 }, 2.1F },
175               { IVec{ 6, 6, 6 }, 6.1F },
176       } }
177 };
178
179 //! Input forces for reduction
180 std::vector<RVec> const c_sampleForcesFull{
181     RVec{ 0.02F, 0.87F, 0.95F }, RVec{ 0.66F, 0.67F, 0.38F }, RVec{ 0.45F, 0.04F, 0.94F },
182     RVec{ 0.54F, 0.76F, 0.58F }, RVec{ 0.83F, 0.31F, 0.73F }, RVec{ 0.71F, 0.06F, 0.35F },
183     RVec{ 0.32F, 0.35F, 0.61F }, RVec{ 0.27F, 0.98F, 0.83F }, RVec{ 0.11F, 0.3F, 0.42F },
184     RVec{ 0.95F, 0.69F, 0.58F }, RVec{ 0.29F, 0.1F, 0.68F },  RVec{ 0.94F, 0.62F, 0.51F },
185     RVec{ 0.47F, 0.04F, 0.47F }, RVec{ 0.34F, 0.71F, 0.52F }
186 };
187
188 /*! \brief A structure for the input atom data, which depends in size on atom count */
189 struct TestSystem
190 {
191     //! Gridline indices
192     GridLineIndicesVector gridLineIndices;
193     //! Charges
194     ChargesVector charges;
195     //! Coordinates
196     CoordinatesVector coordinates;
197 };
198
199 /*! \brief Test systems to use
200  *
201  * The coordinates are intentionally bogus in this test - only the
202  * size matters; the gridline indices are fed directly as inputs */
203 // TODO use NaN?
204 std::map<std::string, TestSystem> c_testSystems = {
205     { "1 atom",
206       { GridLineIndicesVector{ { IVec(4, 2, 6) } },
207         ChargesVector{ 4.95F },
208         CoordinatesVector(1, { 1e6, 1e7, -1e8 }) } },
209     { "2 atoms",
210       { GridLineIndicesVector{ { IVec(1, 4, 10), IVec(0, 6, 6) } },
211         ChargesVector{ { 3.11F, 3.97F } },
212         CoordinatesVector(2, { 1e6, 1e7, -1e8 }) } },
213     { "13 atoms",
214       { GridLineIndicesVector{ {
215                 IVec{ 0, 1, 4 },
216                 IVec{ 6, 3, 0 },
217                 IVec{ 7, 2, 2 },
218                 IVec{ 8, 3, 1 },
219                 IVec{ 4, 0, 3 },
220                 IVec{ 0, 0, 0 },
221                 IVec{ 8, 5, 8 },
222                 IVec{ 4, 4, 2 },
223                 IVec{ 7, 1, 7 },
224                 IVec{ 8, 5, 5 },
225                 IVec{ 2, 6, 5 },
226                 IVec{ 1, 6, 2 },
227                 IVec{ 7, 1, 8 },
228         } },
229         ChargesVector{ { 1.08F, 2.09F, 1.1F, 4.13F, 3.31F, 2.8F, 5.83F, 5.09F, 6.1F, 2.86F, 0.24F, 5.76F, 5.19F } },
230         CoordinatesVector(13, { 1e6, 1e7, -1e8 }) } },
231 };
232
233 /*! \brief Convenience typedef of the test input parameters - unit cell box, PME interpolation
234  * order, grid dimensions, grid values, overwriting/reducing the input forces, atom count.
235  */
236 typedef std::tuple<std::string, int, IVec, std::string, std::string> GatherInputParameters;
237
238 //! Help GoogleTest name our test cases
239 std::string nameOfTest(const testing::TestParamInfo<GatherInputParameters>& info)
240 {
241     std::string testName = formatString(
242             "box_%s_"
243             "order_%d_"
244             "grid_%d_%d_%d_"
245             "gridvalues_%s_"
246             "system_%s",
247             std::get<0>(info.param).c_str(),
248             std::get<1>(info.param),
249             std::get<2>(info.param)[XX],
250             std::get<2>(info.param)[YY],
251             std::get<2>(info.param)[ZZ],
252             std::get<3>(info.param).c_str(),
253             std::get<4>(info.param).c_str());
254
255     // Note that the returned names must be unique and may use only
256     // alphanumeric ASCII characters. It's not supposed to contain
257     // underscores (see the GoogleTest FAQ
258     // why-should-test-suite-names-and-test-names-not-contain-underscore),
259     // but doing so works for now, is likely to remain so, and makes
260     // such test names much more readable.
261     testName = replaceAll(testName, "-", "_");
262     testName = replaceAll(testName, ".", "_");
263     testName = replaceAll(testName, " ", "_");
264     return testName;
265 }
266
267 //! Test fixture
268 class GatherTest : public ::testing::TestWithParam<GatherInputParameters>
269 {
270 public:
271     GatherTest() = default;
272     //! Sets the input atom data references and programs once
273     static void SetUpTestSuite()
274     {
275         s_pmeTestHardwareContexts    = createPmeTestHardwareContextList();
276         g_allowPmeWithSyclForTesting = true; // We support Gather with SYCL
277     }
278
279     static void TearDownTestSuite()
280     {
281         // Revert the value back.
282         g_allowPmeWithSyclForTesting = false;
283     }
284
285     //! The test
286     static void runTest()
287     {
288         /* Getting the input */
289         int         pmeOrder;
290         IVec        gridSize;
291         std::string boxName, gridValuesName, testSystemName;
292         std::tie(boxName, pmeOrder, gridSize, gridValuesName, testSystemName) = GetParam();
293         Matrix3x3                        box               = c_inputBoxes.at(boxName);
294         const SparseRealGridValuesInput& nonZeroGridValues = c_inputGrids.at(gridValuesName);
295         TestSystem                       testSystem        = c_testSystems.at(testSystemName);
296         int                              atomCount         = testSystem.coordinates.size();
297         SplineData                       splineData        = getSplineData(pmeOrder, atomCount);
298
299         /* Storing the input where it's needed, running the test */
300         t_inputrec inputRec;
301         inputRec.nkx         = gridSize[XX];
302         inputRec.nky         = gridSize[YY];
303         inputRec.nkz         = gridSize[ZZ];
304         inputRec.pme_order   = pmeOrder;
305         inputRec.coulombtype = CoulombInteractionType::Pme;
306         inputRec.epsilon_r   = 1.0;
307
308         TestReferenceData refData;
309         for (const auto& pmeTestHardwareContext : s_pmeTestHardwareContexts)
310         {
311             pmeTestHardwareContext->activate();
312             CodePath   codePath       = pmeTestHardwareContext->codePath();
313             const bool supportedInput = pmeSupportsInputForMode(
314                     *getTestHardwareEnvironment()->hwinfo(), &inputRec, codePath);
315             if (!supportedInput)
316             {
317                 /* Testing the failure for the unsupported input */
318                 EXPECT_THROW_GMX(pmeInitWrapper(&inputRec, codePath, nullptr, nullptr, nullptr, box),
319                                  NotImplementedError);
320                 continue;
321             }
322
323             /* Describing the test uniquely */
324             SCOPED_TRACE(
325                     formatString("Testing force gathering on %s for PME grid size %d %d %d"
326                                  ", order %d, %d atoms",
327                                  pmeTestHardwareContext->description().c_str(),
328                                  gridSize[XX],
329                                  gridSize[YY],
330                                  gridSize[ZZ],
331                                  pmeOrder,
332                                  atomCount));
333
334             PmeSafePointer                          pmeSafe = pmeInitWrapper(&inputRec,
335                                                     codePath,
336                                                     pmeTestHardwareContext->deviceContext(),
337                                                     pmeTestHardwareContext->deviceStream(),
338                                                     pmeTestHardwareContext->pmeGpuProgram(),
339                                                     box);
340             std::unique_ptr<StatePropagatorDataGpu> stateGpu =
341                     (codePath == CodePath::GPU)
342                             ? makeStatePropagatorDataGpu(*pmeSafe.get(),
343                                                          pmeTestHardwareContext->deviceContext(),
344                                                          pmeTestHardwareContext->deviceStream())
345                             : nullptr;
346
347             pmeInitAtoms(pmeSafe.get(), stateGpu.get(), codePath, testSystem.coordinates, testSystem.charges);
348
349             /* Setting some more inputs */
350             pmeSetRealGrid(pmeSafe.get(), codePath, nonZeroGridValues);
351             pmeSetGridLineIndices(pmeSafe.get(), codePath, testSystem.gridLineIndices);
352             for (int dimIndex = 0; dimIndex < DIM; dimIndex++)
353             {
354                 pmeSetSplineData(pmeSafe.get(),
355                                  codePath,
356                                  splineData.splineValues[dimIndex],
357                                  PmeSplineDataType::Values,
358                                  dimIndex);
359                 pmeSetSplineData(pmeSafe.get(),
360                                  codePath,
361                                  splineData.splineDerivatives[dimIndex],
362                                  PmeSplineDataType::Derivatives,
363                                  dimIndex);
364             }
365
366             /* Explicitly copying the sample forces to be able to modify them */
367             auto inputForcesFull(c_sampleForcesFull);
368             GMX_RELEASE_ASSERT(ssize(inputForcesFull) >= atomCount, "Bad input forces size");
369             auto forces = ForcesVector(inputForcesFull).subArray(0, atomCount);
370
371             /* Running the force gathering itself */
372             pmePerformGather(pmeSafe.get(), codePath, forces);
373             pmeFinalizeTest(pmeSafe.get(), codePath);
374
375             /* Check the output forces correctness */
376             TestReferenceChecker forceChecker(refData.rootChecker());
377             const auto           ulpTolerance = 3 * pmeOrder;
378             forceChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpTolerance));
379             forceChecker.checkSequence(forces.begin(), forces.end(), "Forces");
380         }
381     }
382
383     static std::vector<std::unique_ptr<PmeTestHardwareContext>> s_pmeTestHardwareContexts;
384 };
385
386 std::vector<std::unique_ptr<PmeTestHardwareContext>> GatherTest::s_pmeTestHardwareContexts;
387
388 //! Test for PME force gathering
389 TEST_P(GatherTest, WorksWith)
390 {
391     EXPECT_NO_THROW_GMX(runTest());
392 }
393
394 //! Moved out from instantiations for readability
395 const auto c_inputBoxNames = ::testing::Values("rect", "tric");
396 //! Moved out from instantiations for readability
397 const auto c_inputGridNames = ::testing::Values("first", "second");
398 //! Moved out from instantiations for readability
399 const auto c_inputTestSystemNames = ::testing::Values("1 atom", "2 atoms", "13 atoms");
400
401 //! Instantiation of the PME gathering test
402 INSTANTIATE_TEST_SUITE_P(Pme,
403                          GatherTest,
404                          ::testing::Combine(c_inputBoxNames,
405                                             ::testing::ValuesIn(c_inputPmeOrders),
406                                             ::testing::ValuesIn(c_inputGridSizes),
407                                             c_inputGridNames,
408                                             c_inputTestSystemNames),
409                          nameOfTest);
410
411 } // namespace
412 } // namespace test
413 } // namespace gmx