2 * This file is part of the GROMACS molecular simulation package.
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.
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 PmeSplineAndSpreadTest() = default;
88 /* Getting the input */
92 CoordinatesVector coordinates;
93 ChargesVector charges;
95 std::tie(box, pmeOrder, gridSize, coordinates, charges) = GetParam();
96 const size_t atomCount = coordinates.size();
98 /* Storing the input where it's needed */
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;
107 TestReferenceData refData;
109 const std::map<PmeSplineAndSpreadOptions, std::string> optionsToTest = {{PmeSplineAndSpreadOptions::SplineAndSpreadUnified, "spline computation and charge spreading (fused)"},
110 {PmeSplineAndSpreadOptions::SplineOnly, "spline computation"},
111 {PmeSplineAndSpreadOptions::SpreadOnly, "charge spreading"}};
113 // There is a subtle problem with multiple comparisons against same reference data:
114 // The subsequent (GPU) spreading runs at one point didn't actually copy the output grid into the proper buffer,
115 // but the reference data was already marked as checked (hasBeenChecked_) by the CPU run, so nothing failed.
116 // For now we will manually track that the count of the grid entries is the same on each run.
117 // This is just a hack for a single specific output though.
118 // What would be much better TODO is to split different codepaths into separate tests,
119 // while making them use the same reference files.
120 bool gridValuesSizeAssigned = false;
121 size_t previousGridValuesSize;
123 for (const auto &context : getPmeTestEnv()->getHardwareContexts())
125 CodePath codePath = context->getCodePath();
126 const bool supportedInput = pmeSupportsInputForMode(*getPmeTestEnv()->hwinfo(), &inputRec, codePath);
129 /* Testing the failure for the unsupported input */
130 EXPECT_THROW(pmeInitAtoms(&inputRec, codePath, nullptr, nullptr, coordinates, charges, box), NotImplementedError);
134 for (const auto &option : optionsToTest)
136 /* Describing the test uniquely in case it fails */
138 SCOPED_TRACE(formatString("Testing %s with %s %sfor PME grid size %d %d %d"
139 ", order %d, %zu atoms",
140 option.second.c_str(), codePathToString(codePath),
141 context->getDescription().c_str(),
142 gridSize[XX], gridSize[YY], gridSize[ZZ],
146 /* Running the test */
148 PmeSafePointer pmeSafe = pmeInitAtoms(&inputRec, codePath, context->getDeviceInfo(),
149 context->getPmeGpuProgram(), coordinates, charges, box);
151 const bool computeSplines = (option.first == PmeSplineAndSpreadOptions::SplineOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
152 const bool spreadCharges = (option.first == PmeSplineAndSpreadOptions::SpreadOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
156 // Here we should set up the results of the spline computation so that the spread can run.
157 // What is lazy and works is running the separate spline so that it will set it up for us:
158 pmePerformSplineAndSpread(pmeSafe.get(), codePath, true, false);
159 // We know that it is tested in another iteration.
160 // TODO: Clean alternative: read and set the reference gridline indices, spline params
163 pmePerformSplineAndSpread(pmeSafe.get(), codePath, computeSplines, spreadCharges);
164 pmeFinalizeTest(pmeSafe.get(), codePath);
166 /* Outputs correctness check */
167 /* All tolerances were picked empirically for single precision on CPU */
169 TestReferenceChecker rootChecker(refData.rootChecker());
171 const auto maxGridSize = std::max(std::max(gridSize[XX], gridSize[YY]), gridSize[ZZ]);
172 const auto ulpToleranceSplineValues = 4 * (pmeOrder - 2) * maxGridSize;
173 /* 4 is a modest estimate for amount of operations; (pmeOrder - 2) is a number of iterations;
174 * maxGridSize is inverse of the smallest positive fractional coordinate (which are interpolated by the splines).
179 const char *dimString[] = { "X", "Y", "Z" };
182 SCOPED_TRACE(formatString("Testing spline values with tolerance of %d", ulpToleranceSplineValues));
183 TestReferenceChecker splineValuesChecker(rootChecker.checkCompound("Splines", "Values"));
184 splineValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineValues));
185 for (int i = 0; i < DIM; i++)
187 auto splineValuesDim = pmeGetSplineData(pmeSafe.get(), codePath, PmeSplineDataType::Values, i);
188 splineValuesChecker.checkSequence(splineValuesDim.begin(), splineValuesDim.end(), dimString[i]);
191 /* Spline derivatives */
192 const auto ulpToleranceSplineDerivatives = 4 * ulpToleranceSplineValues;
193 /* 4 is just a wild guess since the derivatives are deltas of neighbor spline values which could differ greatly */
194 SCOPED_TRACE(formatString("Testing spline derivatives with tolerance of %d", ulpToleranceSplineDerivatives));
195 TestReferenceChecker splineDerivativesChecker(rootChecker.checkCompound("Splines", "Derivatives"));
196 splineDerivativesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineDerivatives));
197 for (int i = 0; i < DIM; i++)
199 auto splineDerivativesDim = pmeGetSplineData(pmeSafe.get(), codePath, PmeSplineDataType::Derivatives, i);
200 splineDerivativesChecker.checkSequence(splineDerivativesDim.begin(), splineDerivativesDim.end(), dimString[i]);
203 /* Particle gridline indices */
204 auto gridLineIndices = pmeGetGridlineIndices(pmeSafe.get(), codePath);
205 rootChecker.checkSequence(gridLineIndices.begin(), gridLineIndices.end(), "Gridline indices");
210 /* The wrapped grid */
211 SparseRealGridValuesOutput nonZeroGridValues = pmeGetRealGrid(pmeSafe.get(), codePath);
212 TestReferenceChecker gridValuesChecker(rootChecker.checkCompound("NonZeroGridValues", "RealSpaceGrid"));
213 const auto ulpToleranceGrid = 2 * ulpToleranceSplineValues * static_cast<int>(ceil(sqrt(atomCount)));
214 /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */
215 SCOPED_TRACE(formatString("Testing grid values with tolerance of %d", ulpToleranceGrid));
216 if (!gridValuesSizeAssigned)
218 previousGridValuesSize = nonZeroGridValues.size();
219 gridValuesSizeAssigned = true;
223 EXPECT_EQ(previousGridValuesSize, nonZeroGridValues.size());
226 gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
227 for (const auto &point : nonZeroGridValues)
229 gridValuesChecker.checkReal(point.second, point.first.c_str());
238 /*! \brief Test for spline parameter computation and charge spreading. */
239 TEST_P(PmeSplineAndSpreadTest, ReproducesOutputs)
241 EXPECT_NO_THROW(runTest());
244 /* Valid input instances */
246 //! A couple of valid inputs for boxes.
247 std::vector<Matrix3x3> const c_sampleBoxes
263 //! A couple of valid inputs for grid sizes.
264 std::vector<IVec> const c_sampleGridSizes
275 std::vector<real> const c_sampleChargesFull
277 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
280 auto const c_sampleCharges1 = ChargesVector(c_sampleChargesFull).subArray(0, 1);
282 auto const c_sampleCharges2 = ChargesVector(c_sampleChargesFull).subArray(1, 2);
284 auto const c_sampleCharges13 = ChargesVector(c_sampleChargesFull).subArray(3, 13);
286 //! Random coordinate vectors
287 CoordinatesVector const c_sampleCoordinatesFull
292 16.0f, 1.02f, 0.22f // 2 box lengths in x
302 16.0001f, 1.52f, 1.19f // > 2 box lengths in x
304 1.43f, 1.1f, 4.1f // > 2 box lengths in z
306 -1.08f, 1.19f, 0.08f // negative x
310 1.32f, -1.48f, 0.16f // negative y
314 0.95f, 7.7f, -0.48f // > 2 box lengths in y, negative z
325 //! 1 coordinate vector
326 CoordinatesVector const c_sampleCoordinates1(c_sampleCoordinatesFull.begin(), c_sampleCoordinatesFull.begin() + 1);
327 //! 2 coordinate vectors
328 CoordinatesVector const c_sampleCoordinates2(c_sampleCoordinatesFull.begin() + 1, c_sampleCoordinatesFull.begin() + 3);
329 //! 13 coordinate vectors
330 CoordinatesVector const c_sampleCoordinates13(c_sampleCoordinatesFull.begin() + 3, c_sampleCoordinatesFull.begin() + 16);
332 //! moved out from instantiantions for readability
333 auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
334 //! moved out from instantiantions for readability
335 auto c_inputPmeOrders = ::testing::Range(3, 5 + 1);
336 //! moved out from instantiantions for readability
337 auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
339 /*! \brief Instantiation of the test with valid input and 1 atom */
340 INSTANTIATE_TEST_CASE_P(SaneInput1, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
341 ::testing::Values(c_sampleCoordinates1),
342 ::testing::Values(c_sampleCharges1)
345 /*! \brief Instantiation of the test with valid input and 2 atoms */
346 INSTANTIATE_TEST_CASE_P(SaneInput2, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
347 ::testing::Values(c_sampleCoordinates2),
348 ::testing::Values(c_sampleCharges2)
350 /*! \brief Instantiation of the test with valid input and 13 atoms */
351 INSTANTIATE_TEST_CASE_P(SaneInput13, PmeSplineAndSpreadTest, ::testing::Combine(c_inputBoxes, c_inputPmeOrders, c_inputGridSizes,
352 ::testing::Values(c_sampleCoordinates13),
353 ::testing::Values(c_sampleCharges13)