Merge branch release-2020 into merge-2020-into-2021
[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, 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 SetUpTestCase() { s_pmeTestHardwareContexts = createPmeTestHardwareContextList(); }
89
90     //! The test
91     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 = eelPME;
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(formatString(
149                         "Testing %s on %s for PME grid size %d %d %d"
150                         ", order %d, %zu atoms",
151                         option.second.c_str(), pmeTestHardwareContext->description().c_str(),
152                         gridSize[XX], gridSize[YY], gridSize[ZZ], pmeOrder, atomCount));
153
154                 /* Running the test */
155
156                 PmeSafePointer pmeSafe =
157                         pmeInitWrapper(&inputRec, codePath, pmeTestHardwareContext->deviceContext(),
158                                        pmeTestHardwareContext->deviceStream(),
159                                        pmeTestHardwareContext->pmeGpuProgram(), box);
160                 std::unique_ptr<StatePropagatorDataGpu> stateGpu =
161                         (codePath == CodePath::GPU)
162                                 ? makeStatePropagatorDataGpu(*pmeSafe.get(),
163                                                              pmeTestHardwareContext->deviceContext(),
164                                                              pmeTestHardwareContext->deviceStream())
165                                 : nullptr;
166
167                 pmeInitAtoms(pmeSafe.get(), stateGpu.get(), codePath, coordinates, charges);
168
169                 const bool computeSplines =
170                         (option.first == PmeSplineAndSpreadOptions::SplineOnly)
171                         || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
172                 const bool spreadCharges =
173                         (option.first == PmeSplineAndSpreadOptions::SpreadOnly)
174                         || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
175
176                 if (!computeSplines)
177                 {
178                     // Here we should set up the results of the spline computation so that the spread can run.
179                     // What is lazy and works is running the separate spline so that it will set it up for us:
180                     pmePerformSplineAndSpread(pmeSafe.get(), codePath, true, false);
181                     // We know that it is tested in another iteration.
182                     // TODO: Clean alternative: read and set the reference gridline indices, spline params
183                 }
184
185                 pmePerformSplineAndSpread(pmeSafe.get(), codePath, computeSplines, spreadCharges);
186                 pmeFinalizeTest(pmeSafe.get(), codePath);
187
188                 /* Outputs correctness check */
189                 /* All tolerances were picked empirically for single precision on CPU */
190
191                 TestReferenceChecker rootChecker(refData.rootChecker());
192
193                 const auto maxGridSize = std::max(std::max(gridSize[XX], gridSize[YY]), gridSize[ZZ]);
194                 const auto ulpToleranceSplineValues = 4 * (pmeOrder - 2) * maxGridSize;
195                 /* 4 is a modest estimate for amount of operations; (pmeOrder - 2) is a number of iterations;
196                  * maxGridSize is inverse of the smallest positive fractional coordinate (which are interpolated by the splines).
197                  */
198
199                 if (computeSplines)
200                 {
201                     const char* dimString[] = { "X", "Y", "Z" };
202
203                     /* Spline values */
204                     SCOPED_TRACE(formatString("Testing spline values with tolerance of %d",
205                                               ulpToleranceSplineValues));
206                     TestReferenceChecker splineValuesChecker(
207                             rootChecker.checkCompound("Splines", "Values"));
208                     splineValuesChecker.setDefaultTolerance(
209                             relativeToleranceAsUlp(1.0, ulpToleranceSplineValues));
210                     for (int i = 0; i < DIM; i++)
211                     {
212                         auto splineValuesDim =
213                                 pmeGetSplineData(pmeSafe.get(), codePath, PmeSplineDataType::Values, i);
214                         splineValuesChecker.checkSequence(splineValuesDim.begin(),
215                                                           splineValuesDim.end(), dimString[i]);
216                     }
217
218                     /* Spline derivatives */
219                     const auto ulpToleranceSplineDerivatives = 4 * ulpToleranceSplineValues;
220                     /* 4 is just a wild guess since the derivatives are deltas of neighbor spline values which could differ greatly */
221                     SCOPED_TRACE(formatString("Testing spline derivatives with tolerance of %d",
222                                               ulpToleranceSplineDerivatives));
223                     TestReferenceChecker splineDerivativesChecker(
224                             rootChecker.checkCompound("Splines", "Derivatives"));
225                     splineDerivativesChecker.setDefaultTolerance(
226                             relativeToleranceAsUlp(1.0, ulpToleranceSplineDerivatives));
227                     for (int i = 0; i < DIM; i++)
228                     {
229                         auto splineDerivativesDim = pmeGetSplineData(
230                                 pmeSafe.get(), codePath, PmeSplineDataType::Derivatives, i);
231                         splineDerivativesChecker.checkSequence(
232                                 splineDerivativesDim.begin(), splineDerivativesDim.end(), dimString[i]);
233                     }
234
235                     /* Particle gridline indices */
236                     auto gridLineIndices = pmeGetGridlineIndices(pmeSafe.get(), codePath);
237                     rootChecker.checkSequence(gridLineIndices.begin(), gridLineIndices.end(),
238                                               "Gridline indices");
239                 }
240
241                 if (spreadCharges)
242                 {
243                     /* The wrapped grid */
244                     SparseRealGridValuesOutput nonZeroGridValues = pmeGetRealGrid(pmeSafe.get(), codePath);
245                     TestReferenceChecker gridValuesChecker(
246                             rootChecker.checkCompound("NonZeroGridValues", "RealSpaceGrid"));
247                     const auto ulpToleranceGrid =
248                             2 * ulpToleranceSplineValues
249                             * static_cast<int>(ceil(sqrt(static_cast<real>(atomCount))));
250                     /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */
251                     SCOPED_TRACE(formatString("Testing grid values with tolerance of %d", ulpToleranceGrid));
252                     if (!gridValuesSizeAssigned)
253                     {
254                         previousGridValuesSize = nonZeroGridValues.size();
255                         gridValuesSizeAssigned = true;
256                     }
257                     else
258                     {
259                         EXPECT_EQ(previousGridValuesSize, nonZeroGridValues.size());
260                     }
261
262                     gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
263                     for (const auto& point : nonZeroGridValues)
264                     {
265                         gridValuesChecker.checkReal(point.second, point.first.c_str());
266                     }
267                 }
268             }
269         }
270     }
271
272     static std::vector<std::unique_ptr<PmeTestHardwareContext>> s_pmeTestHardwareContexts;
273 };
274
275 std::vector<std::unique_ptr<PmeTestHardwareContext>> PmeSplineAndSpreadTest::s_pmeTestHardwareContexts;
276
277
278 /*! \brief Test for spline parameter computation and charge spreading. */
279 TEST_P(PmeSplineAndSpreadTest, ReproducesOutputs)
280 {
281     EXPECT_NO_THROW_GMX(runTest());
282 }
283
284 /* Valid input instances */
285
286 //! A couple of valid inputs for boxes.
287 std::vector<Matrix3x3> const c_sampleBoxes{
288     // normal box
289     Matrix3x3{ { 8.0F, 0.0F, 0.0F, 0.0F, 3.4F, 0.0F, 0.0F, 0.0F, 2.0F } },
290     // triclinic box
291     Matrix3x3{ { 7.0F, 0.0F, 0.0F, 0.0F, 4.1F, 0.0F, 3.5F, 2.0F, 12.2F } },
292 };
293
294 //! A couple of valid inputs for grid sizes.
295 std::vector<IVec> const c_sampleGridSizes{ IVec{ 16, 12, 14 }, IVec{ 19, 17, 11 } };
296
297 //! Random charges
298 std::vector<real> const c_sampleChargesFull{ 4.95F, 3.11F, 3.97F, 1.08F, 2.09F, 1.1F,
299                                              4.13F, 3.31F, 2.8F,  5.83F, 5.09F, 6.1F,
300                                              2.86F, 0.24F, 5.76F, 5.19F, 0.72F };
301 //! 1 charge
302 auto const c_sampleCharges1 = ChargesVector(c_sampleChargesFull).subArray(0, 1);
303 //! 2 charges
304 auto const c_sampleCharges2 = ChargesVector(c_sampleChargesFull).subArray(1, 2);
305 //! 13 charges
306 auto const c_sampleCharges13 = ChargesVector(c_sampleChargesFull).subArray(3, 13);
307
308 //! Random coordinate vectors
309 CoordinatesVector const c_sampleCoordinatesFull{ { 5.59F, 1.37F, 0.95F },
310                                                  {
311                                                          16.0F, 1.02F, 0.22F // 2 box lengths in x
312                                                  },
313                                                  { 0.034F, 1.65F, 0.22F },
314                                                  { 0.33F, 0.92F, 1.56F },
315                                                  { 1.16F, 0.75F, 0.39F },
316                                                  { 0.5F, 1.63F, 1.14F },
317                                                  {
318                                                          16.0001F, 1.52F, 1.19F // > 2 box lengths in x
319                                                  },
320                                                  {
321                                                          1.43F, 1.1F, 4.1F // > 2 box lengths in z
322                                                  },
323                                                  {
324                                                          -1.08F, 1.19F, 0.08F // negative x
325                                                  },
326                                                  { 1.6F, 0.93F, 0.53F },
327                                                  {
328                                                          1.32F, -1.48F, 0.16F // negative y
329                                                  },
330                                                  { 0.87F, 0.0F, 0.33F },
331                                                  {
332                                                          0.95F, 7.7F, -0.48F // > 2 box lengths in y, negative z
333                                                  },
334                                                  { 1.23F, 0.91F, 0.68F },
335                                                  { 0.19F, 1.45F, 0.94F },
336                                                  { 1.28F, 0.46F, 0.38F },
337                                                  { 1.21F, 0.23F, 1.0F } };
338 //! 1 coordinate vector
339 CoordinatesVector const c_sampleCoordinates1(c_sampleCoordinatesFull.begin(),
340                                              c_sampleCoordinatesFull.begin() + 1);
341 //! 2 coordinate vectors
342 CoordinatesVector const c_sampleCoordinates2(c_sampleCoordinatesFull.begin() + 1,
343                                              c_sampleCoordinatesFull.begin() + 3);
344 //! 13 coordinate vectors
345 CoordinatesVector const c_sampleCoordinates13(c_sampleCoordinatesFull.begin() + 3,
346                                               c_sampleCoordinatesFull.begin() + 16);
347
348 //! moved out from instantiantions for readability
349 auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
350 //! moved out from instantiantions for readability
351 auto c_inputPmeOrders = ::testing::Range(3, 5 + 1);
352 //! moved out from instantiantions for readability
353 auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
354
355 /*! \brief Instantiation of the test with valid input and 1 atom */
356 INSTANTIATE_TEST_CASE_P(SaneInput1,
357                         PmeSplineAndSpreadTest,
358                         ::testing::Combine(c_inputBoxes,
359                                            c_inputPmeOrders,
360                                            c_inputGridSizes,
361                                            ::testing::Values(c_sampleCoordinates1),
362                                            ::testing::Values(c_sampleCharges1)));
363
364 /*! \brief Instantiation of the test with valid input and 2 atoms */
365 INSTANTIATE_TEST_CASE_P(SaneInput2,
366                         PmeSplineAndSpreadTest,
367                         ::testing::Combine(c_inputBoxes,
368                                            c_inputPmeOrders,
369                                            c_inputGridSizes,
370                                            ::testing::Values(c_sampleCoordinates2),
371                                            ::testing::Values(c_sampleCharges2)));
372 /*! \brief Instantiation of the test with valid input and 13 atoms */
373 INSTANTIATE_TEST_CASE_P(SaneInput13,
374                         PmeSplineAndSpreadTest,
375                         ::testing::Combine(c_inputBoxes,
376                                            c_inputPmeOrders,
377                                            c_inputGridSizes,
378                                            ::testing::Values(c_sampleCoordinates13),
379                                            ::testing::Values(c_sampleCharges13)));
380 } // namespace
381 } // namespace test
382 } // namespace gmx