Fix #3655
[gromacs.git] / src / programs / mdrun / tests / simulator.cpp
blob31be90d878ebe8ed8cc7a2039986f9e614e1b600
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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 /*! \internal \file
37 * \brief
38 * Tests to compare two simulators which are expected to be identical
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \author Pascal Merz <pascal.merz@me.com>
42 * \ingroup module_mdrun_integration_tests
44 #include "gmxpre.h"
46 #include "config.h"
48 #include "gromacs/topology/ifunc.h"
49 #include "gromacs/utility/stringutil.h"
51 #include "testutils/mpitest.h"
52 #include "testutils/setenv.h"
53 #include "testutils/simulationdatabase.h"
55 #include "moduletest.h"
56 #include "simulatorcomparison.h"
58 namespace gmx
60 namespace test
62 namespace
65 /*! \brief Test fixture base for two equivalent simulators
67 * This test ensures that two simulator code paths (called via different mdp
68 * options and/or environment variables) yield identical coordinate, velocity,
69 * box, force and energy trajectories, up to some (arbitrary) precision.
71 * These tests are useful to check that re-implementations of existing simulators
72 * are correct, and that different code paths expected to yield identical results
73 * are equivalent.
75 using SimulatorComparisonTestParams =
76 std::tuple<std::tuple<std::string, std::string, std::string, std::string>, std::string>;
77 class SimulatorComparisonTest :
78 public MdrunTestFixture,
79 public ::testing::WithParamInterface<SimulatorComparisonTestParams>
83 TEST_P(SimulatorComparisonTest, WithinTolerances)
85 const auto& params = GetParam();
86 const auto& mdpParams = std::get<0>(params);
87 const auto& simulationName = std::get<0>(mdpParams);
88 const auto& integrator = std::get<1>(mdpParams);
89 const auto& tcoupling = std::get<2>(mdpParams);
90 const auto& pcoupling = std::get<3>(mdpParams);
91 const auto& environmentVariable = std::get<1>(params);
93 // TODO At some point we should also test PME-only ranks.
94 const int numRanksAvailable = getNumberOfTestMpiRanks();
95 if (!isNumberOfPpRanksSupported(simulationName, numRanksAvailable))
97 fprintf(stdout,
98 "Test system '%s' cannot run with %d ranks.\n"
99 "The supported numbers are: %s\n",
100 simulationName.c_str(), numRanksAvailable,
101 reportNumbersOfPpRanksSupported(simulationName).c_str());
102 return;
105 if (integrator == "md-vv" && pcoupling == "Parrinello-Rahman")
107 // do_md calls this MTTK, requires Nose-Hoover, and
108 // does not work with constraints or anisotropically
109 return;
112 const auto hasConservedField = !(tcoupling == "no" && pcoupling == "no");
114 SCOPED_TRACE(formatString(
115 "Comparing two simulations of '%s' "
116 "with integrator '%s', '%s' temperature coupling, and '%s' pressure coupling "
117 "switching environment variable '%s'",
118 simulationName.c_str(), integrator.c_str(), tcoupling.c_str(), pcoupling.c_str(),
119 environmentVariable.c_str()));
121 const auto mdpFieldValues = prepareMdpFieldValues(simulationName.c_str(), integrator.c_str(),
122 tcoupling.c_str(), pcoupling.c_str());
124 EnergyTermsToCompare energyTermsToCompare{ {
125 { interaction_function[F_EPOT].longname, relativeToleranceAsPrecisionDependentUlp(10.0, 100, 80) },
126 { interaction_function[F_EKIN].longname, relativeToleranceAsPrecisionDependentUlp(60.0, 100, 80) },
127 { interaction_function[F_PRES].longname,
128 relativeToleranceAsPrecisionDependentFloatingPoint(10.0, 0.01, 0.001) },
129 } };
130 if (hasConservedField)
132 energyTermsToCompare.emplace(interaction_function[F_ECONSERVED].longname,
133 relativeToleranceAsPrecisionDependentUlp(50.0, 100, 80));
136 if (simulationName == "argon12")
138 // Without constraints, we can be more strict
139 energyTermsToCompare = { {
140 { interaction_function[F_EPOT].longname,
141 relativeToleranceAsPrecisionDependentUlp(10.0, 24, 80) },
142 { interaction_function[F_EKIN].longname,
143 relativeToleranceAsPrecisionDependentUlp(10.0, 24, 80) },
144 { interaction_function[F_PRES].longname,
145 relativeToleranceAsPrecisionDependentFloatingPoint(10.0, 0.001, 0.0001) },
146 } };
147 if (hasConservedField)
149 energyTermsToCompare.emplace(interaction_function[F_ECONSERVED].longname,
150 relativeToleranceAsPrecisionDependentUlp(10.0, 24, 80));
154 // Specify how trajectory frame matching must work.
155 const TrajectoryFrameMatchSettings trajectoryMatchSettings{ true,
156 true,
157 true,
158 ComparisonConditions::MustCompare,
159 ComparisonConditions::MustCompare,
160 ComparisonConditions::MustCompare,
161 MaxNumFrames::compareAllFrames() };
162 TrajectoryTolerances trajectoryTolerances = TrajectoryComparison::s_defaultTrajectoryTolerances;
163 if (simulationName != "argon12")
165 trajectoryTolerances.velocities = trajectoryTolerances.coordinates;
168 // Build the functor that will compare reference and test
169 // trajectory frames in the chosen way.
170 const TrajectoryComparison trajectoryComparison{ trajectoryMatchSettings, trajectoryTolerances };
172 // Set file names
173 const auto simulator1TrajectoryFileName = fileManager_.getTemporaryFilePath("sim1.trr");
174 const auto simulator1EdrFileName = fileManager_.getTemporaryFilePath("sim1.edr");
175 const auto simulator2TrajectoryFileName = fileManager_.getTemporaryFilePath("sim2.trr");
176 const auto simulator2EdrFileName = fileManager_.getTemporaryFilePath("sim2.edr");
178 // Run grompp
179 runner_.tprFileName_ = fileManager_.getTemporaryFilePath("sim.tpr");
180 runner_.useTopGroAndNdxFromDatabase(simulationName);
181 runner_.useStringAsMdpFile(prepareMdpFileContents(mdpFieldValues));
182 runGrompp(&runner_);
184 // Backup current state of environment variable and unset it
185 const char* environmentVariableBackup = getenv(environmentVariable.c_str());
186 gmxUnsetenv(environmentVariable.c_str());
188 // Do first mdrun
189 runner_.fullPrecisionTrajectoryFileName_ = simulator1TrajectoryFileName;
190 runner_.edrFileName_ = simulator1EdrFileName;
191 runMdrun(&runner_);
193 // Set environment variable
194 const int overWriteEnvironmentVariable = 1;
195 gmxSetenv(environmentVariable.c_str(), "ON", overWriteEnvironmentVariable);
197 // Do second mdrun
198 runner_.fullPrecisionTrajectoryFileName_ = simulator2TrajectoryFileName;
199 runner_.edrFileName_ = simulator2EdrFileName;
200 runMdrun(&runner_);
202 // Reset or unset environment variable to leave further tests undisturbed
203 if (environmentVariableBackup != nullptr)
205 // set environment variable
206 gmxSetenv(environmentVariable.c_str(), environmentVariableBackup, overWriteEnvironmentVariable);
208 else
210 // unset environment variable
211 gmxUnsetenv(environmentVariable.c_str());
214 // Compare simulation results
215 compareEnergies(simulator1EdrFileName, simulator2EdrFileName, energyTermsToCompare);
216 compareTrajectories(simulator1TrajectoryFileName, simulator2TrajectoryFileName, trajectoryComparison);
219 // TODO: The time for OpenCL kernel compilation means these tests time
220 // out. Once that compilation is cached for the whole process, these
221 // tests can run in such configurations.
222 // These tests are very sensitive, so we only run them in double precision.
223 // As we change call ordering, they might actually become too strict to be useful.
224 #if !GMX_GPU_OPENCL && GMX_DOUBLE
225 INSTANTIATE_TEST_CASE_P(
226 SimulatorsAreEquivalentDefaultModular,
227 SimulatorComparisonTest,
228 ::testing::Combine(::testing::Combine(::testing::Values("argon12", "tip3p5"),
229 ::testing::Values("md-vv"),
230 ::testing::Values("no", "v-rescale", "berendsen"),
231 ::testing::Values("no")),
232 ::testing::Values("GMX_DISABLE_MODULAR_SIMULATOR")));
233 INSTANTIATE_TEST_CASE_P(
234 SimulatorsAreEquivalentDefaultLegacy,
235 SimulatorComparisonTest,
236 ::testing::Combine(::testing::Combine(::testing::Values("argon12", "tip3p5"),
237 ::testing::Values("md"),
238 ::testing::Values("no", "v-rescale", "berendsen"),
239 ::testing::Values("no", "Parrinello-Rahman")),
240 ::testing::Values("GMX_USE_MODULAR_SIMULATOR")));
241 #else
242 INSTANTIATE_TEST_CASE_P(
243 DISABLED_SimulatorsAreEquivalentDefaultModular,
244 SimulatorComparisonTest,
245 ::testing::Combine(::testing::Combine(::testing::Values("argon12", "tip3p5"),
246 ::testing::Values("md-vv"),
247 ::testing::Values("no", "v-rescale", "berendsen"),
248 ::testing::Values("no")),
249 ::testing::Values("GMX_DISABLE_MODULAR_SIMULATOR")));
250 INSTANTIATE_TEST_CASE_P(
251 DISABLED_SimulatorsAreEquivalentDefaultLegacy,
252 SimulatorComparisonTest,
253 ::testing::Combine(::testing::Combine(::testing::Values("argon12", "tip3p5"),
254 ::testing::Values("md"),
255 ::testing::Values("no", "v-rescale", "berendsen"),
256 ::testing::Values("no", "Parrinello-Rahman")),
257 ::testing::Values("GMX_USE_MODULAR_SIMULATOR")));
258 #endif
260 } // namespace
261 } // namespace test
262 } // namespace gmx