2 * This file is part of the GROMACS molecular simulation package.
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.
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"
56 #include "testhardwarecontexts.h"
65 //! PME spline and spread code path being tested
66 enum class PmeSplineAndSpreadOptions
70 SplineAndSpreadUnified
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
77 typedef std::tuple<Matrix3x3, int, IVec, CoordinatesVector, ChargesVector> SplineAndSpreadInputParameters;
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.
82 class PmeSplineAndSpreadTest : public ::testing::TestWithParam<SplineAndSpreadInputParameters>
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<PmeSplineAndSpreadOptions, std::string> optionsToTest = {
111 { PmeSplineAndSpreadOptions::SplineAndSpreadUnified,
112 "spline computation and charge spreading (fused)" },
113 { PmeSplineAndSpreadOptions::SplineOnly, "spline computation" },
114 { 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
119 // into the proper buffer, but the reference data was already marked as checked
120 // (hasBeenChecked_) by the CPU run, so nothing failed. For now we will manually track that
121 // the count of the grid entries is the same on each run. This is just a hack for a single
122 // specific output though. What would be much better TODO is to split different codepaths
123 // into separate tests, while making them use the same reference files.
124 bool gridValuesSizeAssigned = false;
125 size_t previousGridValuesSize;
127 for (const auto& context : getPmeTestEnv()->getHardwareContexts())
129 CodePath codePath = context->getCodePath();
130 const bool supportedInput =
131 pmeSupportsInputForMode(*getPmeTestEnv()->hwinfo(), &inputRec, codePath);
134 /* Testing the failure for the unsupported input */
135 EXPECT_THROW_GMX(pmeInitWrapper(&inputRec, codePath, nullptr, nullptr, box),
136 NotImplementedError);
140 for (const auto& option : optionsToTest)
142 /* Describing the test uniquely in case it fails */
145 formatString("Testing %s with %s %sfor PME grid size %d %d %d"
146 ", order %d, %zu atoms",
147 option.second.c_str(), codePathToString(codePath),
148 context->getDescription().c_str(), gridSize[XX], gridSize[YY],
149 gridSize[ZZ], pmeOrder, atomCount));
151 /* Running the test */
153 PmeSafePointer pmeSafe = pmeInitWrapper(&inputRec, codePath, context->getDeviceInfo(),
154 context->getPmeGpuProgram(), box);
155 std::unique_ptr<StatePropagatorDataGpu> stateGpu =
156 (codePath == CodePath::GPU)
157 ? makeStatePropagatorDataGpu(*pmeSafe.get(), context->deviceContext())
160 pmeInitAtoms(pmeSafe.get(), stateGpu.get(), codePath, coordinates, charges);
162 const bool computeSplines =
163 (option.first == PmeSplineAndSpreadOptions::SplineOnly)
164 || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
165 const bool spreadCharges =
166 (option.first == PmeSplineAndSpreadOptions::SpreadOnly)
167 || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
171 // Here we should set up the results of the spline computation so that the spread can run.
172 // What is lazy and works is running the separate spline so that it will set it up for us:
173 pmePerformSplineAndSpread(pmeSafe.get(), codePath, true, false);
174 // We know that it is tested in another iteration.
175 // TODO: Clean alternative: read and set the reference gridline indices, spline params
178 pmePerformSplineAndSpread(pmeSafe.get(), codePath, computeSplines, spreadCharges);
179 pmeFinalizeTest(pmeSafe.get(), codePath);
181 /* Outputs correctness check */
182 /* All tolerances were picked empirically for single precision on CPU */
184 TestReferenceChecker rootChecker(refData.rootChecker());
186 const auto maxGridSize = std::max(std::max(gridSize[XX], gridSize[YY]), gridSize[ZZ]);
187 const auto ulpToleranceSplineValues = 4 * (pmeOrder - 2) * maxGridSize;
188 /* 4 is a modest estimate for amount of operations; (pmeOrder - 2) is a number of iterations;
189 * maxGridSize is inverse of the smallest positive fractional coordinate (which are interpolated by the splines).
194 const char* dimString[] = { "X", "Y", "Z" };
197 SCOPED_TRACE(formatString("Testing spline values with tolerance of %d",
198 ulpToleranceSplineValues));
199 TestReferenceChecker splineValuesChecker(
200 rootChecker.checkCompound("Splines", "Values"));
201 splineValuesChecker.setDefaultTolerance(
202 relativeToleranceAsUlp(1.0, ulpToleranceSplineValues));
203 for (int i = 0; i < DIM; i++)
205 auto splineValuesDim =
206 pmeGetSplineData(pmeSafe.get(), codePath, PmeSplineDataType::Values, i);
207 splineValuesChecker.checkSequence(splineValuesDim.begin(),
208 splineValuesDim.end(), dimString[i]);
211 /* Spline derivatives */
212 const auto ulpToleranceSplineDerivatives = 4 * ulpToleranceSplineValues;
213 /* 4 is just a wild guess since the derivatives are deltas of neighbor spline values which could differ greatly */
214 SCOPED_TRACE(formatString("Testing spline derivatives with tolerance of %d",
215 ulpToleranceSplineDerivatives));
216 TestReferenceChecker splineDerivativesChecker(
217 rootChecker.checkCompound("Splines", "Derivatives"));
218 splineDerivativesChecker.setDefaultTolerance(
219 relativeToleranceAsUlp(1.0, ulpToleranceSplineDerivatives));
220 for (int i = 0; i < DIM; i++)
222 auto splineDerivativesDim = pmeGetSplineData(
223 pmeSafe.get(), codePath, PmeSplineDataType::Derivatives, i);
224 splineDerivativesChecker.checkSequence(
225 splineDerivativesDim.begin(), splineDerivativesDim.end(), dimString[i]);
228 /* Particle gridline indices */
229 auto gridLineIndices = pmeGetGridlineIndices(pmeSafe.get(), codePath);
230 rootChecker.checkSequence(gridLineIndices.begin(), gridLineIndices.end(),
236 /* The wrapped grid */
237 SparseRealGridValuesOutput nonZeroGridValues = pmeGetRealGrid(pmeSafe.get(), codePath);
238 TestReferenceChecker gridValuesChecker(
239 rootChecker.checkCompound("NonZeroGridValues", "RealSpaceGrid"));
240 const auto ulpToleranceGrid =
241 2 * ulpToleranceSplineValues
242 * static_cast<int>(ceil(sqrt(static_cast<real>(atomCount))));
243 /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */
244 SCOPED_TRACE(formatString("Testing grid values with tolerance of %d", ulpToleranceGrid));
245 if (!gridValuesSizeAssigned)
247 previousGridValuesSize = nonZeroGridValues.size();
248 gridValuesSizeAssigned = true;
252 EXPECT_EQ(previousGridValuesSize, nonZeroGridValues.size());
255 gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
256 for (const auto& point : nonZeroGridValues)
258 gridValuesChecker.checkReal(point.second, point.first.c_str());
267 /*! \brief Test for spline parameter computation and charge spreading. */
268 TEST_P(PmeSplineAndSpreadTest, ReproducesOutputs)
270 EXPECT_NO_THROW_GMX(runTest());
273 /* Valid input instances */
275 //! A couple of valid inputs for boxes.
276 std::vector<Matrix3x3> const c_sampleBoxes{
278 Matrix3x3{ { 8.0F, 0.0F, 0.0F, 0.0F, 3.4F, 0.0F, 0.0F, 0.0F, 2.0F } },
280 Matrix3x3{ { 7.0F, 0.0F, 0.0F, 0.0F, 4.1F, 0.0F, 3.5F, 2.0F, 12.2F } },
283 //! A couple of valid inputs for grid sizes.
284 std::vector<IVec> const c_sampleGridSizes{ IVec{ 16, 12, 14 }, IVec{ 19, 17, 11 } };
287 std::vector<real> const c_sampleChargesFull{ 4.95F, 3.11F, 3.97F, 1.08F, 2.09F, 1.1F,
288 4.13F, 3.31F, 2.8F, 5.83F, 5.09F, 6.1F,
289 2.86F, 0.24F, 5.76F, 5.19F, 0.72F };
291 auto const c_sampleCharges1 = ChargesVector(c_sampleChargesFull).subArray(0, 1);
293 auto const c_sampleCharges2 = ChargesVector(c_sampleChargesFull).subArray(1, 2);
295 auto const c_sampleCharges13 = ChargesVector(c_sampleChargesFull).subArray(3, 13);
297 //! Random coordinate vectors
298 CoordinatesVector const c_sampleCoordinatesFull{ { 5.59F, 1.37F, 0.95F },
300 16.0F, 1.02F, 0.22F // 2 box lengths in x
302 { 0.034F, 1.65F, 0.22F },
303 { 0.33F, 0.92F, 1.56F },
304 { 1.16F, 0.75F, 0.39F },
305 { 0.5F, 1.63F, 1.14F },
307 16.0001F, 1.52F, 1.19F // > 2 box lengths in x
310 1.43F, 1.1F, 4.1F // > 2 box lengths in z
313 -1.08F, 1.19F, 0.08F // negative x
315 { 1.6F, 0.93F, 0.53F },
317 1.32F, -1.48F, 0.16F // negative y
319 { 0.87F, 0.0F, 0.33F },
321 0.95F, 7.7F, -0.48F // > 2 box lengths in y, negative z
323 { 1.23F, 0.91F, 0.68F },
324 { 0.19F, 1.45F, 0.94F },
325 { 1.28F, 0.46F, 0.38F },
326 { 1.21F, 0.23F, 1.0F } };
327 //! 1 coordinate vector
328 CoordinatesVector const c_sampleCoordinates1(c_sampleCoordinatesFull.begin(),
329 c_sampleCoordinatesFull.begin() + 1);
330 //! 2 coordinate vectors
331 CoordinatesVector const c_sampleCoordinates2(c_sampleCoordinatesFull.begin() + 1,
332 c_sampleCoordinatesFull.begin() + 3);
333 //! 13 coordinate vectors
334 CoordinatesVector const c_sampleCoordinates13(c_sampleCoordinatesFull.begin() + 3,
335 c_sampleCoordinatesFull.begin() + 16);
337 //! moved out from instantiantions for readability
338 auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
339 //! moved out from instantiantions for readability
340 auto c_inputPmeOrders = ::testing::Range(3, 5 + 1);
341 //! moved out from instantiantions for readability
342 auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
344 /*! \brief Instantiation of the test with valid input and 1 atom */
345 INSTANTIATE_TEST_CASE_P(SaneInput1,
346 PmeSplineAndSpreadTest,
347 ::testing::Combine(c_inputBoxes,
350 ::testing::Values(c_sampleCoordinates1),
351 ::testing::Values(c_sampleCharges1)));
353 /*! \brief Instantiation of the test with valid input and 2 atoms */
354 INSTANTIATE_TEST_CASE_P(SaneInput2,
355 PmeSplineAndSpreadTest,
356 ::testing::Combine(c_inputBoxes,
359 ::testing::Values(c_sampleCoordinates2),
360 ::testing::Values(c_sampleCharges2)));
361 /*! \brief Instantiation of the test with valid input and 13 atoms */
362 INSTANTIATE_TEST_CASE_P(SaneInput13,
363 PmeSplineAndSpreadTest,
364 ::testing::Combine(c_inputBoxes,
367 ::testing::Values(c_sampleCoordinates13),
368 ::testing::Values(c_sampleCharges13)));