Run simulator equivalence tests in double precision only
[gromacs.git] / src / programs / mdrun / tests / simulator.cpp
blob682865eaaf5ff8bd896c8e8a3028cb719c48a4e5
1 /*
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 /*! \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/trajectory/energyframe.h"
50 #include "gromacs/utility/stringutil.h"
52 #include "simulatorcomparison.h"
54 namespace gmx
56 namespace test
58 namespace
61 /*! \brief Test fixture base for two equivalent simulators
63 * This test ensures that two simulator code paths (called via different mdp
64 * options and/or environment variables) yield identical coordinate, velocity,
65 * box, force and energy trajectories, up to some (arbitrary) precision.
67 * These tests are useful to check that re-implementations of existing simulators
68 * are correct, and that different code paths expected to yield identical results
69 * are equivalent.
71 using SimulatorComparisonTestParams = std::tuple<
72 std::tuple < std::string, std::string, std::string>,
73 std::string>;
74 class SimulatorComparisonTest :
75 public MdrunTestFixture,
76 public ::testing::WithParamInterface < SimulatorComparisonTestParams >
80 TEST_P(SimulatorComparisonTest, WithinTolerances)
82 auto params = GetParam();
83 auto mdpParams = std::get<0>(params);
84 auto simulationName = std::get<0>(mdpParams);
85 auto integrator = std::get<1>(mdpParams);
86 auto tcoupling = std::get<2>(mdpParams);
87 auto envVariable = std::get<1>(params);
88 SCOPED_TRACE(formatString("Comparing two simulations of '%s' "
89 "with integrator '%s' and '%s' temperature coupling, "
90 "switching environment variable '%s'",
91 simulationName.c_str(), integrator.c_str(),
92 tcoupling.c_str(), envVariable.c_str()));
94 auto mdpFieldValues = prepareMdpFieldValues(simulationName.c_str(),
95 integrator.c_str(),
96 tcoupling.c_str(), "no");
98 EnergyTermsToCompare energyTermsToCompare
101 interaction_function[F_EPOT].longname,
102 relativeToleranceAsPrecisionDependentUlp(10.0, 100, 40)
105 interaction_function[F_EKIN].longname,
106 relativeToleranceAsPrecisionDependentUlp(60.0, 100, 40)
109 interaction_function[F_PRES].longname,
110 relativeToleranceAsPrecisionDependentFloatingPoint(10.0, 0.01, 0.001)
114 if (simulationName == "argon12")
116 // Without constraints, we can be more strict
117 energyTermsToCompare =
120 interaction_function[F_EPOT].longname,
121 relativeToleranceAsPrecisionDependentUlp(10.0, 24, 40)
124 interaction_function[F_EKIN].longname,
125 relativeToleranceAsPrecisionDependentUlp(10.0, 24, 40)
128 interaction_function[F_PRES].longname,
129 relativeToleranceAsPrecisionDependentFloatingPoint(10.0, 0.001, 0.0001)
134 // Specify how trajectory frame matching must work.
135 TrajectoryFrameMatchSettings trajectoryMatchSettings {
136 true, true, true,
137 ComparisonConditions::MustCompare,
138 ComparisonConditions::MustCompare,
139 ComparisonConditions::MustCompare
141 TrajectoryTolerances trajectoryTolerances = TrajectoryComparison::s_defaultTrajectoryTolerances;
142 if (simulationName != "argon12")
144 trajectoryTolerances.velocities = trajectoryTolerances.coordinates;
147 // Build the functor that will compare reference and test
148 // trajectory frames in the chosen way.
149 TrajectoryComparison trajectoryComparison {
150 trajectoryMatchSettings, trajectoryTolerances
153 int numWarningsToTolerate = 0;
154 executeSimulatorComparisonTest(
155 envVariable,
156 &fileManager_, &runner_,
157 simulationName, numWarningsToTolerate,
158 mdpFieldValues,
159 energyTermsToCompare,
160 trajectoryComparison);
163 // TODO: The time for OpenCL kernel compilation means these tests time
164 // out. Once that compilation is cached for the whole process, these
165 // tests can run in such configurations.
166 // These tests are very sensitive, so we only run them in double precision.
167 // As we change call ordering, they might actually become too strict to be useful.
168 #if GMX_GPU != GMX_GPU_OPENCL && GMX_DOUBLE
169 INSTANTIATE_TEST_CASE_P(SimulatorsAreEquivalentDefaultModular, SimulatorComparisonTest,
170 ::testing::Combine(
171 ::testing::Combine(::testing::Values("argon12", "tip3p5"),
172 ::testing::Values("md-vv"),
173 ::testing::Values("no", "v-rescale")),
174 ::testing::Values("GMX_DISABLE_MODULAR_SIMULATOR")));
175 INSTANTIATE_TEST_CASE_P(SimulatorsAreEquivalentDefaultLegacy, SimulatorComparisonTest,
176 ::testing::Combine(
177 ::testing::Combine(::testing::Values("argon12", "tip3p5"),
178 ::testing::Values("md"),
179 ::testing::Values("no", "v-rescale")),
180 ::testing::Values("GMX_USE_MODULAR_SIMULATOR")));
181 #else
182 INSTANTIATE_TEST_CASE_P(DISABLED_SimulatorsAreEquivalentDefaultModular, SimulatorComparisonTest,
183 ::testing::Combine(
184 ::testing::Combine(::testing::Values("argon12", "tip3p5"),
185 ::testing::Values("md-vv"),
186 ::testing::Values("no", "v-rescale")),
187 ::testing::Values("GMX_DISABLE_MODULAR_SIMULATOR")));
188 INSTANTIATE_TEST_CASE_P(DISABLED_SimulatorsAreEquivalentDefaultLegacy, SimulatorComparisonTest,
189 ::testing::Combine(
190 ::testing::Combine(::testing::Values("argon12", "tip3p5"),
191 ::testing::Values("md"),
192 ::testing::Values("no", "v-rescale")),
193 ::testing::Values("GMX_USE_MODULAR_SIMULATOR")));
194 #endif
196 } // namespace
197 } // namespace test
198 } // namespace gmx