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 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 //! Default constructor
76 PmeSolveTest() = default;
81 /* Getting the input */
84 SparseComplexGridValuesInput nonZeroGridValues;
88 PmeSolveAlgorithm method;
89 std::tie(box, gridSize, nonZeroGridValues, epsilon_r, ewaldCoeff_q, ewaldCoeff_lj, method) = GetParam();
91 /* Storing the input where it's needed, running the test */
93 inputRec.nkx = gridSize[XX];
94 inputRec.nky = gridSize[YY];
95 inputRec.nkz = gridSize[ZZ];
96 inputRec.pme_order = 4;
97 inputRec.coulombtype = eelPME;
98 inputRec.epsilon_r = epsilon_r;
101 case PmeSolveAlgorithm::Coulomb:
104 case PmeSolveAlgorithm::LennardJones:
105 inputRec.vdwtype = evdwPME;
109 GMX_THROW(InternalError("Unknown PME solver"));
112 TestReferenceData refData;
113 const std::map<CodePath, std::string> modesToTest = {{CodePath::CPU, "CPU"},
114 {CodePath::CUDA, "CUDA"}};
115 for (const auto &mode : modesToTest)
117 const bool supportedInput = pmeSupportsInputForMode(&inputRec, mode.first);
120 /* Testing the failure for the unsupported input */
121 EXPECT_THROW(pmeInitEmpty(&inputRec, mode.first, nullptr, box, ewaldCoeff_q, ewaldCoeff_lj), NotImplementedError);
125 std::map<GridOrdering, std::string> gridOrderingsToTest = {{GridOrdering::YZX, "YZX"}};
126 if (mode.first == CodePath::CUDA)
128 gridOrderingsToTest[GridOrdering::XYZ] = "XYZ";
130 const auto contextsToTest = getPmeTestEnv()->getHardwareContexts(mode.first);
131 for (const auto &gridOrdering : gridOrderingsToTest)
133 for (const auto &context : contextsToTest)
135 for (bool computeEnergyAndVirial : {false, true})
137 /* Describing the test*/
138 SCOPED_TRACE(formatString("Testing solving (%s, %s, %s energy/virial) with %s %sfor PME grid size %d %d %d, Ewald coefficients %g %g",
139 (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
140 gridOrdering.second.c_str(),
141 computeEnergyAndVirial ? "with" : "without",
143 context.getDescription().c_str(),
144 gridSize[XX], gridSize[YY], gridSize[ZZ],
145 ewaldCoeff_q, ewaldCoeff_lj
148 /* Running the test */
149 PmeSafePointer pmeSafe = pmeInitEmpty(&inputRec, mode.first, context.getDeviceInfo(), box, ewaldCoeff_q, ewaldCoeff_lj);
150 pmeSetComplexGrid(pmeSafe.get(), mode.first, gridOrdering.first, nonZeroGridValues);
151 const real cellVolume = box[0] * box[4] * box[8];
152 //FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
153 pmePerformSolve(pmeSafe.get(), mode.first, method, cellVolume, gridOrdering.first, computeEnergyAndVirial);
154 pmeFinalizeTest(pmeSafe.get(), mode.first);
156 /* Check the outputs */
157 TestReferenceChecker checker(refData.rootChecker());
159 SparseComplexGridValuesOutput nonZeroGridValuesOutput = pmeGetComplexGrid(pmeSafe.get(), mode.first, gridOrdering.first);
160 /* Transformed grid */
161 TestReferenceChecker gridValuesChecker(checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
163 real gridValuesMagnitude = 1.0;
164 for (const auto &point : nonZeroGridValuesOutput)
166 gridValuesMagnitude = std::max(std::fabs(point.second.re), gridValuesMagnitude);
167 gridValuesMagnitude = std::max(std::fabs(point.second.im), gridValuesMagnitude);
169 // Spline moduli participate 3 times in the computation; 2 is an additional factor for SIMD exp() precision
170 gmx_uint64_t gridUlpToleranceFactor = DIM * 2;
171 if (method == PmeSolveAlgorithm::LennardJones)
173 // Lennard Jones is more complex and also uses erfc(), relax more
174 gridUlpToleranceFactor *= 2;
176 const gmx_uint64_t splineModuliDoublePrecisionUlps
177 = getSplineModuliDoublePrecisionUlps(inputRec.pme_order + 1);
179 = relativeToleranceAsPrecisionDependentUlp(gridValuesMagnitude,
180 gridUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
181 gridUlpToleranceFactor * splineModuliDoublePrecisionUlps);
182 gridValuesChecker.setDefaultTolerance(gridTolerance);
184 for (const auto &point : nonZeroGridValuesOutput)
186 // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
187 // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
188 if (fabs(point.second.re) >= GMX_FLOAT_MIN)
190 gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
192 if (fabs(point.second.im) >= GMX_FLOAT_MIN)
194 gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
198 if (computeEnergyAndVirial)
200 // Extract the energy and virial
203 std::tie(energy, virial) = pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), mode.first, method);
205 // These quantities are computed based on the grid values, so must have
206 // checking relative tolerances at least as large. Virial needs more flops
207 // than energy, so needs a larger tolerance.
210 double energyMagnitude = 10.0;
211 // TODO This factor is arbitrary, do a proper error-propagation analysis
212 gmx_uint64_t energyUlpToleranceFactor = gridUlpToleranceFactor * 2;
214 = relativeToleranceAsPrecisionDependentUlp(energyMagnitude,
215 energyUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
216 energyUlpToleranceFactor * splineModuliDoublePrecisionUlps);
217 TestReferenceChecker energyChecker(checker);
218 energyChecker.setDefaultTolerance(energyTolerance);
219 energyChecker.checkReal(energy, "Energy");
222 double virialMagnitude = 1000.0;
223 // TODO This factor is arbitrary, do a proper error-propagation analysis
224 gmx_uint64_t virialUlpToleranceFactor = energyUlpToleranceFactor * 2;
226 = relativeToleranceAsPrecisionDependentUlp(virialMagnitude,
227 virialUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
228 virialUlpToleranceFactor * splineModuliDoublePrecisionUlps);
229 TestReferenceChecker virialChecker(checker.checkCompound("Matrix", "Virial"));
230 virialChecker.setDefaultTolerance(virialTolerance);
231 for (int i = 0; i < DIM; i++)
233 for (int j = 0; j <= i; j++)
235 std::string valueId = formatString("Cell %d %d", i, j);
236 virialChecker.checkReal(virial[i * DIM + j], valueId.c_str());
248 /*! \brief Test for PME solving */
249 TEST_P(PmeSolveTest, ReproducesOutputs)
251 EXPECT_NO_THROW(runTest());
254 /* Valid input instances */
256 //! A couple of valid inputs for boxes.
257 static std::vector<Matrix3x3> const c_sampleBoxes
273 //! A couple of valid inputs for grid sizes
274 static std::vector<IVec> const c_sampleGridSizes
284 //! Moved out from instantiations for readability
285 const auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
286 //! Moved out from instantiations for readability
287 const auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
289 //! 2 sample complex grids - only non-zero values have to be listed
290 static std::vector<SparseComplexGridValuesInput> const c_sampleGrids
292 SparseComplexGridValuesInput {{
323 SparseComplexGridValuesInput {{
350 //! Moved out from instantiations for readability
351 const auto c_inputGrids = ::testing::ValuesIn(c_sampleGrids);
352 //! Moved out from instantiations for readability
353 const auto c_inputEpsilon_r = ::testing::Values(1.2);
354 //! Moved out from instantiations for readability
355 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
356 //! Moved out from instantiations for readability
357 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
358 //! Moved out from instantiations for readability
359 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
361 //! Instantiation of the PME solving test
362 INSTANTIATE_TEST_CASE_P(SaneInput, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
363 c_inputEpsilon_r, c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj, c_inputMethods));
365 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
366 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffQ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
367 c_inputEpsilon_r, ::testing::Values(0.4), c_inputEwaldCoeff_lj,
368 ::testing::Values(PmeSolveAlgorithm::Coulomb)));
370 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
371 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
372 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
373 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffLJ, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
374 c_inputEpsilon_r, c_inputEwaldCoeff_q, ::testing::Values(2.35),
375 ::testing::Values(PmeSolveAlgorithm::LennardJones)));
377 //! A few more instances to check that different epsilon_r actually affects results of all solvers
378 INSTANTIATE_TEST_CASE_P(DifferentEpsilonR, PmeSolveTest, ::testing::Combine(c_inputBoxes, c_inputGridSizes, c_inputGrids,
379 testing::Values(1.9), c_inputEwaldCoeff_q, c_inputEwaldCoeff_lj,