6d1286f79429badb3781c6324933137cf627e969
[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, 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  * \ingroup module_ewald
41  */
42
43 #include "gmxpre.h"
44
45 #include <string>
46
47 #include <gmock/gmock.h>
48
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/utility/stringutil.h"
51
52 #include "testutils/refdata.h"
53 #include "testutils/testasserts.h"
54
55 #include "pmetestcommon.h"
56
57 namespace gmx
58 {
59 namespace test
60 {
61 namespace
62 {
63
64 /* Valid input instances */
65
66 //! A couple of valid inputs for boxes.
67 static std::vector<Matrix3x3> const c_sampleBoxes
68 {
69     // normal box
70     Matrix3x3 {{
71                    8.0f, 0.0f, 0.0f,
72                    0.0f, 3.4f, 0.0f,
73                    0.0f, 0.0f, 2.0f
74                }},
75     // triclinic box
76     Matrix3x3 {{
77                    7.0f, 0.0f, 0.0f,
78                    0.0f, 4.1f, 0.0f,
79                    3.5f, 2.0f, 12.2f
80                }},
81 };
82
83 //! A couple of valid inputs for grid sizes
84 static std::vector<IVec> const c_sampleGridSizes
85 {
86     IVec {
87         16, 12, 14
88     },
89     IVec {
90         13, 15, 11
91     }
92 };
93 //! Random charges
94 static std::vector<real> const c_sampleChargesFull
95 {
96     4.95f, 3.11f, 3.97f, 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, 0.72f
97 };
98
99 //! All the input atom gridline indices
100 static std::vector<IVec> const c_sampleGridLineIndicesFull
101 {
102     IVec {
103         4, 2, 6
104     },
105     IVec {
106         1, 4, 10
107     },
108     IVec {
109         0, 6, 6
110     },
111     IVec {
112         0, 1, 4
113     },
114     IVec {
115         6, 3, 0
116     },
117     IVec {
118         7, 2, 2
119     },
120     IVec {
121         8, 3, 1
122     },
123     IVec {
124         4, 0, 3
125     },
126     IVec {
127         0, 0, 0
128     },
129     IVec {
130         8, 5, 8
131     },
132     IVec {
133         4, 4, 2
134     },
135     IVec {
136         7, 1, 7
137     },
138     IVec {
139         8, 5, 5
140     },
141     IVec {
142         2, 6, 5
143     },
144     IVec {
145         1, 6, 2
146     },
147     IVec {
148         7, 1, 8
149     },
150     IVec {
151         3, 5, 1
152     },
153 };
154
155 // Spline values/derivatives below are also generated randomly, so they are bogus,
156 // but that should not affect the reproducibility, which we're after
157
158 //! A lot of bogus input spline values - should have at list (max PME order = 5) * (DIM = 3) * (total unique atom number in all test cases = 16) values
159 static std::vector<real> const c_sampleSplineValuesFull
160 {
161     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, 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, 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, 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, 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, 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, 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, 0.86f, 0.54f, 0.86f, 0.94f, 0.37f, 0.35f, 0.81f, 0.89f, 0.48f,
162     0.34f, 0.18f, 0.11f, 0.02f, 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, 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, 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, 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, 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, 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, 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, 0.21f, 0.43f, 0.32f, 0.44f, 0.72f,
163     0.21f, 0.4f, 0.93f, 0.07f, 0.11f, 0.41f, 0.24f, 0.04f, 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, 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, 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, 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, 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, 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, 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, 0.3f,
164     0.0f, 0.6f, 0.99f, 0.69f,
165 };
166
167 //! A lot of bogus input spline derivatives - should have at list (max PME order = 5) * (DIM = 3) * (total unique atom number in all test cases = 16) values
168 static std::vector<real> const c_sampleSplineDerivativesFull
169 {
170     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, 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, 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, 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, 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, 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, 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, 0.59f, 0.33f, 0.2f, 0.49f, 0.41f, 0.25f, 0.4f, 0.16f, 0.83f,
171     0.44f, 0.82f, 0.21f, 0.95f, 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, 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, 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, 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, 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, 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, 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, 0.48f, 0.37f, 0.93f, 0.5f, 0.2f,
172     0.92f, 0.09f, 0.74f, 0.55f, 0.44f, 0.05f, 0.13f, 0.17f, 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, 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, 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, 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, 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, 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, 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, 0.85f,
173     0.09f, 0.83f, 0.02f, 0.91f,
174 };
175
176 //! 2 c_sample grids - only non-zero values have to be listed
177 static std::vector<SparseRealGridValuesInput> const c_sampleGrids
178 {
179     SparseRealGridValuesInput {{
180                                    IVec {
181                                        0, 0, 0
182                                    }, 3.5f
183                                }, {
184                                    IVec {
185                                        7, 0, 0
186                                    }, -2.5f
187                                },
188                                {
189                                    IVec {
190                                        3, 5, 7
191                                    }, -0.006f
192                                }, {
193                                    IVec {
194                                        1, 6, 7
195                                    }, -2.76f
196                                }, {
197                                    IVec {
198                                        3, 1, 2
199                                    }, 0.6f
200                                },  {
201                                    IVec {
202                                        6, 2, 4
203                                    }, 7.1f
204                                }, {
205                                    IVec {
206                                        4, 5, 6
207                                    }, 4.1f
208                                }, {
209                                    IVec {
210                                        4, 4, 6
211                                    }, -3.7f
212                                }, },
213     SparseRealGridValuesInput {{
214                                    IVec {
215                                        0, 4, 0
216                                    }, 6.f
217                                }, {
218                                    IVec {
219                                        4, 2, 7
220                                    }, 13.76f
221                                }, {
222                                    IVec {
223                                        0, 6, 7
224                                    }, 3.6f
225                                }, {
226                                    IVec {
227                                        3, 2, 8
228                                    }, 0.61f
229                                },
230                                {
231                                    IVec {
232                                        5, 4, 3
233                                    }, 2.1f
234                                },
235                                {
236                                    IVec {
237                                        2, 5, 10
238                                    }, 3.6f
239                                }, {
240                                    IVec {
241                                        5, 3, 6
242                                    }, 2.1f
243                                }, {
244                                    IVec {
245                                        6, 6, 6
246                                    }, 6.1f
247                                }, }
248 };
249
250 //! Input forces for reduction
251 static std::vector<RVec> const c_sampleForcesFull {
252     RVec {
253         0.02f, 0.87f, 0.95f
254     }, RVec {
255         0.66f, 0.67f, 0.38f
256     }, RVec {
257         0.45f, 0.04f, 0.94f
258     }, RVec {
259         0.54f, 0.76f, 0.58f
260     },
261     RVec {
262         0.83f, 0.31f, 0.73f
263     }, RVec {
264         0.71f, 0.06f, 0.35f
265     }, RVec {
266         0.32f, 0.35f, 0.61f
267     }, RVec {
268         0.27f, 0.98f, 0.83f
269     },
270     RVec {
271         0.11f, 0.3f, 0.42f
272     }, RVec {
273         0.95f, 0.69f, 0.58f
274     }, RVec {
275         0.29f, 0.1f, 0.68f
276     }, RVec {
277         0.94f, 0.62f, 0.51f
278     },
279     RVec {
280         0.47f, 0.04f, 0.47f
281     }, RVec {
282         0.34f, 0.71f, 0.52f
283     }
284 };
285
286 //! PME orders to test
287 static std::vector<int> const pmeOrders {
288     3, 4, 5
289 };
290 //! Atom counts to test
291 static std::vector<size_t> const atomCounts {
292     1, 2, 13
293 };
294
295 /* Helper structures for test input */
296
297 //! A structure for all the spline data which depends in size both on the PME order and atom count
298 struct AtomAndPmeOrderSizedData
299 {
300     //! Spline values
301     SplineParamsVector splineValues;
302     //! Spline derivatives
303     SplineParamsVector splineDerivatives;
304 };
305
306 //! A structure for all the input atom data which depends in size on atom count - including a range of spline data for different PME orders
307 struct AtomSizedData
308 {
309     //! Gridline indices
310     GridLineIndicesVector                   gridLineIndices;
311     //! Charges
312     ChargesVector                           charges;
313     //! Coordinates
314     CoordinatesVector                       coordinates;
315     //! Spline data for different orders
316     std::map<int, AtomAndPmeOrderSizedData> splineDataByPmeOrder;
317 };
318
319 //! A range of the test input data sets, uniquely identified by the atom count
320 typedef std::map<size_t, AtomSizedData> InputDataByAtomCount;
321
322 /*! \brief Convenience typedef of the test input parameters - unit cell box, PME interpolation order, grid dimensions,
323  *  grid values, overwriting/reducing the input forces, atom count.
324  *
325  * The rest of the atom-related input data - gridline indices, spline theta values, spline dtheta values, atom charges -
326  * is looked up in the inputAtomDataSets_ test fixture variable.
327  */
328 typedef std::tuple<Matrix3x3, int, IVec, SparseRealGridValuesInput, PmeForceOutputHandling, size_t> GatherInputParameters;
329
330 //! Test fixture
331 class PmeGatherTest : public ::testing::TestWithParam<GatherInputParameters>
332 {
333     private:
334         //! Storage of all the input atom datasets
335         static InputDataByAtomCount s_inputAtomDataSets_;
336
337     public:
338         //! Default constructor
339         PmeGatherTest()  = default;
340         //! Sets the input atom data references once
341         static void SetUpTestCase()
342         {
343             size_t start = 0;
344             for (auto atomCount : atomCounts)
345             {
346                 AtomSizedData atomData;
347                 atomData.gridLineIndices = GridLineIndicesVector(c_sampleGridLineIndicesFull).subArray(start, atomCount);
348                 atomData.charges         = ChargesVector(c_sampleChargesFull).subArray(start, atomCount);
349                 start                   += atomCount;
350                 atomData.coordinates.resize(atomCount, RVec {1e6, 1e7, -1e8});
351                 /* The coordinates are intentionally bogus in this test - only the size matters; the gridline indices are fed directly as inputs */
352                 for (auto pmeOrder : pmeOrders)
353                 {
354                     AtomAndPmeOrderSizedData splineData;
355                     const size_t             dimSize = atomCount * pmeOrder;
356                     for (int dimIndex = 0; dimIndex < DIM; dimIndex++)
357                     {
358                         splineData.splineValues[dimIndex]      = SplineParamsDimVector(c_sampleSplineValuesFull).subArray(dimIndex * dimSize, dimSize);
359                         splineData.splineDerivatives[dimIndex] = SplineParamsDimVector(c_sampleSplineDerivativesFull).subArray(dimIndex * dimSize, dimSize);
360                     }
361                     atomData.splineDataByPmeOrder[pmeOrder] = splineData;
362                 }
363                 s_inputAtomDataSets_[atomCount] = atomData;
364             }
365         }
366         //! The test
367         void runTest()
368         {
369             /* Getting the input */
370             Matrix3x3                 box;
371             int                       pmeOrder;
372             IVec                      gridSize;
373             size_t                    atomCount;
374             SparseRealGridValuesInput nonZeroGridValues;
375             PmeForceOutputHandling    inputForceTreatment;
376             std::tie(box, pmeOrder, gridSize, nonZeroGridValues, inputForceTreatment, atomCount) = GetParam();
377             auto inputAtomData       = s_inputAtomDataSets_[atomCount];
378             auto inputAtomSplineData = inputAtomData.splineDataByPmeOrder[pmeOrder];
379
380             /* Storing the input where it's needed, running the test */
381             t_inputrec inputRec;
382             inputRec.nkx         = gridSize[XX];
383             inputRec.nky         = gridSize[YY];
384             inputRec.nkz         = gridSize[ZZ];
385             inputRec.pme_order   = pmeOrder;
386             inputRec.coulombtype = eelPME;
387             inputRec.epsilon_r   = 1.0;
388
389             TestReferenceData refData;
390             for (const auto &context : getPmeTestEnv()->getHardwareContexts())
391             {
392                 CodePath   codePath       = context->getCodePath();
393                 const bool supportedInput = pmeSupportsInputForMode(&inputRec, codePath);
394                 if (!supportedInput)
395                 {
396                     /* Testing the failure for the unsupported input */
397                     EXPECT_THROW(pmeInitAtoms(&inputRec, codePath, nullptr, nullptr, inputAtomData.coordinates, inputAtomData.charges, box), NotImplementedError);
398                     continue;
399                 }
400
401                 /* Describing the test uniquely */
402                 SCOPED_TRACE(formatString("Testing force gathering with %s %sfor PME grid size %d %d %d"
403                                           ", order %d, %zu atoms, %s",
404                                           codePathToString(codePath), context->getDescription().c_str(),
405                                           gridSize[XX], gridSize[YY], gridSize[ZZ],
406                                           pmeOrder,
407                                           atomCount,
408                                           (inputForceTreatment == PmeForceOutputHandling::ReduceWithInput) ? "with reduction" : "without reduction"
409                                           ));
410
411                 PmeSafePointer pmeSafe = pmeInitAtoms(&inputRec, codePath, context->getDeviceInfo(),
412                                                       context->getPmeGpuProgram(), inputAtomData.coordinates, inputAtomData.charges, box);
413
414                 /* Setting some more inputs */
415                 pmeSetRealGrid(pmeSafe.get(), codePath, nonZeroGridValues);
416                 pmeSetGridLineIndices(pmeSafe.get(), codePath, inputAtomData.gridLineIndices);
417                 for (int dimIndex = 0; dimIndex < DIM; dimIndex++)
418                 {
419                     pmeSetSplineData(pmeSafe.get(), codePath, inputAtomSplineData.splineValues[dimIndex], PmeSplineDataType::Values, dimIndex);
420                     pmeSetSplineData(pmeSafe.get(), codePath, inputAtomSplineData.splineDerivatives[dimIndex], PmeSplineDataType::Derivatives, dimIndex);
421                 }
422
423                 /* Explicitly copying the sample forces to be able to modify them */
424                 auto inputForcesFull(c_sampleForcesFull);
425                 GMX_RELEASE_ASSERT(inputForcesFull.size() >= atomCount, "Bad input forces size");
426                 auto forces = ForcesVector(inputForcesFull).subArray(0, atomCount);
427
428                 /* Running the force gathering itself */
429                 pmePerformGather(pmeSafe.get(), codePath, inputForceTreatment, forces);
430                 pmeFinalizeTest(pmeSafe.get(), codePath);
431
432                 /* Check the output forces correctness */
433                 TestReferenceChecker forceChecker(refData.rootChecker());
434                 const auto           ulpTolerance = 3 * pmeOrder;
435                 forceChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpTolerance));
436                 forceChecker.checkSequence(forces.begin(), forces.end(), "Forces");
437             }
438         }
439 };
440
441 // An instance of static atom data
442 InputDataByAtomCount PmeGatherTest::s_inputAtomDataSets_;
443
444 //! Test for PME force gathering
445 TEST_P(PmeGatherTest, ReproducesOutputs)
446 {
447     EXPECT_NO_THROW(runTest());
448 }
449
450 //! Instantiation of the PME gathering test
451 INSTANTIATE_TEST_CASE_P(SaneInput, PmeGatherTest, ::testing::Combine(::testing::ValuesIn(c_sampleBoxes),
452                                                                          ::testing::ValuesIn(pmeOrders),
453                                                                          ::testing::ValuesIn(c_sampleGridSizes),
454                                                                          ::testing::ValuesIn(c_sampleGrids),
455                                                                          ::testing::Values(PmeForceOutputHandling::Set, PmeForceOutputHandling::ReduceWithInput),
456                                                                          ::testing::ValuesIn(atomCounts)));
457
458 }
459 }
460 }