Update bundled GoogleTest to current HEAD
[alexxy/gromacs.git] / src / programs / mdrun / tests / rerun.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  * \brief
39  * Tests for the mdrun -rerun functionality
40  *
41  * \author Mark Abraham <mark.j.abraham@gmail.com>
42  * \ingroup module_mdrun_integration_tests
43  */
44 #include "gmxpre.h"
45
46 #include "config.h"
47
48 #include "gromacs/topology/ifunc.h"
49 #include "gromacs/utility/strconvert.h"
50 #include "gromacs/utility/stringutil.h"
51
52 #include "testutils/mpitest.h"
53 #include "testutils/simulationdatabase.h"
54
55 #include "moduletest.h"
56 #include "simulatorcomparison.h"
57
58 namespace gmx
59 {
60 namespace test
61 {
62 namespace
63 {
64 /*! \brief Test fixture base for mdrun -rerun
65  *
66  * This test ensures mdrun can run a simulation, writing a trajectory
67  * and matching energies, and reproduce the same energies from a rerun
68  * to within a tight tolerance. It says nothing about whether a rerun
69  * can reproduce energies from a trajectory generated with older code,
70  * since that is not a useful property. Whether mdrun produced correct
71  * energies then and now needs different kinds of testing, but if
72  * true, this test ensures the rerun has the expected property.
73  *
74  * The limitations of mdrun and its output means that reproducing the
75  * same energies is currently only meaningful for integration without
76  * thermostats or barostats, however the present form of the test
77  * infrastructure has in-principle support for such, if that is ever
78  * needed/useful.
79  *
80  * We should also not compare pressure, because with constraints the
81  * non-search steps need a much larger tolerance, and per Issue #1868
82  * we should stop computing pressure in reruns anyway.
83  *
84  * Similarly, per 1868, in the present implementation the kinetic
85  * energy quantities are not generally reproducible, either.
86  *
87  * The choices for tolerance are arbitrary but sufficient.  Rerun does
88  * pair search every frame, so it cannot in general exactly reproduce
89  * quantities from a normal run, because the accumulation order
90  * differs. (Nor does it reproduce pair-search frames exactly,
91  * either). */
92 class MdrunRerunTest :
93     public MdrunTestFixture,
94     public ::testing::WithParamInterface<std::tuple<std::string, std::string>>
95 {
96 public:
97     //! Trajectory components to compare
98     static const TrajectoryFrameMatchSettings trajectoryMatchSettings;
99 };
100
101 // Compare box, positions and forces, but not velocities
102 // (velocities are ignored in reruns)
103 const TrajectoryFrameMatchSettings MdrunRerunTest::trajectoryMatchSettings = {
104     true,
105     true,
106     true,
107     ComparisonConditions::MustCompare,
108     ComparisonConditions::NoComparison,
109     ComparisonConditions::MustCompare,
110     MaxNumFrames::compareAllFrames()
111 };
112
113 void executeRerunTest(TestFileManager*            fileManager,
114                       SimulationRunner*           runner,
115                       const std::string&          simulationName,
116                       int                         numWarningsToTolerate,
117                       const MdpFieldValues&       mdpFieldValues,
118                       const EnergyTermsToCompare& energyTermsToCompare,
119                       const TrajectoryComparison& trajectoryComparison)
120 {
121     // TODO At some point we should also test PME-only ranks.
122     int numRanksAvailable = getNumberOfTestMpiRanks();
123     if (!isNumberOfPpRanksSupported(simulationName, numRanksAvailable))
124     {
125         fprintf(stdout,
126                 "Test system '%s' cannot run with %d ranks.\n"
127                 "The supported numbers are: %s\n",
128                 simulationName.c_str(),
129                 numRanksAvailable,
130                 reportNumbersOfPpRanksSupported(simulationName).c_str());
131         return;
132     }
133
134     // Set file names
135     auto simulator1TrajectoryFileName = fileManager->getTemporaryFilePath("sim1.trr");
136     auto simulator1EdrFileName        = fileManager->getTemporaryFilePath("sim1.edr");
137     auto simulator2TrajectoryFileName = fileManager->getTemporaryFilePath("sim2.trr");
138     auto simulator2EdrFileName        = fileManager->getTemporaryFilePath("sim2.edr");
139
140     // Run grompp
141     runner->tprFileName_ = fileManager->getTemporaryFilePath("sim.tpr");
142     runner->useTopGroAndNdxFromDatabase(simulationName);
143     runner->useStringAsMdpFile(prepareMdpFileContents(mdpFieldValues));
144     auto options = std::vector<SimulationOptionTuple>();
145     if (numWarningsToTolerate > 0)
146     {
147         options.emplace_back(SimulationOptionTuple("-maxwarn", std::to_string(numWarningsToTolerate)));
148     }
149     runGrompp(runner, options);
150
151     // Do first mdrun
152     runner->fullPrecisionTrajectoryFileName_ = simulator1TrajectoryFileName;
153     runner->edrFileName_                     = simulator1EdrFileName;
154     runMdrun(runner);
155
156     // Do second mdrun
157     runner->fullPrecisionTrajectoryFileName_ = simulator2TrajectoryFileName;
158     runner->edrFileName_                     = simulator2EdrFileName;
159     runMdrun(runner, { SimulationOptionTuple("-rerun", simulator1TrajectoryFileName) });
160
161     // Compare simulation results
162     compareEnergies(simulator1EdrFileName, simulator2EdrFileName, energyTermsToCompare);
163     compareTrajectories(simulator1TrajectoryFileName, simulator2TrajectoryFileName, trajectoryComparison);
164 }
165
166 TEST_P(MdrunRerunTest, WithinTolerances)
167 {
168     auto params         = GetParam();
169     auto simulationName = std::get<0>(params);
170     auto integrator     = std::get<1>(params);
171     SCOPED_TRACE(
172             formatString("Comparing normal and rerun of simulation '%s' "
173                          "with integrator '%s'",
174                          simulationName.c_str(),
175                          integrator.c_str()));
176
177     auto mdpFieldValues =
178             prepareMdpFieldValues(simulationName.c_str(), integrator.c_str(), "no", "no");
179
180     // bd is much less reproducible in a rerun than the other integrators
181     const int            toleranceScaleFactor = (integrator == "bd") ? 2 : 1;
182     EnergyTermsToCompare energyTermsToCompare{ {
183             { interaction_function[F_EPOT].longname,
184               relativeToleranceAsPrecisionDependentUlp(
185                       10.0, 24 * toleranceScaleFactor, 40 * toleranceScaleFactor) },
186     } };
187
188     // Specify how trajectory frame matching must work
189     TrajectoryComparison trajectoryComparison{ trajectoryMatchSettings,
190                                                TrajectoryComparison::s_defaultTrajectoryTolerances };
191
192     int numWarningsToTolerate = 0;
193     executeRerunTest(&fileManager_,
194                      &runner_,
195                      simulationName,
196                      numWarningsToTolerate,
197                      mdpFieldValues,
198                      energyTermsToCompare,
199                      trajectoryComparison);
200 }
201
202 // TODO The time for OpenCL kernel compilation means these tests time
203 // out. Once that compilation is cached for the whole process, these
204 // tests can run in such configurations.
205 #if !GMX_GPU_OPENCL
206 INSTANTIATE_TEST_SUITE_P(
207         NormalMdrunIsReproduced,
208         MdrunRerunTest,
209         ::testing::Combine(::testing::Values("argon12", "tip3p5", "alanine_vsite_vacuo"),
210                            ::testing::Values("md", "md-vv", "bd", "sd")));
211 #else
212 INSTANTIATE_TEST_SUITE_P(
213         DISABLED_NormalMdrunIsReproduced,
214         MdrunRerunTest,
215         ::testing::Combine(::testing::Values("argon12", "tip3p5", "alanine_vsite_vacuo"),
216                            ::testing::Values("md", "md-vv", "bd", "sd")));
217 #endif
218
219 class MdrunRerunFreeEnergyTest :
220     public MdrunTestFixture,
221     public ::testing::WithParamInterface<std::tuple<std::string, std::string, int>>
222 {
223 };
224
225 TEST_P(MdrunRerunFreeEnergyTest, WithinTolerances)
226 {
227     auto params          = GetParam();
228     auto simulationName  = std::get<0>(params);
229     auto integrator      = std::get<1>(params);
230     auto initLambdaState = std::get<2>(params);
231     SCOPED_TRACE(
232             formatString("Comparing normal and rerun of simulation '%s' "
233                          "with integrator '%s' for initial lambda state %d",
234                          simulationName.c_str(),
235                          integrator.c_str(),
236                          initLambdaState));
237
238     auto mdpFieldValues =
239             prepareMdpFieldValues(simulationName.c_str(), integrator.c_str(), "no", "no");
240     mdpFieldValues["init-lambda-state"] = toString(initLambdaState);
241
242     EnergyTermsToCompare energyTermsToCompare{
243         { { interaction_function[F_EPOT].longname, relativeToleranceAsPrecisionDependentUlp(10.0, 24, 32) },
244           { interaction_function[F_DVDL_COUL].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) },
245           { interaction_function[F_DVDL_VDW].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) },
246           { interaction_function[F_DVDL_BONDED].longname,
247             relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) },
248           { interaction_function[F_DVDL_RESTRAINT].longname,
249             relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) } }
250     };
251
252     // Specify how trajectory frame matching must work
253     TrajectoryComparison trajectoryComparison{ MdrunRerunTest::trajectoryMatchSettings,
254                                                TrajectoryComparison::s_defaultTrajectoryTolerances };
255
256     // The md integrator triggers a warning for nearly decoupled
257     // states, which we need to suppress. TODO sometimes?
258     int numWarningsToTolerate = (integrator == "md") ? 1 : 0;
259     executeRerunTest(&fileManager_,
260                      &runner_,
261                      simulationName,
262                      numWarningsToTolerate,
263                      mdpFieldValues,
264                      energyTermsToCompare,
265                      trajectoryComparison);
266 }
267
268 // TODO The time for OpenCL kernel compilation means these tests time
269 // out. Once that compilation is cached for the whole process, these
270 // tests can run in such configurations.
271 #if !GMX_GPU_OPENCL
272 INSTANTIATE_TEST_SUITE_P(MdrunIsReproduced,
273                          MdrunRerunFreeEnergyTest,
274                          ::testing::Combine(::testing::Values("nonanol_vacuo"),
275                                             ::testing::Values("md", "md-vv", "sd"),
276                                             ::testing::Range(0, 11)));
277 #else
278 INSTANTIATE_TEST_SUITE_P(DISABLED_MdrunIsReproduced,
279                          MdrunRerunFreeEnergyTest,
280                          ::testing::Combine(::testing::Values("nonanol_vacuo"),
281                                             ::testing::Values("md", "md-vv", "sd"),
282                                             ::testing::Range(0, 11)));
283 #endif
284
285 } // namespace
286 } // namespace test
287 } // namespace gmx