2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
36 * \brief Tests for the Leap-Frog integrator
38 * The test creates a system of independent particles exerting constant
39 * external forces and makes several numerical integration timesteps.
40 * The results are compared with the analytical solution (for the systems
41 * without the temperature coupling) and with the pre-computed reference
42 * values. The tests use runners that are created for each available
43 * implementation of the tested algorithm.
45 * \todo Add tests for integrators with pressure control.
46 * \todo Add PBC handling test.
48 * \author Artem Zhmurov <zhmurov@gmail.com>
49 * \ingroup module_mdlib
61 #include <unordered_map>
64 #include <gtest/gtest.h>
66 #include "gromacs/gpu_utils/gpu_testutils.h"
67 #include "gromacs/math/vec.h"
68 #include "gromacs/math/vectypes.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/utility/stringutil.h"
72 #include "testutils/refdata.h"
73 #include "testutils/testasserts.h"
75 #include "leapfrogtestdata.h"
76 #include "leapfrogtestrunners.h"
85 /*! \brief The parameters for the test.
87 * The test will run for combinations of:
92 * 4. Velocity components
94 * 6. Number of temperature coupling groups
96 struct LeapFrogTestParameters
98 //! Total number of atoms
102 //! Number of integration steps
108 //! Number of temperature coupling group (zero for no temperature coupling)
109 int numTCoupleGroups;
110 //! Number of steps between pressure coupling steps (zero for no pressure coupling).
114 //! The set of parameters combinations to run the test on
115 const LeapFrogTestParameters parametersSets[] = {
116 { 1, 0.001, 1, { 0.0, 0.0, 0.0 }, { 0.0, 0.0, 0.0 }, 0, 0 }, // Zero velocity and force
117 { 1, 0.001, 1, { 0.0, 0.0, 0.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // Zero velocity
118 { 1, 0.001, 1, { 1.0, -2.0, 3.0 }, { 0.0, 0.0, 0.0 }, 0, 0 }, // Zero force
119 { 1, 0.001, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // 1 particle
120 { 10, 0.001, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // 10 particles
121 { 100, 0.001, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // 100 particles
122 { 300, 0.001, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // 300 particles
123 { 1, 0.0005, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // 0.0005 ps timestep
124 { 1, 0.001, 10, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // 10 step
125 { 1, 0.001, 100, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 0 }, // 100 steps
126 { 100, 0.001, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 1, 0 }, // 1 temperature couple group
127 { 100, 0.001, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 2, 0 }, // 2 temperature couple groups
128 { 100, 0.001, 1, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 10, 0 }, // 10 temperature couple groups
129 { 100, 0.001, 10, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 1 }, // With pressure coupling
130 { 100, 0.001, 10, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 2, 1 }, // With both temperature and pressure coupling
131 { 100, 0.001, 10, { 1.0, -2.0, 3.0 }, { -3.0, 2.0, -1.0 }, 0, 3 }
132 }; // Do pressure coupling not on every step
135 /*! \brief Test fixture for LeapFrog integrator.
137 class LeapFrogTest : public ::testing::TestWithParam<LeapFrogTestParameters>
140 //! Availiable runners (CPU and GPU versions of the Leap-Frog)
141 static std::unordered_map<std::string, void (*)(LeapFrogTestData* testData, const int numSteps)> s_runners_;
143 TestReferenceData refData_;
144 //! Checker for reference data
145 TestReferenceChecker checker_;
147 LeapFrogTest() : checker_(refData_.rootChecker()) {}
149 //! Setup the runners one for all parameters sets
150 static void SetUpTestCase()
153 // All runners should be registered here under appropriate conditions
155 s_runners_["LeapFrogSimple"] = integrateLeapFrogSimple;
156 if (GMX_GPU == GMX_GPU_CUDA && canComputeOnGpu())
158 s_runners_["LeapFrogGpu"] = integrateLeapFrogGpu;
162 /*! \brief Test the numerical integrator against analytical solution for simple constant force case.
164 * \param[in] tolerance Tolerance
165 * \param[in] testData Test data object
166 * \param[in] totalTime Total numerical integration time
168 void testAgainstAnalyticalSolution(FloatingPointTolerance tolerance,
169 const LeapFrogTestData& testData,
170 const real totalTime)
172 for (int i = 0; i < testData.numAtoms_; i++)
176 for (int d = 0; d < DIM; d++)
178 // Analytical solution for constant-force particle movement
179 real x0 = testData.x0_[i][d];
180 real v0 = testData.v0_[i][d];
181 real f = testData.f_[i][d];
182 real im = testData.inverseMasses_[i];
184 xAnalytical[d] = x0 + v0 * totalTime + 0.5 * f * totalTime * totalTime * im;
185 vAnalytical[d] = v0 + f * totalTime * im;
187 EXPECT_REAL_EQ_TOL(xAnalytical[d], testData.xPrime_[i][d], tolerance) << gmx::formatString(
188 "Coordinate %d of atom %d is different from analytical solution.", d, i);
190 EXPECT_REAL_EQ_TOL(vAnalytical[d], testData.v_[i][d], tolerance) << gmx::formatString(
191 "Velocity component %d of atom %d is different from analytical solution.", d, i);
196 /*! \brief Test the numerical integrator against pre-computed reference values.
198 * \param[in] testData Test data object
200 void testAgainstReferenceData(const LeapFrogTestData& testData)
202 TestReferenceChecker finalPositionsRef(
203 checker_.checkSequenceCompound("FinalPositions", testData.numAtoms_));
204 for (int i = 0; i < testData.numAtoms_; i++)
206 const gmx::RVec& xPrime = testData.xPrime_[i];
207 TestReferenceChecker xPrimeRef(finalPositionsRef.checkCompound("Atom", nullptr));
208 xPrimeRef.checkReal(xPrime[XX], "XX");
209 xPrimeRef.checkReal(xPrime[YY], "YY");
210 xPrimeRef.checkReal(xPrime[ZZ], "ZZ");
213 TestReferenceChecker finalVelocitiesRef(
214 checker_.checkSequenceCompound("FinalVelocities", testData.numAtoms_));
215 for (int i = 0; i < testData.numAtoms_; i++)
217 const gmx::RVec& v = testData.v_[i];
218 TestReferenceChecker vRef(finalVelocitiesRef.checkCompound("Atom", nullptr));
219 vRef.checkReal(v[XX], "XX");
220 vRef.checkReal(v[YY], "YY");
221 vRef.checkReal(v[ZZ], "ZZ");
226 std::unordered_map<std::string, void (*)(LeapFrogTestData* testData, const int numSteps)> LeapFrogTest::s_runners_;
228 TEST_P(LeapFrogTest, SimpleIntegration)
230 // Cycle through all available runners
231 for (const auto& runner : s_runners_)
233 std::string runnerName = runner.first;
235 LeapFrogTestParameters parameters = GetParam();
237 std::string testDescription = formatString(
238 "Testing %s with %d atoms for %d timesteps with %d temperature coupling groups and "
239 "%s pressure coupling (dt = %f, v0=(%f, %f, %f), f0=(%f, %f, %f), nstpcouple = %d)",
240 runnerName.c_str(), parameters.numAtoms, parameters.numSteps,
241 parameters.numTCoupleGroups, parameters.nstpcouple == 0 ? "without" : "with",
242 parameters.timestep, parameters.v[XX], parameters.v[YY], parameters.v[ZZ],
243 parameters.f[XX], parameters.f[YY], parameters.f[ZZ], parameters.nstpcouple);
244 SCOPED_TRACE(testDescription);
246 std::unique_ptr<LeapFrogTestData> testData = std::make_unique<LeapFrogTestData>(
247 parameters.numAtoms, parameters.timestep, parameters.v, parameters.f,
248 parameters.numTCoupleGroups, parameters.nstpcouple);
250 runner.second(testData.get(), parameters.numSteps);
252 real totalTime = parameters.numSteps * parameters.timestep;
253 // TODO For the case of constant force, the numerical scheme is exact and
254 // the only source of errors is floating point arithmetic. Hence,
255 // the tolerance can be calculated.
256 FloatingPointTolerance tolerance = absoluteTolerance(parameters.numSteps * 0.000005);
258 // Test against the analytical solution (without temperature coupling)
259 if (parameters.numTCoupleGroups == 0 && parameters.nstpcouple == 0)
261 testAgainstAnalyticalSolution(tolerance, *testData, totalTime);
264 checker_.setDefaultTolerance(tolerance);
265 testAgainstReferenceData(*testData);
269 INSTANTIATE_TEST_CASE_P(WithParameters, LeapFrogTest, ::testing::ValuesIn(parametersSets));