clang-tidy: google tests applicable
[alexxy/gromacs.git] / src / gromacs / ewald / tests / pmesplinespreadtest.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 spline computation and charge spreading 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 //! PME spline and spread code path being tested
65 enum class PmeSplineAndSpreadOptions
66 {
67     SplineOnly,
68     SpreadOnly,
69     SplineAndSpreadUnified
70 };
71
72 /*! \brief Convenience typedef of input parameters - unit cell box, PME interpolation order, grid dimensions,
73  * particle coordinates, particle charges
74  * TODO: consider inclusion of local grid offsets/sizes or PME nodes counts to test the PME DD
75  */
76 typedef std::tuple<Matrix3x3, int, IVec, CoordinatesVector, ChargesVector> SplineAndSpreadInputParameters;
77
78 /*! \brief Test fixture for testing both atom spline parameter computation and charge spreading.
79  * These 2 stages of PME are tightly coupled in the code.
80  */
81 class PmeSplineAndSpreadTest : public ::testing::TestWithParam<SplineAndSpreadInputParameters>
82 {
83     public:
84         //! Default constructor
85         PmeSplineAndSpreadTest() = default;
86         //! The test
87         void runTest()
88         {
89             /* Getting the input */
90             Matrix3x3               box;
91             int                     pmeOrder;
92             IVec                    gridSize;
93             CoordinatesVector       coordinates;
94             ChargesVector           charges;
95
96             std::tie(box, pmeOrder, gridSize, coordinates, charges) = GetParam();
97             const size_t atomCount = coordinates.size();
98
99             /* Storing the input where it's needed */
100             t_inputrec inputRec;
101             inputRec.nkx         = gridSize[XX];
102             inputRec.nky         = gridSize[YY];
103             inputRec.nkz         = gridSize[ZZ];
104             inputRec.pme_order   = pmeOrder;
105             inputRec.coulombtype = eelPME;
106             inputRec.epsilon_r   = 1.0;
107
108             TestReferenceData                                      refData;
109
110             const std::map<PmeSplineAndSpreadOptions, std::string> optionsToTest = {{PmeSplineAndSpreadOptions::SplineAndSpreadUnified, "spline computation and charge spreading (fused)"},
111                                                                                     {PmeSplineAndSpreadOptions::SplineOnly, "spline computation"},
112                                                                                     {PmeSplineAndSpreadOptions::SpreadOnly, "charge spreading"}};
113
114             // There is a subtle problem with multiple comparisons against same reference data:
115             // The subsequent (GPU) spreading runs at one point didn't actually copy the output grid into the proper buffer,
116             // but the reference data was already marked as checked (hasBeenChecked_) by the CPU run, so nothing failed.
117             // For now we will manually track that the count of the grid entries is the same on each run.
118             // This is just a hack for a single specific output though.
119             // What would be much better TODO is to split different codepaths into separate tests,
120             // while making them use the same reference files.
121             bool       gridValuesSizeAssigned = false;
122             size_t     previousGridValuesSize;
123
124             for (const auto &context : getPmeTestEnv()->getHardwareContexts())
125             {
126                 CodePath   codePath       = context->getCodePath();
127                 const bool supportedInput = pmeSupportsInputForMode(&inputRec, codePath);
128                 if (!supportedInput)
129                 {
130                     /* Testing the failure for the unsupported input */
131                     EXPECT_THROW(pmeInitAtoms(&inputRec, codePath, nullptr, nullptr, coordinates, charges, box), NotImplementedError);
132                     continue;
133                 }
134
135                 for (const auto &option : optionsToTest)
136                 {
137                     /* Describing the test uniquely in case it fails */
138
139                     SCOPED_TRACE(formatString("Testing %s with %s %sfor PME grid size %d %d %d"
140                                               ", order %d, %zu atoms",
141                                               option.second.c_str(), codePathToString(codePath),
142                                               context->getDescription().c_str(),
143                                               gridSize[XX], gridSize[YY], gridSize[ZZ],
144                                               pmeOrder,
145                                               atomCount));
146
147                     /* Running the test */
148
149                     PmeSafePointer pmeSafe = pmeInitAtoms(&inputRec, codePath, context->getDeviceInfo(),
150                                                           context->getPmeGpuProgram(), coordinates, charges, box);
151
152                     const bool     computeSplines = (option.first == PmeSplineAndSpreadOptions::SplineOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
153                     const bool     spreadCharges  = (option.first == PmeSplineAndSpreadOptions::SpreadOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
154
155                     if (!computeSplines)
156                     {
157                         // Here we should set up the results of the spline computation so that the spread can run.
158                         // What is lazy and works is running the separate spline so that it will set it up for us:
159                         pmePerformSplineAndSpread(pmeSafe.get(), codePath, true, false);
160                         // We know that it is tested in another iteration.
161                         // TODO: Clean alternative: read and set the reference gridline indices, spline params
162                     }
163
164                     pmePerformSplineAndSpread(pmeSafe.get(), codePath, computeSplines, spreadCharges);
165                     pmeFinalizeTest(pmeSafe.get(), codePath);
166
167                     /* Outputs correctness check */
168                     /* All tolerances were picked empirically for single precision on CPU */
169
170                     TestReferenceChecker rootChecker(refData.rootChecker());
171
172                     const auto           maxGridSize              = std::max(std::max(gridSize[XX], gridSize[YY]), gridSize[ZZ]);
173                     const auto           ulpToleranceSplineValues = 4 * (pmeOrder - 2) * maxGridSize;
174                     /* 4 is a modest estimate for amount of operations; (pmeOrder - 2) is a number of iterations;
175                      * maxGridSize is inverse of the smallest positive fractional coordinate (which are interpolated by the splines).
176                      */
177
178                     if (computeSplines)
179                     {
180                         const char *dimString[] = { "X", "Y", "Z" };
181
182                         /* Spline values */
183                         SCOPED_TRACE(formatString("Testing spline values with tolerance of %ld", ulpToleranceSplineValues));
184                         TestReferenceChecker splineValuesChecker(rootChecker.checkCompound("Splines", "Values"));
185                         splineValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineValues));
186                         for (int i = 0; i < DIM; i++)
187                         {
188                             auto splineValuesDim = pmeGetSplineData(pmeSafe.get(), codePath, PmeSplineDataType::Values, i);
189                             splineValuesChecker.checkSequence(splineValuesDim.begin(), splineValuesDim.end(), dimString[i]);
190                         }
191
192                         /* Spline derivatives */
193                         const auto ulpToleranceSplineDerivatives = 4 * ulpToleranceSplineValues;
194                         /* 4 is just a wild guess since the derivatives are deltas of neighbor spline values which could differ greatly */
195                         SCOPED_TRACE(formatString("Testing spline derivatives with tolerance of %ld", ulpToleranceSplineDerivatives));
196                         TestReferenceChecker splineDerivativesChecker(rootChecker.checkCompound("Splines", "Derivatives"));
197                         splineDerivativesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineDerivatives));
198                         for (int i = 0; i < DIM; i++)
199                         {
200                             auto splineDerivativesDim = pmeGetSplineData(pmeSafe.get(), codePath, PmeSplineDataType::Derivatives, i);
201                             splineDerivativesChecker.checkSequence(splineDerivativesDim.begin(), splineDerivativesDim.end(), dimString[i]);
202                         }
203
204                         /* Particle gridline indices */
205                         auto gridLineIndices = pmeGetGridlineIndices(pmeSafe.get(), codePath);
206                         rootChecker.checkSequence(gridLineIndices.begin(), gridLineIndices.end(), "Gridline indices");
207                     }
208
209                     if (spreadCharges)
210                     {
211                         /* The wrapped grid */
212                         SparseRealGridValuesOutput nonZeroGridValues = pmeGetRealGrid(pmeSafe.get(), codePath);
213                         TestReferenceChecker       gridValuesChecker(rootChecker.checkCompound("NonZeroGridValues", "RealSpaceGrid"));
214                         const auto                 ulpToleranceGrid = 2 * ulpToleranceSplineValues * static_cast<int>(ceil(sqrt(atomCount)));
215                         /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */
216                         SCOPED_TRACE(formatString("Testing grid values with tolerance of %ld", ulpToleranceGrid));
217                         if (!gridValuesSizeAssigned)
218                         {
219                             previousGridValuesSize = nonZeroGridValues.size();
220                             gridValuesSizeAssigned = true;
221                         }
222                         else
223                         {
224                             EXPECT_EQ(previousGridValuesSize, nonZeroGridValues.size());
225                         }
226
227                         gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
228                         for (const auto &point : nonZeroGridValues)
229                         {
230                             gridValuesChecker.checkReal(point.second, point.first.c_str());
231                         }
232                     }
233                 }
234             }
235         }
236 };
237
238
239 /*! \brief Test for spline parameter computation and charge spreading. */
240 TEST_P(PmeSplineAndSpreadTest, ReproducesOutputs)
241 {
242     EXPECT_NO_THROW(runTest());
243 }
244
245 /* Valid input instances */
246
247 //! A couple of valid inputs for boxes.
248 std::vector<Matrix3x3> const c_sampleBoxes
249 {
250     // normal box
251     Matrix3x3 {{
252                    8.0f, 0.0f, 0.0f,
253                    0.0f, 3.4f, 0.0f,
254                    0.0f, 0.0f, 2.0f
255                }},
256     // triclinic box
257     Matrix3x3 {{
258                    7.0f, 0.0f, 0.0f,
259                    0.0f, 4.1f, 0.0f,
260                    3.5f, 2.0f, 12.2f
261                }},
262 };
263
264 //! A couple of valid inputs for grid sizes.
265 std::vector<IVec> const c_sampleGridSizes
266 {
267     IVec {
268         16, 12, 14
269     },
270     IVec {
271         19, 17, 11
272     }
273 };
274
275 //! Random charges
276 std::vector<real> const c_sampleChargesFull
277 {
278     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
279 };
280 //! 1 charge
281 auto const c_sampleCharges1 = ChargesVector(c_sampleChargesFull).subArray(0, 1);
282 //! 2 charges
283 auto const c_sampleCharges2 = ChargesVector(c_sampleChargesFull).subArray(1, 2);
284 //! 13 charges
285 auto const c_sampleCharges13 = ChargesVector(c_sampleChargesFull).subArray(3, 13);
286
287 //! Random coordinate vectors
288 CoordinatesVector const c_sampleCoordinatesFull
289 {
290     {
291         5.59f, 1.37f, 0.95f
292     }, {
293         16.0f, 1.02f, 0.22f      // 2 box lengths in x
294     }, {
295         0.034f, 1.65f, 0.22f
296     }, {
297         0.33f, 0.92f, 1.56f
298     }, {
299         1.16f, 0.75f, 0.39f
300     }, {
301         0.5f, 1.63f, 1.14f
302     }, {
303         16.0001f, 1.52f, 1.19f  // > 2 box lengths in x
304     }, {
305         1.43f, 1.1f, 4.1f       // > 2 box lengths in z
306     }, {
307         -1.08f, 1.19f, 0.08f    // negative x
308     }, {
309         1.6f, 0.93f, 0.53f
310     }, {
311         1.32f, -1.48f, 0.16f    // negative y
312     }, {
313         0.87f, 0.0f, 0.33f
314     }, {
315         0.95f, 7.7f, -0.48f     // > 2 box lengths in y, negative z
316     }, {
317         1.23f, 0.91f, 0.68f
318     }, {
319         0.19f, 1.45f, 0.94f
320     }, {
321         1.28f, 0.46f, 0.38f
322     }, {
323         1.21f, 0.23f, 1.0f
324     }
325 };
326 //! 1 coordinate vector
327 CoordinatesVector const c_sampleCoordinates1(c_sampleCoordinatesFull.begin(), c_sampleCoordinatesFull.begin() + 1);
328 //! 2 coordinate vectors
329 CoordinatesVector const c_sampleCoordinates2(c_sampleCoordinatesFull.begin() + 1, c_sampleCoordinatesFull.begin() + 3);
330 //! 13 coordinate vectors
331 CoordinatesVector const c_sampleCoordinates13(c_sampleCoordinatesFull.begin() + 3, c_sampleCoordinatesFull.begin() + 16);
332
333 //! moved out from instantiantions for readability
334 auto c_inputBoxes     = ::testing::ValuesIn(c_sampleBoxes);
335 //! moved out from instantiantions for readability
336 auto c_inputPmeOrders = ::testing::Range(3, 5 + 1);
337 //! moved out from instantiantions for readability
338 auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
339
340 /*! \brief Instantiation of the test with valid input and 1 atom */
341 INSTANTIATE_TEST_CASE_P(SaneInput1, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
342                                                                                    ::testing::Values(c_sampleCoordinates1),
343                                                                                    ::testing::Values(c_sampleCharges1)
344                                                                                ));
345
346 /*! \brief Instantiation of the test with valid input and 2 atoms */
347 INSTANTIATE_TEST_CASE_P(SaneInput2, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
348                                                                                    ::testing::Values(c_sampleCoordinates2),
349                                                                                    ::testing::Values(c_sampleCharges2)
350                                                                                ));
351 /*! \brief Instantiation of the test with valid input and 13 atoms */
352 INSTANTIATE_TEST_CASE_P(SaneInput13, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
353                                                                                     ::testing::Values(c_sampleCoordinates13),
354                                                                                     ::testing::Values(c_sampleCharges13)
355                                                                                 ));
356 }  // namespace
357 }  // namespace test
358 }  // namespace gmx