2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019,2020,2021, 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 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \ingroup module_ewald
48 #include <gmock/gmock.h>
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/utility/stringutil.h"
53 #include "testutils/refdata.h"
54 #include "testutils/test_hardware_environment.h"
55 #include "testutils/testasserts.h"
57 #include "pmetestcommon.h"
66 //! A couple of valid inputs for grid sizes
67 std::vector<IVec> const c_inputGridSizes{ IVec{ 16, 12, 28 }, IVec{ 9, 7, 23 } };
69 //! Two input complex grids - only non-zero values have to be listed
70 const std::map<std::string, SparseComplexGridValuesInput> c_inputGridValues = {
72 SparseComplexGridValuesInput{
73 { IVec{ 0, 0, 0 }, t_complex{ 3.5F, 6.7F } },
74 { IVec{ 7, 0, 0 }, t_complex{ -2.5F, -0.7F } },
75 { IVec{ 3, 5, 7 }, t_complex{ -0.006F, 1e-8F } },
76 { IVec{ 3, 1, 2 }, t_complex{ 0.6F, 7.9F } },
77 { IVec{ 6, 2, 4 }, t_complex{ 30.1F, 2.45F } },
80 SparseComplexGridValuesInput{
81 { IVec{ 0, 4, 0 }, t_complex{ 0.0F, 0.3F } },
82 { IVec{ 4, 2, 7 }, t_complex{ 13.76F, -40.0F } },
83 { IVec{ 0, 6, 7 }, t_complex{ 3.6F, 0.0F } },
84 { IVec{ 2, 5, 10 }, t_complex{ 3.6F, 10.65F } },
88 /*! \brief Convenience typedef of the test input parameters - unit cell box, complex grid dimensions, complex grid values,
89 * electrostatic constant epsilon_r, Ewald splitting parameters ewaldcoeff_q and ewaldcoeff_lj, solver type
90 * Output: transformed local grid (Fourier space); optionally reciprocal energy and virial matrix.
92 * Implement and test Lorentz-Berthelot
94 typedef std::tuple<std::string, IVec, std::string, double, double, double, PmeSolveAlgorithm> SolveInputParameters;
96 const char* enumValueToString(PmeSolveAlgorithm enumValue)
98 static constexpr gmx::EnumerationArray<PmeSolveAlgorithm, const char*> s_pmeSolveAlgorithmNames = {
101 return s_pmeSolveAlgorithmNames[enumValue];
104 //! Help GoogleTest name our test cases
105 std::string nameOfTest(const testing::TestParamInfo<SolveInputParameters>& info)
107 std::string testName = formatString(
115 std::get<0>(info.param).c_str(),
116 std::get<1>(info.param)[XX],
117 std::get<1>(info.param)[YY],
118 std::get<1>(info.param)[ZZ],
119 std::get<2>(info.param).c_str(),
120 std::get<3>(info.param),
121 std::get<4>(info.param),
122 std::get<5>(info.param),
123 enumValueToString(std::get<6>(info.param)));
125 // Note that the returned names must be unique and may use only
126 // alphanumeric ASCII characters. It's not supposed to contain
127 // underscores (see the GoogleTest FAQ
128 // why-should-test-suite-names-and-test-names-not-contain-underscore),
129 // but doing so works for now, is likely to remain so, and makes
130 // such test names much more readable.
131 testName = replaceAll(testName, "-", "_");
132 testName = replaceAll(testName, ".", "_");
133 testName = replaceAll(testName, " ", "_");
138 class SolveTest : public ::testing::TestWithParam<SolveInputParameters>
141 SolveTest() = default;
143 //! Sets the programs once
144 static void SetUpTestSuite()
146 s_pmeTestHardwareContexts = createPmeTestHardwareContextList();
147 g_allowPmeWithSyclForTesting = true; // We support Solve with SYCL
150 static void TearDownTestSuite()
152 // Revert the value back.
153 g_allowPmeWithSyclForTesting = false;
157 static void runTest()
159 /* Getting the input */
163 double ewaldCoeff_lj;
164 PmeSolveAlgorithm method;
165 std::string boxName, gridValuesName;
166 std::tie(boxName, gridSize, gridValuesName, epsilon_r, ewaldCoeff_q, ewaldCoeff_lj, method) =
168 Matrix3x3 box = c_inputBoxes.at(boxName);
169 const SparseComplexGridValuesInput& nonZeroGridValues = c_inputGridValues.at(gridValuesName);
171 /* Storing the input where it's needed, running the test */
173 inputRec.nkx = gridSize[XX];
174 inputRec.nky = gridSize[YY];
175 inputRec.nkz = gridSize[ZZ];
176 inputRec.pme_order = 4;
177 inputRec.coulombtype = CoulombInteractionType::Pme;
178 inputRec.epsilon_r = epsilon_r;
181 case PmeSolveAlgorithm::Coulomb: break;
183 case PmeSolveAlgorithm::LennardJones: inputRec.vdwtype = VanDerWaalsType::Pme; break;
185 default: GMX_THROW(InternalError("Unknown PME solver"));
188 TestReferenceData refData;
189 for (const auto& pmeTestHardwareContext : s_pmeTestHardwareContexts)
191 pmeTestHardwareContext->activate();
192 CodePath codePath = pmeTestHardwareContext->codePath();
193 const bool supportedInput = pmeSupportsInputForMode(
194 *getTestHardwareEnvironment()->hwinfo(), &inputRec, codePath);
197 /* Testing the failure for the unsupported input */
199 pmeInitWrapper(&inputRec, codePath, nullptr, nullptr, nullptr, box, ewaldCoeff_q, ewaldCoeff_lj),
200 NotImplementedError);
204 std::map<GridOrdering, std::string> gridOrderingsToTest = { { GridOrdering::YZX,
206 if (codePath == CodePath::GPU)
208 gridOrderingsToTest[GridOrdering::XYZ] = "XYZ";
210 for (const auto& gridOrdering : gridOrderingsToTest)
212 for (bool computeEnergyAndVirial : { false, true })
214 /* Describing the test*/
215 SCOPED_TRACE(formatString(
216 "Testing solving (%s, %s, %s energy/virial) on %s for PME grid "
217 "size %d %d %d, Ewald coefficients %g %g",
218 (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
219 gridOrdering.second.c_str(),
220 computeEnergyAndVirial ? "with" : "without",
221 pmeTestHardwareContext->description().c_str(),
228 /* Running the test */
229 PmeSafePointer pmeSafe = pmeInitWrapper(&inputRec,
231 pmeTestHardwareContext->deviceContext(),
232 pmeTestHardwareContext->deviceStream(),
233 pmeTestHardwareContext->pmeGpuProgram(),
237 pmeSetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first, nonZeroGridValues);
238 const real cellVolume = box[0] * box[4] * box[8];
239 // FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
240 pmePerformSolve(pmeSafe.get(), codePath, method, cellVolume, gridOrdering.first, computeEnergyAndVirial);
241 pmeFinalizeTest(pmeSafe.get(), codePath);
243 /* Check the outputs */
244 TestReferenceChecker checker(refData.rootChecker());
246 SparseComplexGridValuesOutput nonZeroGridValuesOutput =
247 pmeGetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first);
248 /* Transformed grid */
249 TestReferenceChecker gridValuesChecker(
250 checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
252 real gridValuesMagnitude = 1.0;
253 for (const auto& point : nonZeroGridValuesOutput)
255 gridValuesMagnitude = std::max(std::fabs(point.second.re), gridValuesMagnitude);
256 gridValuesMagnitude = std::max(std::fabs(point.second.im), gridValuesMagnitude);
258 // Spline moduli participate 3 times in the computation; 2 is an additional factor for SIMD exp() precision
259 uint64_t gridUlpToleranceFactor = DIM * 2;
260 if (method == PmeSolveAlgorithm::LennardJones)
262 // Lennard Jones is more complex and also uses erfc(), relax more
263 gridUlpToleranceFactor *= 2;
265 const uint64_t splineModuliDoublePrecisionUlps =
266 getSplineModuliDoublePrecisionUlps(inputRec.pme_order + 1);
267 auto gridTolerance = relativeToleranceAsPrecisionDependentUlp(
269 gridUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
270 gridUlpToleranceFactor * splineModuliDoublePrecisionUlps);
271 gridValuesChecker.setDefaultTolerance(gridTolerance);
273 for (const auto& point : nonZeroGridValuesOutput)
275 // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
276 // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
277 if (fabs(point.second.re) >= GMX_FLOAT_MIN)
279 gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
281 if (fabs(point.second.im) >= GMX_FLOAT_MIN)
283 gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
287 if (computeEnergyAndVirial)
289 // Extract the energy and virial
291 pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), codePath, method);
292 const auto& energy = (method == PmeSolveAlgorithm::Coulomb)
293 ? output.coulombEnergy_
294 : output.lennardJonesEnergy_;
295 const auto& virial = (method == PmeSolveAlgorithm::Coulomb)
296 ? output.coulombVirial_
297 : output.lennardJonesVirial_;
299 // These quantities are computed based on the grid values, so must have
300 // checking relative tolerances at least as large. Virial needs more flops
301 // than energy, so needs a larger tolerance.
304 double energyMagnitude = 10.0;
305 // TODO This factor is arbitrary, do a proper error-propagation analysis
306 uint64_t energyUlpToleranceFactor = gridUlpToleranceFactor * 2;
307 auto energyTolerance = relativeToleranceAsPrecisionDependentUlp(
309 energyUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
310 energyUlpToleranceFactor * splineModuliDoublePrecisionUlps);
311 TestReferenceChecker energyChecker(checker);
312 energyChecker.setDefaultTolerance(energyTolerance);
313 energyChecker.checkReal(energy, "Energy");
316 double virialMagnitude = 1000.0;
317 // TODO This factor is arbitrary, do a proper error-propagation analysis
318 uint64_t virialUlpToleranceFactor = energyUlpToleranceFactor * 2;
319 auto virialTolerance = relativeToleranceAsPrecisionDependentUlp(
321 virialUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
322 virialUlpToleranceFactor * splineModuliDoublePrecisionUlps);
323 TestReferenceChecker virialChecker(
324 checker.checkCompound("Matrix", "Virial"));
325 virialChecker.setDefaultTolerance(virialTolerance);
326 for (int i = 0; i < DIM; i++)
328 for (int j = 0; j <= i; j++)
330 std::string valueId = formatString("Cell %d %d", i, j);
331 virialChecker.checkReal(virial[i][j], valueId.c_str());
340 static std::vector<std::unique_ptr<PmeTestHardwareContext>> s_pmeTestHardwareContexts;
343 std::vector<std::unique_ptr<PmeTestHardwareContext>> SolveTest::s_pmeTestHardwareContexts;
345 /*! \brief Test for PME solving */
346 TEST_P(SolveTest, WorksWith)
348 EXPECT_NO_THROW_GMX(runTest());
351 //! Moved out from instantiations for readability
352 const auto c_inputBoxNames = ::testing::Values("rect", "tric");
353 //! Moved out from instantiations for readability
354 const auto c_inputGridNames = ::testing::Values("first", "second");
355 //! Moved out from instantiations for readability
356 const auto c_inputEpsilon_r = ::testing::Values(1.2);
357 //! Moved out from instantiations for readability
358 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
359 //! Moved out from instantiations for readability
360 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
361 //! Moved out from instantiations for readability
362 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
364 //! Instantiation of the PME solving test
365 INSTANTIATE_TEST_SUITE_P(Pme,
367 ::testing::Combine(c_inputBoxNames,
368 ::testing::ValuesIn(c_inputGridSizes),
372 c_inputEwaldCoeff_lj,
376 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
377 INSTANTIATE_TEST_SUITE_P(PmeDiffEwaldQ,
379 ::testing::Combine(c_inputBoxNames,
380 ::testing::ValuesIn(c_inputGridSizes),
383 ::testing::Values(0.4),
384 c_inputEwaldCoeff_lj,
385 ::testing::Values(PmeSolveAlgorithm::Coulomb)),
388 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
389 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
390 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
391 INSTANTIATE_TEST_SUITE_P(PmeDiffEwaldLJ,
393 ::testing::Combine(c_inputBoxNames,
394 ::testing::ValuesIn(c_inputGridSizes),
398 ::testing::Values(2.35),
399 ::testing::Values(PmeSolveAlgorithm::LennardJones)),
402 //! A few more instances to check that different epsilon_r actually affects results of all solvers
403 INSTANTIATE_TEST_SUITE_P(PmeDiffEps,
405 ::testing::Combine(c_inputBoxNames,
406 ::testing::ValuesIn(c_inputGridSizes),
408 testing::Values(1.9),
410 c_inputEwaldCoeff_lj,