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