2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019, 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 solving 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"
63 /*! \brief Convenience typedef of the test input parameters - unit cell box, complex grid dimensions, complex grid values,
64 * electrostatic constant epsilon_r, Ewald splitting parameters ewaldcoeff_q and ewaldcoeff_lj, solver type
65 * Output: transformed local grid (Fourier space); optionally reciprocal energy and virial matrix.
67 * Implement and test Lorentz-Berthelot
69 typedef std::tuple<Matrix3x3, IVec, SparseComplexGridValuesInput, double, double, double, PmeSolveAlgorithm> SolveInputParameters;
72 class PmeSolveTest : public ::testing::TestWithParam<SolveInputParameters>
75 PmeSolveTest() = default;
80 /* Getting the input */
83 SparseComplexGridValuesInput nonZeroGridValues;
87 PmeSolveAlgorithm method;
88 std::tie(box, gridSize, nonZeroGridValues, epsilon_r, ewaldCoeff_q, ewaldCoeff_lj, method) = GetParam();
90 /* Storing the input where it's needed, running the test */
92 inputRec.nkx = gridSize[XX];
93 inputRec.nky = gridSize[YY];
94 inputRec.nkz = gridSize[ZZ];
95 inputRec.pme_order = 4;
96 inputRec.coulombtype = eelPME;
97 inputRec.epsilon_r = epsilon_r;
100 case PmeSolveAlgorithm::Coulomb:
103 case PmeSolveAlgorithm::LennardJones:
104 inputRec.vdwtype = evdwPME;
108 GMX_THROW(InternalError("Unknown PME solver"));
111 TestReferenceData refData;
112 for (const auto &context : getPmeTestEnv()->getHardwareContexts())
114 CodePath codePath = context->getCodePath();
115 const bool supportedInput = pmeSupportsInputForMode(*getPmeTestEnv()->hwinfo(), &inputRec, codePath);
118 /* Testing the failure for the unsupported input */
119 EXPECT_THROW(pmeInitEmpty(&inputRec, codePath, nullptr, nullptr, box, ewaldCoeff_q, ewaldCoeff_lj), NotImplementedError);
123 std::map<GridOrdering, std::string> gridOrderingsToTest = {{GridOrdering::YZX, "YZX"}};
124 if (codePath == CodePath::GPU)
126 gridOrderingsToTest[GridOrdering::XYZ] = "XYZ";
128 for (const auto &gridOrdering : gridOrderingsToTest)
130 for (bool computeEnergyAndVirial : {false, true})
132 /* Describing the test*/
133 SCOPED_TRACE(formatString("Testing solving (%s, %s, %s energy/virial) with %s %sfor PME grid size %d %d %d, Ewald coefficients %g %g",
134 (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
135 gridOrdering.second.c_str(),
136 computeEnergyAndVirial ? "with" : "without",
137 codePathToString(codePath),
138 context->getDescription().c_str(),
139 gridSize[XX], gridSize[YY], gridSize[ZZ],
140 ewaldCoeff_q, ewaldCoeff_lj
143 /* Running the test */
144 PmeSafePointer pmeSafe = pmeInitEmpty(&inputRec, codePath, context->getDeviceInfo(),
145 context->getPmeGpuProgram(), box, ewaldCoeff_q, ewaldCoeff_lj);
146 pmeSetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first, nonZeroGridValues);
147 const real cellVolume = box[0] * box[4] * box[8];
148 //FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
149 pmePerformSolve(pmeSafe.get(), codePath, method, cellVolume, gridOrdering.first, computeEnergyAndVirial);
150 pmeFinalizeTest(pmeSafe.get(), codePath);
152 /* Check the outputs */
153 TestReferenceChecker checker(refData.rootChecker());
155 SparseComplexGridValuesOutput nonZeroGridValuesOutput = pmeGetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first);
156 /* Transformed grid */
157 TestReferenceChecker gridValuesChecker(checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
159 real gridValuesMagnitude = 1.0;
160 for (const auto &point : nonZeroGridValuesOutput)
162 gridValuesMagnitude = std::max(std::fabs(point.second.re), gridValuesMagnitude);
163 gridValuesMagnitude = std::max(std::fabs(point.second.im), gridValuesMagnitude);
165 // Spline moduli participate 3 times in the computation; 2 is an additional factor for SIMD exp() precision
166 uint64_t gridUlpToleranceFactor = DIM * 2;
167 if (method == PmeSolveAlgorithm::LennardJones)
169 // Lennard Jones is more complex and also uses erfc(), relax more
170 gridUlpToleranceFactor *= 2;
172 const uint64_t splineModuliDoublePrecisionUlps
173 = getSplineModuliDoublePrecisionUlps(inputRec.pme_order + 1);
175 = relativeToleranceAsPrecisionDependentUlp(gridValuesMagnitude,
176 gridUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
177 gridUlpToleranceFactor * splineModuliDoublePrecisionUlps);
178 gridValuesChecker.setDefaultTolerance(gridTolerance);
180 for (const auto &point : nonZeroGridValuesOutput)
182 // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
183 // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
184 if (fabs(point.second.re) >= GMX_FLOAT_MIN)
186 gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
188 if (fabs(point.second.im) >= GMX_FLOAT_MIN)
190 gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
194 if (computeEnergyAndVirial)
196 // Extract the energy and virial
197 const auto output = pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), codePath, method);
198 const auto &energy = (method == PmeSolveAlgorithm::Coulomb) ? output.coulombEnergy_ : output.lennardJonesEnergy_;
199 const auto &virial = (method == PmeSolveAlgorithm::Coulomb) ? output.coulombVirial_ : output.lennardJonesVirial_;
201 // These quantities are computed based on the grid values, so must have
202 // checking relative tolerances at least as large. Virial needs more flops
203 // than energy, so needs a larger tolerance.
206 double energyMagnitude = 10.0;
207 // TODO This factor is arbitrary, do a proper error-propagation analysis
208 uint64_t energyUlpToleranceFactor = gridUlpToleranceFactor * 2;
210 = relativeToleranceAsPrecisionDependentUlp(energyMagnitude,
211 energyUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
212 energyUlpToleranceFactor * splineModuliDoublePrecisionUlps);
213 TestReferenceChecker energyChecker(checker);
214 energyChecker.setDefaultTolerance(energyTolerance);
215 energyChecker.checkReal(energy, "Energy");
218 double virialMagnitude = 1000.0;
219 // TODO This factor is arbitrary, do a proper error-propagation analysis
220 uint64_t virialUlpToleranceFactor = energyUlpToleranceFactor * 2;
222 = relativeToleranceAsPrecisionDependentUlp(virialMagnitude,
223 virialUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
224 virialUlpToleranceFactor * splineModuliDoublePrecisionUlps);
225 TestReferenceChecker virialChecker(checker.checkCompound("Matrix", "Virial"));
226 virialChecker.setDefaultTolerance(virialTolerance);
227 for (int i = 0; i < DIM; i++)
229 for (int j = 0; j <= i; j++)
231 std::string valueId = formatString("Cell %d %d", i, j);
232 virialChecker.checkReal(virial[i][j], valueId.c_str());
243 /*! \brief Test for PME solving */
244 TEST_P(PmeSolveTest, ReproducesOutputs)
246 EXPECT_NO_THROW(runTest());
249 /* Valid input instances */
251 //! A couple of valid inputs for boxes.
252 std::vector<Matrix3x3> const c_sampleBoxes
268 //! A couple of valid inputs for grid sizes
269 std::vector<IVec> const c_sampleGridSizes
279 //! Moved out from instantiations for readability
280 const auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
281 //! Moved out from instantiations for readability
282 const auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
284 //! 2 sample complex grids - only non-zero values have to be listed
285 std::vector<SparseComplexGridValuesInput> const c_sampleGrids
287 SparseComplexGridValuesInput {{
318 SparseComplexGridValuesInput {{
345 //! Moved out from instantiations for readability
346 const auto c_inputGrids = ::testing::ValuesIn(c_sampleGrids);
347 //! Moved out from instantiations for readability
348 const auto c_inputEpsilon_r = ::testing::Values(1.2);
349 //! Moved out from instantiations for readability
350 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
351 //! Moved out from instantiations for readability
352 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
353 //! Moved out from instantiations for readability
354 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
356 //! Instantiation of the PME solving test
357 INSTANTIATE_TEST_CASE_P(SaneInput, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
358 c_inputEpsilon_r, c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj, c_inputMethods));
360 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
361 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffQ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
362 c_inputEpsilon_r, ::testing::Values(0.4), c_inputEwaldCoeff_lj,
363 ::testing::Values(PmeSolveAlgorithm::Coulomb)));
365 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
366 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
367 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
368 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffLJ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
369 c_inputEpsilon_r, c_inputEwaldCoeff_q, ::testing::Values(2.35),
370 ::testing::Values(PmeSolveAlgorithm::LennardJones)));
372 //! A few more instances to check that different epsilon_r actually affects results of all solvers
373 INSTANTIATE_TEST_CASE_P(DifferentEpsilonR, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
374 testing::Values(1.9), c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj,