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