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 * \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/test_hardware_environment.h"
54 #include "testutils/testasserts.h"
56 #include "pmetestcommon.h"
64 /*! \brief Convenience typedef of the test input parameters - unit cell box, complex grid dimensions, complex grid values,
65 * electrostatic constant epsilon_r, Ewald splitting parameters ewaldcoeff_q and ewaldcoeff_lj, solver type
66 * Output: transformed local grid (Fourier space); optionally reciprocal energy and virial matrix.
68 * Implement and test Lorentz-Berthelot
70 typedef std::tuple<Matrix3x3, IVec, SparseComplexGridValuesInput, double, double, double, PmeSolveAlgorithm> SolveInputParameters;
73 class PmeSolveTest : public ::testing::TestWithParam<SolveInputParameters>
76 PmeSolveTest() = default;
78 //! Sets the programs once
79 static void SetUpTestSuite() { s_pmeTestHardwareContexts = createPmeTestHardwareContextList(); }
84 /* Getting the input */
87 SparseComplexGridValuesInput nonZeroGridValues;
91 PmeSolveAlgorithm method;
92 std::tie(box, gridSize, nonZeroGridValues, epsilon_r, ewaldCoeff_q, ewaldCoeff_lj, method) =
95 /* Storing the input where it's needed, running the test */
97 inputRec.nkx = gridSize[XX];
98 inputRec.nky = gridSize[YY];
99 inputRec.nkz = gridSize[ZZ];
100 inputRec.pme_order = 4;
101 inputRec.coulombtype = CoulombInteractionType::Pme;
102 inputRec.epsilon_r = epsilon_r;
105 case PmeSolveAlgorithm::Coulomb: break;
107 case PmeSolveAlgorithm::LennardJones: inputRec.vdwtype = VanDerWaalsType::Pme; break;
109 default: GMX_THROW(InternalError("Unknown PME solver"));
112 TestReferenceData refData;
113 for (const auto& pmeTestHardwareContext : s_pmeTestHardwareContexts)
115 pmeTestHardwareContext->activate();
116 CodePath codePath = pmeTestHardwareContext->codePath();
117 const bool supportedInput = pmeSupportsInputForMode(
118 *getTestHardwareEnvironment()->hwinfo(), &inputRec, codePath);
121 /* Testing the failure for the unsupported input */
123 pmeInitWrapper(&inputRec, codePath, nullptr, nullptr, nullptr, box, ewaldCoeff_q, ewaldCoeff_lj),
124 NotImplementedError);
128 std::map<GridOrdering, std::string> gridOrderingsToTest = { { GridOrdering::YZX,
130 if (codePath == CodePath::GPU)
132 gridOrderingsToTest[GridOrdering::XYZ] = "XYZ";
134 for (const auto& gridOrdering : gridOrderingsToTest)
136 for (bool computeEnergyAndVirial : { false, true })
138 /* Describing the test*/
139 SCOPED_TRACE(formatString(
140 "Testing solving (%s, %s, %s energy/virial) on %s for PME grid "
141 "size %d %d %d, Ewald coefficients %g %g",
142 (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
143 gridOrdering.second.c_str(),
144 computeEnergyAndVirial ? "with" : "without",
145 pmeTestHardwareContext->description().c_str(),
152 /* Running the test */
153 PmeSafePointer pmeSafe = pmeInitWrapper(&inputRec,
155 pmeTestHardwareContext->deviceContext(),
156 pmeTestHardwareContext->deviceStream(),
157 pmeTestHardwareContext->pmeGpuProgram(),
161 pmeSetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first, nonZeroGridValues);
162 const real cellVolume = box[0] * box[4] * box[8];
163 // FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
164 pmePerformSolve(pmeSafe.get(), codePath, method, cellVolume, gridOrdering.first, computeEnergyAndVirial);
165 pmeFinalizeTest(pmeSafe.get(), codePath);
167 /* Check the outputs */
168 TestReferenceChecker checker(refData.rootChecker());
170 SparseComplexGridValuesOutput nonZeroGridValuesOutput =
171 pmeGetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first);
172 /* Transformed grid */
173 TestReferenceChecker gridValuesChecker(
174 checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
176 real gridValuesMagnitude = 1.0;
177 for (const auto& point : nonZeroGridValuesOutput)
179 gridValuesMagnitude = std::max(std::fabs(point.second.re), gridValuesMagnitude);
180 gridValuesMagnitude = std::max(std::fabs(point.second.im), gridValuesMagnitude);
182 // Spline moduli participate 3 times in the computation; 2 is an additional factor for SIMD exp() precision
183 uint64_t gridUlpToleranceFactor = DIM * 2;
184 if (method == PmeSolveAlgorithm::LennardJones)
186 // Lennard Jones is more complex and also uses erfc(), relax more
187 gridUlpToleranceFactor *= 2;
189 const uint64_t splineModuliDoublePrecisionUlps =
190 getSplineModuliDoublePrecisionUlps(inputRec.pme_order + 1);
191 auto gridTolerance = relativeToleranceAsPrecisionDependentUlp(
193 gridUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
194 gridUlpToleranceFactor * splineModuliDoublePrecisionUlps);
195 gridValuesChecker.setDefaultTolerance(gridTolerance);
197 for (const auto& point : nonZeroGridValuesOutput)
199 // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
200 // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
201 if (fabs(point.second.re) >= GMX_FLOAT_MIN)
203 gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
205 if (fabs(point.second.im) >= GMX_FLOAT_MIN)
207 gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
211 if (computeEnergyAndVirial)
213 // Extract the energy and virial
215 pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), codePath, method);
216 const auto& energy = (method == PmeSolveAlgorithm::Coulomb)
217 ? output.coulombEnergy_
218 : output.lennardJonesEnergy_;
219 const auto& virial = (method == PmeSolveAlgorithm::Coulomb)
220 ? output.coulombVirial_
221 : output.lennardJonesVirial_;
223 // These quantities are computed based on the grid values, so must have
224 // checking relative tolerances at least as large. Virial needs more flops
225 // than energy, so needs a larger tolerance.
228 double energyMagnitude = 10.0;
229 // TODO This factor is arbitrary, do a proper error-propagation analysis
230 uint64_t energyUlpToleranceFactor = gridUlpToleranceFactor * 2;
231 auto energyTolerance = relativeToleranceAsPrecisionDependentUlp(
233 energyUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
234 energyUlpToleranceFactor * splineModuliDoublePrecisionUlps);
235 TestReferenceChecker energyChecker(checker);
236 energyChecker.setDefaultTolerance(energyTolerance);
237 energyChecker.checkReal(energy, "Energy");
240 double virialMagnitude = 1000.0;
241 // TODO This factor is arbitrary, do a proper error-propagation analysis
242 uint64_t virialUlpToleranceFactor = energyUlpToleranceFactor * 2;
243 auto virialTolerance = relativeToleranceAsPrecisionDependentUlp(
245 virialUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
246 virialUlpToleranceFactor * splineModuliDoublePrecisionUlps);
247 TestReferenceChecker virialChecker(
248 checker.checkCompound("Matrix", "Virial"));
249 virialChecker.setDefaultTolerance(virialTolerance);
250 for (int i = 0; i < DIM; i++)
252 for (int j = 0; j <= i; j++)
254 std::string valueId = formatString("Cell %d %d", i, j);
255 virialChecker.checkReal(virial[i][j], valueId.c_str());
264 static std::vector<std::unique_ptr<PmeTestHardwareContext>> s_pmeTestHardwareContexts;
267 std::vector<std::unique_ptr<PmeTestHardwareContext>> PmeSolveTest::s_pmeTestHardwareContexts;
269 /*! \brief Test for PME solving */
270 TEST_P(PmeSolveTest, ReproducesOutputs)
272 EXPECT_NO_THROW_GMX(runTest());
275 /* Valid input instances */
277 //! A couple of valid inputs for boxes.
278 std::vector<Matrix3x3> const c_sampleBoxes{
280 Matrix3x3{ { 8.0F, 0.0F, 0.0F, 0.0F, 3.4F, 0.0F, 0.0F, 0.0F, 2.0F } },
282 Matrix3x3{ { 7.0F, 0.0F, 0.0F, 0.0F, 4.1F, 0.0F, 3.5F, 2.0F, 12.2F } },
285 //! A couple of valid inputs for grid sizes
286 std::vector<IVec> const c_sampleGridSizes{ IVec{ 16, 12, 28 }, IVec{ 9, 7, 23 } };
288 //! Moved out from instantiations for readability
289 const auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
290 //! Moved out from instantiations for readability
291 const auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
293 //! 2 sample complex grids - only non-zero values have to be listed
294 std::vector<SparseComplexGridValuesInput> const c_sampleGrids{
295 SparseComplexGridValuesInput{
296 { IVec{ 0, 0, 0 }, t_complex{ 3.5F, 6.7F } },
297 { IVec{ 7, 0, 0 }, t_complex{ -2.5F, -0.7F } },
298 { IVec{ 3, 5, 7 }, t_complex{ -0.006F, 1e-8F } },
299 { IVec{ 3, 1, 2 }, t_complex{ 0.6F, 7.9F } },
300 { IVec{ 6, 2, 4 }, t_complex{ 30.1F, 2.45F } },
302 SparseComplexGridValuesInput{
303 { IVec{ 0, 4, 0 }, t_complex{ 0.0F, 0.3F } },
304 { IVec{ 4, 2, 7 }, t_complex{ 13.76F, -40.0F } },
305 { IVec{ 0, 6, 7 }, t_complex{ 3.6F, 0.0F } },
306 { IVec{ 2, 5, 10 }, t_complex{ 3.6F, 10.65F } },
310 //! Moved out from instantiations for readability
311 const auto c_inputGrids = ::testing::ValuesIn(c_sampleGrids);
312 //! Moved out from instantiations for readability
313 const auto c_inputEpsilon_r = ::testing::Values(1.2);
314 //! Moved out from instantiations for readability
315 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
316 //! Moved out from instantiations for readability
317 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
318 //! Moved out from instantiations for readability
319 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
321 //! Instantiation of the PME solving test
322 INSTANTIATE_TEST_SUITE_P(SaneInput,
324 ::testing::Combine(c_inputBoxes,
329 c_inputEwaldCoeff_lj,
332 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
333 INSTANTIATE_TEST_SUITE_P(DifferentEwaldCoeffQ,
335 ::testing::Combine(c_inputBoxes,
339 ::testing::Values(0.4),
340 c_inputEwaldCoeff_lj,
341 ::testing::Values(PmeSolveAlgorithm::Coulomb)));
343 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
344 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
345 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
346 INSTANTIATE_TEST_SUITE_P(DifferentEwaldCoeffLJ,
348 ::testing::Combine(c_inputBoxes,
353 ::testing::Values(2.35),
354 ::testing::Values(PmeSolveAlgorithm::LennardJones)));
356 //! A few more instances to check that different epsilon_r actually affects results of all solvers
357 INSTANTIATE_TEST_SUITE_P(DifferentEpsilonR,
359 ::testing::Combine(c_inputBoxes,
362 testing::Values(1.9),
364 c_inputEwaldCoeff_lj,