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/mdrunutility/mdmodules.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/utility/stringutil.h"
53 #include "testutils/refdata.h"
54 #include "testutils/testasserts.h"
56 #include "pmetestcommon.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 dimensions,
74 * 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 //! Environment for getting the t_inputrec structure easily
89 //! Default constructor
90 PmeSplineAndSpreadTest() = default;
94 /* Getting the input */
98 CoordinatesVector coordinates;
99 ChargesVector charges;
101 std::tie(box, pmeOrder, gridSize, coordinates, charges) = GetParam();
102 const size_t atomCount = coordinates.size();
104 /* Storing the input where it's needed */
105 t_inputrec *inputRec = mdModules_.inputrec();
106 inputRec->nkx = gridSize[XX];
107 inputRec->nky = gridSize[YY];
108 inputRec->nkz = gridSize[ZZ];
109 inputRec->pme_order = pmeOrder;
110 inputRec->coulombtype = eelPME;
112 TestReferenceData refData;
114 const std::map<CodePath, std::string> modesToTest = {{CodePath::CPU, "CPU"}};
115 const std::map<PmeSplineAndSpreadOptions, std::string> optionsToTest = {{PmeSplineAndSpreadOptions::SplineAndSpreadUnified, "spline computation and charge spreading (fused)"},
116 {PmeSplineAndSpreadOptions::SplineOnly, "spline computation"},
117 {PmeSplineAndSpreadOptions::SpreadOnly, "charge spreading"}};
118 for (const auto &mode : modesToTest)
120 for (const auto &option : optionsToTest)
122 /* Describing the test uniquely in case it fails */
124 SCOPED_TRACE(formatString("Testing %s with %s for PME grid size %d %d %d"
125 ", order %d, %zu atoms",
126 option.second.c_str(), mode.second.c_str(),
127 gridSize[XX], gridSize[YY], gridSize[ZZ],
131 /* Running the test */
133 PmeSafePointer pmeSafe = pmeInitWithAtoms(inputRec, coordinates, charges, box);
135 const bool computeSplines = (option.first == PmeSplineAndSpreadOptions::SplineOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
136 const bool spreadCharges = (option.first == PmeSplineAndSpreadOptions::SpreadOnly) || (option.first == PmeSplineAndSpreadOptions::SplineAndSpreadUnified);
140 // Here we should set up the results of the spline computation so that the spread can run.
141 // What is lazy and works is running the separate spline so that it will set it up for us:
142 pmePerformSplineAndSpread(pmeSafe.get(), mode.first, true, false);
143 // We know that it is tested in another iteration.
144 // TODO: Clean alternative: read and set the reference gridline indices, spline params
147 pmePerformSplineAndSpread(pmeSafe.get(), mode.first, computeSplines, spreadCharges);
149 /* Outputs correctness check */
150 /* All tolerances were picked empirically for single precision on CPU */
152 TestReferenceChecker rootChecker(refData.rootChecker());
154 const auto maxGridSize = std::max(std::max(gridSize[XX], gridSize[YY]), gridSize[ZZ]);
155 const auto ulpToleranceSplineValues = 2 * (pmeOrder - 2) * maxGridSize;
156 /* 2 is empiric, the rest follows from the amount of operations */
160 const char *dimString[] = { "X", "Y", "Z" };
163 SCOPED_TRACE(formatString("Testing spline values with tolerance of %ld", ulpToleranceSplineValues));
164 TestReferenceChecker splineValuesChecker(rootChecker.checkCompound("Splines", "Values"));
165 splineValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineValues));
166 for (int i = 0; i < DIM; i++)
168 auto splineValuesDim = pmeGetSplineData(pmeSafe.get(), mode.first, PmeSplineDataType::Values, i);
169 splineValuesChecker.checkSequence(splineValuesDim.begin(), splineValuesDim.end(), dimString[i]);
172 /* Spline derivatives */
173 const auto ulpToleranceSplineDerivatives = 4 * ulpToleranceSplineValues;
174 /* 4 is just a wild guess since the derivatives are deltas of neighbor spline values which could differ greatly */
175 SCOPED_TRACE(formatString("Testing spline derivatives with tolerance of %ld", ulpToleranceSplineDerivatives));
176 TestReferenceChecker splineDerivativesChecker(rootChecker.checkCompound("Splines", "Derivatives"));
177 splineDerivativesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceSplineDerivatives));
178 for (int i = 0; i < DIM; i++)
180 auto splineDerivativesDim = pmeGetSplineData(pmeSafe.get(), mode.first, PmeSplineDataType::Derivatives, i);
181 splineDerivativesChecker.checkSequence(splineDerivativesDim.begin(), splineDerivativesDim.end(), dimString[i]);
184 /* Particle gridline indices */
185 auto gridLineIndices = pmeGetGridlineIndices(pmeSafe.get(), mode.first);
186 rootChecker.checkSequence(gridLineIndices.begin(), gridLineIndices.end(), "Gridline indices");
191 /* The wrapped grid */
192 SparseRealGridValuesOutput nonZeroGridValues = pmeGetRealGrid(pmeSafe.get(), mode.first);
193 TestReferenceChecker gridValuesChecker(rootChecker.checkCompound("NonZeroGridValues", "RealSpaceGrid"));
194 const auto ulpToleranceGrid = 2 * ulpToleranceSplineValues * (int)(ceil(sqrt(atomCount)));
195 /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */
196 SCOPED_TRACE(formatString("Testing grid values with tolerance of %ld", ulpToleranceGrid));
197 gridValuesChecker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpToleranceGrid));
198 for (const auto &point : nonZeroGridValues)
200 gridValuesChecker.checkReal(point.second, point.first.c_str());
209 /*! \brief Test for PME B-spline moduli computation */
210 TEST_P(PmeSplineAndSpreadTest, ReproducesOutputs)
212 EXPECT_NO_THROW(runTest());
215 /* Valid input instances */
217 //! A couple of valid inputs for boxes.
218 static std::vector<Matrix3x3> const sampleBoxes
234 //! A couple of valid inputs for grid sizes.
235 static std::vector<IVec> const sampleGridSizes
246 static std::vector<real> const sampleChargesFull
248 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
251 static auto const sampleCharges1 = ChargesVector::fromVector(sampleChargesFull.begin(), sampleChargesFull.begin() + 1);
253 static auto const sampleCharges2 = ChargesVector::fromVector(sampleChargesFull.begin() + 1, sampleChargesFull.begin() + 3);
255 static auto const sampleCharges13 = ChargesVector::fromVector(sampleChargesFull.begin() + 3, sampleChargesFull.begin() + 16);
257 //! Random coordinate vectors
258 static CoordinatesVector const sampleCoordinatesFull
263 16.0f, 1.02f, 0.22f // 2 box lengths in x
273 16.0001f, 1.52f, 1.19f // > 2 box lengths in x
275 1.43f, 1.1f, 4.1f // > 2 box lengths in z
277 -1.08f, 1.19f, 0.08f // negative x
281 1.32f, -1.48f, 0.16f // negative y
285 0.95f, 7.7f, -0.48f // > 2 box lengths in y, negative z
296 //! 1 coordinate vector
297 static CoordinatesVector const sampleCoordinates1(sampleCoordinatesFull.begin(), sampleCoordinatesFull.begin() + 1);
298 //! 2 coordinate vectors
299 static CoordinatesVector const sampleCoordinates2(sampleCoordinatesFull.begin() + 1, sampleCoordinatesFull.begin() + 3);
300 //! 13 coordinate vectors
301 static CoordinatesVector const sampleCoordinates13(sampleCoordinatesFull.begin() + 3, sampleCoordinatesFull.begin() + 16);
303 //! moved out from instantiantions for readability
304 auto inputBoxes = ::testing::ValuesIn(sampleBoxes);
305 //! moved out from instantiantions for readability
306 auto inputPmeOrders = ::testing::Range(3, 5 + 1);
307 //! moved out from instantiantions for readability
308 auto inputGridSizes = ::testing::ValuesIn(sampleGridSizes);
310 /*! \brief Instantiation of the PME spline computation test with valid input and 1 atom */
311 INSTANTIATE_TEST_CASE_P(SaneInput1, PmeSplineAndSpreadTest, ::testing::Combine(inputBoxes, inputPmeOrders, inputGridSizes,
312 ::testing::Values(sampleCoordinates1),
313 ::testing::Values(sampleCharges1)
315 /*! \brief Instantiation of the PME spline computation test with valid input and 2 atoms */
316 INSTANTIATE_TEST_CASE_P(SaneInput2, PmeSplineAndSpreadTest, ::testing::Combine(inputBoxes, inputPmeOrders, inputGridSizes,
317 ::testing::Values(sampleCoordinates2),
318 ::testing::Values(sampleCharges2)
320 /*! \brief Instantiation of the PME spline computation test with valid input and 13 atoms */
321 INSTANTIATE_TEST_CASE_P(SaneInput13, PmeSplineAndSpreadTest, ::testing::Combine(inputBoxes, inputPmeOrders, inputGridSizes,
322 ::testing::Values(sampleCoordinates13),
323 ::testing::Values(sampleCharges13)