2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018, 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 for (const auto &context
: getPmeTestEnv()->getHardwareContexts())
115 CodePath codePath
= context
->getCodePath();
116 const bool supportedInput
= pmeSupportsInputForMode(&inputRec
, codePath
);
119 /* Testing the failure for the unsupported input */
120 EXPECT_THROW(pmeInitEmpty(&inputRec
, codePath
, nullptr, nullptr, box
, ewaldCoeff_q
, ewaldCoeff_lj
), NotImplementedError
);
124 std::map
<GridOrdering
, std::string
> gridOrderingsToTest
= {{GridOrdering::YZX
, "YZX"}};
125 if (codePath
== CodePath::CUDA
)
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("Testing solving (%s, %s, %s energy/virial) with %s %sfor PME grid size %d %d %d, Ewald coefficients %g %g",
135 (method
== PmeSolveAlgorithm::LennardJones
) ? "Lennard-Jones" : "Coulomb",
136 gridOrdering
.second
.c_str(),
137 computeEnergyAndVirial
? "with" : "without",
138 codePathToString(codePath
),
139 context
->getDescription().c_str(),
140 gridSize
[XX
], gridSize
[YY
], gridSize
[ZZ
],
141 ewaldCoeff_q
, ewaldCoeff_lj
144 /* Running the test */
145 PmeSafePointer pmeSafe
= pmeInitEmpty(&inputRec
, codePath
, context
->getDeviceInfo(),
146 context
->getPmeGpuProgram(), box
, ewaldCoeff_q
, ewaldCoeff_lj
);
147 pmeSetComplexGrid(pmeSafe
.get(), codePath
, gridOrdering
.first
, nonZeroGridValues
);
148 const real cellVolume
= box
[0] * box
[4] * box
[8];
149 //FIXME - this is box[XX][XX] * box[YY][YY] * box[ZZ][ZZ], should be stored in the PME structure
150 pmePerformSolve(pmeSafe
.get(), codePath
, method
, cellVolume
, gridOrdering
.first
, computeEnergyAndVirial
);
151 pmeFinalizeTest(pmeSafe
.get(), codePath
);
153 /* Check the outputs */
154 TestReferenceChecker
checker(refData
.rootChecker());
156 SparseComplexGridValuesOutput nonZeroGridValuesOutput
= pmeGetComplexGrid(pmeSafe
.get(), codePath
, gridOrdering
.first
);
157 /* Transformed grid */
158 TestReferenceChecker
gridValuesChecker(checker
.checkCompound("NonZeroGridValues", "ComplexSpaceGrid"));
160 real gridValuesMagnitude
= 1.0;
161 for (const auto &point
: nonZeroGridValuesOutput
)
163 gridValuesMagnitude
= std::max(std::fabs(point
.second
.re
), gridValuesMagnitude
);
164 gridValuesMagnitude
= std::max(std::fabs(point
.second
.im
), gridValuesMagnitude
);
166 // Spline moduli participate 3 times in the computation; 2 is an additional factor for SIMD exp() precision
167 gmx_uint64_t gridUlpToleranceFactor
= DIM
* 2;
168 if (method
== PmeSolveAlgorithm::LennardJones
)
170 // Lennard Jones is more complex and also uses erfc(), relax more
171 gridUlpToleranceFactor
*= 2;
173 const gmx_uint64_t splineModuliDoublePrecisionUlps
174 = getSplineModuliDoublePrecisionUlps(inputRec
.pme_order
+ 1);
176 = relativeToleranceAsPrecisionDependentUlp(gridValuesMagnitude
,
177 gridUlpToleranceFactor
* c_splineModuliSinglePrecisionUlps
,
178 gridUlpToleranceFactor
* splineModuliDoublePrecisionUlps
);
179 gridValuesChecker
.setDefaultTolerance(gridTolerance
);
181 for (const auto &point
: nonZeroGridValuesOutput
)
183 // we want an additional safeguard for denormal numbers as they cause an exception in string conversion;
184 // however, using GMX_REAL_MIN causes an "unused item warning" for single precision builds
185 if (fabs(point
.second
.re
) >= GMX_FLOAT_MIN
)
187 gridValuesChecker
.checkReal(point
.second
.re
, (point
.first
+ " re").c_str());
189 if (fabs(point
.second
.im
) >= GMX_FLOAT_MIN
)
191 gridValuesChecker
.checkReal(point
.second
.im
, (point
.first
+ " im").c_str());
195 if (computeEnergyAndVirial
)
197 // Extract the energy and virial
200 std::tie(energy
, virial
) = pmeGetReciprocalEnergyAndVirial(pmeSafe
.get(), codePath
, method
);
202 // These quantities are computed based on the grid values, so must have
203 // checking relative tolerances at least as large. Virial needs more flops
204 // than energy, so needs a larger tolerance.
207 double energyMagnitude
= 10.0;
208 // TODO This factor is arbitrary, do a proper error-propagation analysis
209 gmx_uint64_t energyUlpToleranceFactor
= gridUlpToleranceFactor
* 2;
211 = relativeToleranceAsPrecisionDependentUlp(energyMagnitude
,
212 energyUlpToleranceFactor
* c_splineModuliSinglePrecisionUlps
,
213 energyUlpToleranceFactor
* splineModuliDoublePrecisionUlps
);
214 TestReferenceChecker
energyChecker(checker
);
215 energyChecker
.setDefaultTolerance(energyTolerance
);
216 energyChecker
.checkReal(energy
, "Energy");
219 double virialMagnitude
= 1000.0;
220 // TODO This factor is arbitrary, do a proper error-propagation analysis
221 gmx_uint64_t virialUlpToleranceFactor
= energyUlpToleranceFactor
* 2;
223 = relativeToleranceAsPrecisionDependentUlp(virialMagnitude
,
224 virialUlpToleranceFactor
* c_splineModuliSinglePrecisionUlps
,
225 virialUlpToleranceFactor
* splineModuliDoublePrecisionUlps
);
226 TestReferenceChecker
virialChecker(checker
.checkCompound("Matrix", "Virial"));
227 virialChecker
.setDefaultTolerance(virialTolerance
);
228 for (int i
= 0; i
< DIM
; i
++)
230 for (int j
= 0; j
<= i
; j
++)
232 std::string valueId
= formatString("Cell %d %d", i
, j
);
233 virialChecker
.checkReal(virial
[i
* DIM
+ j
], valueId
.c_str());
244 /*! \brief Test for PME solving */
245 TEST_P(PmeSolveTest
, ReproducesOutputs
)
247 EXPECT_NO_THROW(runTest());
250 /* Valid input instances */
252 //! A couple of valid inputs for boxes.
253 static std::vector
<Matrix3x3
> const c_sampleBoxes
269 //! A couple of valid inputs for grid sizes
270 static std::vector
<IVec
> const c_sampleGridSizes
280 //! Moved out from instantiations for readability
281 const auto c_inputBoxes
= ::testing::ValuesIn(c_sampleBoxes
);
282 //! Moved out from instantiations for readability
283 const auto c_inputGridSizes
= ::testing::ValuesIn(c_sampleGridSizes
);
285 //! 2 sample complex grids - only non-zero values have to be listed
286 static std::vector
<SparseComplexGridValuesInput
> const c_sampleGrids
288 SparseComplexGridValuesInput
{{
319 SparseComplexGridValuesInput
{{
346 //! Moved out from instantiations for readability
347 const auto c_inputGrids
= ::testing::ValuesIn(c_sampleGrids
);
348 //! Moved out from instantiations for readability
349 const auto c_inputEpsilon_r
= ::testing::Values(1.2);
350 //! Moved out from instantiations for readability
351 const auto c_inputEwaldCoeff_q
= ::testing::Values(2.0);
352 //! Moved out from instantiations for readability
353 const auto c_inputEwaldCoeff_lj
= ::testing::Values(0.7);
354 //! Moved out from instantiations for readability
355 const auto c_inputMethods
= ::testing::Values(PmeSolveAlgorithm::Coulomb
, PmeSolveAlgorithm::LennardJones
);
357 //! Instantiation of the PME solving test
358 INSTANTIATE_TEST_CASE_P(SaneInput
, PmeSolveTest
, ::testing::Combine(c_inputBoxes
, c_inputGridSizes
, c_inputGrids
,
359 c_inputEpsilon_r
, c_inputEwaldCoeff_q
, c_inputEwaldCoeff_lj
, c_inputMethods
));
361 //! A few more instances to check that different ewaldCoeff_q actually affects results of the Coulomb solver
362 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffQ
, PmeSolveTest
, ::testing::Combine(c_inputBoxes
, c_inputGridSizes
, c_inputGrids
,
363 c_inputEpsilon_r
, ::testing::Values(0.4), c_inputEwaldCoeff_lj
,
364 ::testing::Values(PmeSolveAlgorithm::Coulomb
)));
366 //! A few more instances to check that different ewaldCoeff_lj actually affects results of the Lennard-Jones solver.
367 //! The value has to be approximately larger than 1 / (box dimensions) to have a meaningful output grid.
368 //! Previous value of 0.3 caused one of the grid cells to be less or greater than GMX_FLOAT_MIN, depending on the architecture.
369 INSTANTIATE_TEST_CASE_P(DifferentEwaldCoeffLJ
, PmeSolveTest
, ::testing::Combine(c_inputBoxes
, c_inputGridSizes
, c_inputGrids
,
370 c_inputEpsilon_r
, c_inputEwaldCoeff_q
, ::testing::Values(2.35),
371 ::testing::Values(PmeSolveAlgorithm::LennardJones
)));
373 //! A few more instances to check that different epsilon_r actually affects results of all solvers
374 INSTANTIATE_TEST_CASE_P(DifferentEpsilonR
, PmeSolveTest
, ::testing::Combine(c_inputBoxes
, c_inputGridSizes
, c_inputGrids
,
375 testing::Values(1.9), c_inputEwaldCoeff_q
, c_inputEwaldCoeff_lj
,