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) =
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: break;
103 case PmeSolveAlgorithm::LennardJones: inputRec.vdwtype = evdwPME; break;
105 default: GMX_THROW(InternalError("Unknown PME solver"));
108 TestReferenceData refData;
109 for (const auto& context : getPmeTestEnv()->getHardwareContexts())
111 CodePath codePath = context->getCodePath();
112 const bool supportedInput =
113 pmeSupportsInputForMode(*getPmeTestEnv()->hwinfo(), &inputRec, codePath);
116 /* Testing the failure for the unsupported input */
117 EXPECT_THROW(pmeInitEmpty(&inputRec, codePath, nullptr, nullptr, box, ewaldCoeff_q,
119 NotImplementedError);
123 std::map<GridOrdering, std::string> gridOrderingsToTest = { { GridOrdering::YZX,
125 if (codePath == CodePath::GPU)
127 gridOrderingsToTest[GridOrdering::XYZ] = "XYZ";
129 for (const auto& gridOrdering : gridOrderingsToTest)
131 for (bool computeEnergyAndVirial : { false, true })
133 /* Describing the test*/
134 SCOPED_TRACE(formatString(
135 "Testing solving (%s, %s, %s energy/virial) with %s %sfor PME grid "
136 "size %d %d %d, Ewald coefficients %g %g",
137 (method == PmeSolveAlgorithm::LennardJones) ? "Lennard-Jones" : "Coulomb",
138 gridOrdering.second.c_str(), computeEnergyAndVirial ? "with" : "without",
139 codePathToString(codePath), context->getDescription().c_str(),
140 gridSize[XX], gridSize[YY], gridSize[ZZ], ewaldCoeff_q, ewaldCoeff_lj));
142 /* Running the test */
143 PmeSafePointer pmeSafe =
144 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,
150 computeEnergyAndVirial);
151 pmeFinalizeTest(pmeSafe.get(), codePath);
153 /* Check the outputs */
154 TestReferenceChecker checker(refData.rootChecker());
156 SparseComplexGridValuesOutput nonZeroGridValuesOutput =
157 pmeGetComplexGrid(pmeSafe.get(), codePath, gridOrdering.first);
158 /* Transformed grid */
159 TestReferenceChecker gridValuesChecker(
160 checker.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
162 real gridValuesMagnitude = 1.0;
163 for (const auto& point : nonZeroGridValuesOutput)
165 gridValuesMagnitude = std::max(std::fabs(point.second.re), gridValuesMagnitude);
166 gridValuesMagnitude = std::max(std::fabs(point.second.im), gridValuesMagnitude);
168 // Spline moduli participate 3 times in the computation; 2 is an additional factor for SIMD exp() precision
169 uint64_t gridUlpToleranceFactor = DIM * 2;
170 if (method == PmeSolveAlgorithm::LennardJones)
172 // Lennard Jones is more complex and also uses erfc(), relax more
173 gridUlpToleranceFactor *= 2;
175 const uint64_t splineModuliDoublePrecisionUlps =
176 getSplineModuliDoublePrecisionUlps(inputRec.pme_order + 1);
177 auto gridTolerance = relativeToleranceAsPrecisionDependentUlp(
178 gridValuesMagnitude, gridUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
179 gridUlpToleranceFactor * splineModuliDoublePrecisionUlps);
180 gridValuesChecker.setDefaultTolerance(gridTolerance);
182 for (const auto& point : nonZeroGridValuesOutput)
184 // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
185 // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
186 if (fabs(point.second.re) >= GMX_FLOAT_MIN)
188 gridValuesChecker.checkReal(point.second.re, (point.first + " re").c_str());
190 if (fabs(point.second.im) >= GMX_FLOAT_MIN)
192 gridValuesChecker.checkReal(point.second.im, (point.first + " im").c_str());
196 if (computeEnergyAndVirial)
198 // Extract the energy and virial
200 pmeGetReciprocalEnergyAndVirial(pmeSafe.get(), codePath, method);
201 const auto& energy = (method == PmeSolveAlgorithm::Coulomb)
202 ? output.coulombEnergy_
203 : output.lennardJonesEnergy_;
204 const auto& virial = (method == PmeSolveAlgorithm::Coulomb)
205 ? output.coulombVirial_
206 : output.lennardJonesVirial_;
208 // These quantities are computed based on the grid values, so must have
209 // checking relative tolerances at least as large. Virial needs more flops
210 // than energy, so needs a larger tolerance.
213 double energyMagnitude = 10.0;
214 // TODO This factor is arbitrary, do a proper error-propagation analysis
215 uint64_t energyUlpToleranceFactor = gridUlpToleranceFactor * 2;
216 auto energyTolerance = relativeToleranceAsPrecisionDependentUlp(
217 energyMagnitude, energyUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
218 energyUlpToleranceFactor * splineModuliDoublePrecisionUlps);
219 TestReferenceChecker energyChecker(checker);
220 energyChecker.setDefaultTolerance(energyTolerance);
221 energyChecker.checkReal(energy, "Energy");
224 double virialMagnitude = 1000.0;
225 // TODO This factor is arbitrary, do a proper error-propagation analysis
226 uint64_t virialUlpToleranceFactor = energyUlpToleranceFactor * 2;
227 auto virialTolerance = relativeToleranceAsPrecisionDependentUlp(
228 virialMagnitude, virialUlpToleranceFactor * c_splineModuliSinglePrecisionUlps,
229 virialUlpToleranceFactor * splineModuliDoublePrecisionUlps);
230 TestReferenceChecker virialChecker(
231 checker.checkCompound("Matrix", "Virial"));
232 virialChecker.setDefaultTolerance(virialTolerance);
233 for (int i = 0; i < DIM; i++)
235 for (int j = 0; j <= i; j++)
237 std::string valueId = formatString("Cell %d %d", i, j);
238 virialChecker.checkReal(virial[i][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 std::vector<Matrix3x3> const c_sampleBoxes{
259 Matrix3x3{ { 8.0F, 0.0F, 0.0F, 0.0F, 3.4F, 0.0F, 0.0F, 0.0F, 2.0F } },
261 Matrix3x3{ { 7.0F, 0.0F, 0.0F, 0.0F, 4.1F, 0.0F, 3.5F, 2.0F, 12.2F } },
264 //! A couple of valid inputs for grid sizes
265 std::vector<IVec> const c_sampleGridSizes{ IVec{ 16, 12, 28 }, IVec{ 9, 7, 23 } };
267 //! Moved out from instantiations for readability
268 const auto c_inputBoxes = ::testing::ValuesIn(c_sampleBoxes);
269 //! Moved out from instantiations for readability
270 const auto c_inputGridSizes = ::testing::ValuesIn(c_sampleGridSizes);
272 //! 2 sample complex grids - only non-zero values have to be listed
273 std::vector<SparseComplexGridValuesInput> const c_sampleGrids{
274 SparseComplexGridValuesInput{
275 { IVec{ 0, 0, 0 }, t_complex{ 3.5F, 6.7F } },
276 { IVec{ 7, 0, 0 }, t_complex{ -2.5F, -0.7F } },
277 { IVec{ 3, 5, 7 }, t_complex{ -0.006F, 1e-8F } },
278 { IVec{ 3, 1, 2 }, t_complex{ 0.6F, 7.9F } },
279 { IVec{ 6, 2, 4 }, t_complex{ 30.1F, 2.45F } },
281 SparseComplexGridValuesInput{
282 { IVec{ 0, 4, 0 }, t_complex{ 0.0F, 0.3F } },
283 { IVec{ 4, 2, 7 }, t_complex{ 13.76F, -40.0F } },
284 { IVec{ 0, 6, 7 }, t_complex{ 3.6F, 0.0F } },
285 { IVec{ 2, 5, 10 }, t_complex{ 3.6F, 10.65F } },
289 //! Moved out from instantiations for readability
290 const auto c_inputGrids = ::testing::ValuesIn(c_sampleGrids);
291 //! Moved out from instantiations for readability
292 const auto c_inputEpsilon_r = ::testing::Values(1.2);
293 //! Moved out from instantiations for readability
294 const auto c_inputEwaldCoeff_q = ::testing::Values(2.0);
295 //! Moved out from instantiations for readability
296 const auto c_inputEwaldCoeff_lj = ::testing::Values(0.7);
297 //! Moved out from instantiations for readability
298 const auto c_inputMethods = ::testing::Values(PmeSolveAlgorithm::Coulomb, PmeSolveAlgorithm::LennardJones);
300 //! Instantiation of the PME solving test
301 INSTANTIATE_TEST_CASE_P(SaneInput,
303 ::testing::Combine(c_inputBoxes,
308 c_inputEwaldCoeff_lj,
311 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
312 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffQ,
314 ::testing::Combine(c_inputBoxes,
318 ::testing::Values(0.4),
319 c_inputEwaldCoeff_lj,
320 ::testing::Values(PmeSolveAlgorithm::Coulomb)));
322 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
323 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
324 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
325 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffLJ,
327 ::testing::Combine(c_inputBoxes,
332 ::testing::Values(2.35),
333 ::testing::Values(PmeSolveAlgorithm::LennardJones)));
335 //! A few more instances to check that different epsilon_r actually affects results of all solvers
336 INSTANTIATE_TEST_CASE_P(DifferentEpsilonR,
338 ::testing::Combine(c_inputBoxes,
341 testing::Values(1.9),
343 c_inputEwaldCoeff_lj,