2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017, 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.
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.
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.
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.
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.
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.
37 * Implements PME spline computation and charge spreading tests.
39 * \author Aleksei Iupinov <a.yupinov@gmail.com>
40 * \ingroup module_ewald
47 #include <gmock/gmock.h>
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/utility/stringutil.h"
52 #include "testutils/refdata.h"
53 #include "testutils/testasserts.h"
55 #include "pmetestcommon.h"
64 //! PME spline and spread code path being tested
65 enum class PmeSplineAndSpreadOptions
69 SplineAndSpreadUnified
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
76 typedef std::tuple<Matrix3x3, int, IVec, CoordinatesVector, ChargesVector> SplineAndSpreadInputParameters;
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.
81 class PmeSplineAndSpreadTest : public ::testing::TestWithParam<SplineAndSpreadInputParameters>
84 //! Default constructor
85 PmeSplineAndSpreadTest() = default;
89 /* Getting the input */
93 CoordinatesVector coordinates;
94 ChargesVector charges;
96 std::tie(box, pmeOrder, gridSize, coordinates, charges) = GetParam();
97 const size_t atomCount = coordinates.size();
99 /* Storing the input where it's needed */
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;
108 TestReferenceData refData;
110 const std::map<CodePath, std::string> modesToTest = {{CodePath::CPU, "CPU"},
111 {CodePath::CUDA, "CUDA"}};
113 const std::map<PmeSplineAndSpreadOptions, std::string> optionsToTest = {{PmeSplineAndSpreadOptions::SplineAndSpreadUnified, "spline computation and charge spreading (fused)"},
114 {PmeSplineAndSpreadOptions::SplineOnly, "spline computation"},
115 {PmeSplineAndSpreadOptions::SpreadOnly, "charge spreading"}};
117 // There is a subtle problem with multiple comparisons against same reference data:
118 // The subsequent (GPU) spreading runs at one point didn't actually copy the output grid into the proper buffer,
119 // but the reference data was already marked as checked (hasBeenChecked_) by the CPU run, so nothing failed.
120 // For now we will manually track that the count of the grid entries is the same on each run.
121 // This is just a hack for a single specific output though.
122 // What would be much better TODO is to split different codepaths into separate tests,
123 // while making them use the same reference files.
124 bool gridValuesSizeAssigned = false;
125 size_t previousGridValuesSize;
127 for (const auto &mode : modesToTest)
129 const bool supportedInput = pmeSupportsInputForMode(&inputRec, mode.first);
132 /* Testing the failure for the unsupported input */
133 EXPECT_THROW(pmeInitAtoms(&inputRec, mode.first, nullptr, coordinates, charges, box), NotImplementedError);
137 const auto contextsToTest = getPmeTestEnv()->getHardwareContexts(mode.first);
138 for (const auto &context : contextsToTest)
140 for (const auto &option : optionsToTest)
142 /* Describing the test uniquely in case it fails */
144 SCOPED_TRACE(formatString("Testing %s with %s %sfor PME grid size %d %d %d"
145 ", order %d, %zu atoms",
146 option.second.c_str(), mode.second.c_str(),
147 context.getDescription().c_str(),
148 gridSize[XX], gridSize[YY], gridSize[ZZ],
152 /* Running the test */
154 PmeSafePointer pmeSafe = pmeInitAtoms(&inputRec, mode.first, context.getDeviceInfo(), coordinates, charges, box);
156 const bool computeSplines = (option.first == PmeSplineAndSpreadOptions::SplineOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
157 const bool spreadCharges = (option.first == PmeSplineAndSpreadOptions::SpreadOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
161 // Here we should set up the results of the spline computation so that the spread can run.
162 // What is lazy and works is running the separate spline so that it will set it up for us:
163 pmePerformSplineAndSpread(pmeSafe.get(), mode.first, true, false);
164 // We know that it is tested in another iteration.
165 // TODO: Clean alternative: read and set the reference gridline indices, spline params
168 pmePerformSplineAndSpread(pmeSafe.get(), mode.first, computeSplines, spreadCharges);
169 pmeFinalizeTest(pmeSafe.get(), mode.first);
171 /* Outputs correctness check */
172 /* All tolerances were picked empirically for single precision on CPU */
174 TestReferenceChecker rootChecker(refData.rootChecker());
176 const auto maxGridSize = std::max(std::max(gridSize[XX], gridSize[YY]), gridSize[ZZ]);
177 const auto ulpToleranceSplineValues = 4 * (pmeOrder - 2) * maxGridSize;
178 /* 4 is a modest estimate for amount of operations; (pmeOrder - 2) is a number of iterations;
179 * maxGridSize is inverse of the smallest positive fractional coordinate (which are interpolated by the splines).
184 const char *dimString[] = { "X", "Y", "Z" };
187 SCOPED_TRACE(formatString("Testing spline values with tolerance of %ld", ulpToleranceSplineValues));
188 TestReferenceChecker splineValuesChecker(rootChecker.checkCompound("Splines", "Values"));
189 splineValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineValues));
190 for (int i = 0; i < DIM; i++)
192 auto splineValuesDim = pmeGetSplineData(pmeSafe.get(), mode.first, PmeSplineDataType::Values, i);
193 splineValuesChecker.checkSequence(splineValuesDim.begin(), splineValuesDim.end(), dimString[i]);
196 /* Spline derivatives */
197 const auto ulpToleranceSplineDerivatives = 4 * ulpToleranceSplineValues;
198 /* 4 is just a wild guess since the derivatives are deltas of neighbor spline values which could differ greatly */
199 SCOPED_TRACE(formatString("Testing spline derivatives with tolerance of %ld", ulpToleranceSplineDerivatives));
200 TestReferenceChecker splineDerivativesChecker(rootChecker.checkCompound("Splines", "Derivatives"));
201 splineDerivativesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineDerivatives));
202 for (int i = 0; i < DIM; i++)
204 auto splineDerivativesDim = pmeGetSplineData(pmeSafe.get(), mode.first, PmeSplineDataType::Derivatives, i);
205 splineDerivativesChecker.checkSequence(splineDerivativesDim.begin(), splineDerivativesDim.end(), dimString[i]);
208 /* Particle gridline indices */
209 auto gridLineIndices = pmeGetGridlineIndices(pmeSafe.get(), mode.first);
210 rootChecker.checkSequence(gridLineIndices.begin(), gridLineIndices.end(), "Gridline indices");
215 /* The wrapped grid */
216 SparseRealGridValuesOutput nonZeroGridValues = pmeGetRealGrid(pmeSafe.get(), mode.first);
217 TestReferenceChecker gridValuesChecker(rootChecker.checkCompound("NonZeroGridValues", "RealSpaceGrid"));
218 const auto ulpToleranceGrid = 2 * ulpToleranceSplineValues * (int)(ceil(sqrt(atomCount)));
219 /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */
220 SCOPED_TRACE(formatString("Testing grid values with tolerance of %ld", ulpToleranceGrid));
221 if (!gridValuesSizeAssigned)
223 previousGridValuesSize = nonZeroGridValues.size();
224 gridValuesSizeAssigned = true;
228 EXPECT_EQ(previousGridValuesSize, nonZeroGridValues.size());
231 gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
232 for (const auto &point : nonZeroGridValues)
234 gridValuesChecker.checkReal(point.second, point.first.c_str());
244 /*! \brief Test for spline parameter computation and charge spreading. */
245 TEST_P(PmeSplineAndSpreadTest, ReproducesOutputs)
247 EXPECT_NO_THROW(runTest());
250 /* Valid input instances */
252 //! A couple of valid inputs for boxes.
253 static std::vector<Matrix3x3> const c_sampleBoxes
269 //! A couple of valid inputs for grid sizes.
270 static std::vector<IVec> const c_sampleGridSizes
281 static std::vector<real> const c_sampleChargesFull
283 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
286 static auto const c_sampleCharges1 = ChargesVector(c_sampleChargesFull).subArray(0, 1);
288 static auto const c_sampleCharges2 = ChargesVector(c_sampleChargesFull).subArray(1, 2);
290 static auto const c_sampleCharges13 = ChargesVector(c_sampleChargesFull).subArray(3, 13);
292 //! Random coordinate vectors
293 static CoordinatesVector const c_sampleCoordinatesFull
298 16.0f, 1.02f, 0.22f // 2 box lengths in x
308 16.0001f, 1.52f, 1.19f // > 2 box lengths in x
310 1.43f, 1.1f, 4.1f // > 2 box lengths in z
312 -1.08f, 1.19f, 0.08f // negative x
316 1.32f, -1.48f, 0.16f // negative y
320 0.95f, 7.7f, -0.48f // > 2 box lengths in y, negative z
331 //! 1 coordinate vector
332 static CoordinatesVector const c_sampleCoordinates1(c_sampleCoordinatesFull.begin(), c_sampleCoordinatesFull.begin() + 1);
333 //! 2 coordinate vectors
334 static CoordinatesVector const c_sampleCoordinates2(c_sampleCoordinatesFull.begin() + 1, c_sampleCoordinatesFull.begin() + 3);
335 //! 13 coordinate vectors
336 static CoordinatesVector const c_sampleCoordinates13(c_sampleCoordinatesFull.begin() + 3, c_sampleCoordinatesFull.begin() + 16);
338 //! moved out from instantiantions for readability
339 auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
340 //! moved out from instantiantions for readability
341 auto c_inputPmeOrders = ::testing::Range(3, 5 + 1);
342 //! moved out from instantiantions for readability
343 auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
345 /*! \brief Instantiation of the test with valid input and 1 atom */
346 INSTANTIATE_TEST_CASE_P(SaneInput1, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
347 ::testing::Values(c_sampleCoordinates1),
348 ::testing::Values(c_sampleCharges1)
351 /*! \brief Instantiation of the test with valid input and 2 atoms */
352 INSTANTIATE_TEST_CASE_P(SaneInput2, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
353 ::testing::Values(c_sampleCoordinates2),
354 ::testing::Values(c_sampleCharges2)
356 /*! \brief Instantiation of the test with valid input and 13 atoms */
357 INSTANTIATE_TEST_CASE_P(SaneInput13, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
358 ::testing::Values(c_sampleCoordinates13),
359 ::testing::Values(c_sampleCharges13)