Fix clang-tidy warnings for clang 7
[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,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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  * \brief
38  * Tests for the mdrun -rerun functionality
39  *
40  * \author Mark Abraham <mark.j.abraham@gmail.com>
41  * \ingroup module_mdrun_integration_tests
42  */
43 #include "gmxpre.h"
44
45 #include "config.h"
46
47 #include <map>
48 #include <memory>
49 #include <string>
50 #include <tuple>
51 #include <vector>
52
53 #include <gtest/gtest.h>
54
55 #include "gromacs/compat/make_unique.h"
56 #include "gromacs/options/filenameoption.h"
57 #include "gromacs/topology/idef.h"
58 #include "gromacs/topology/ifunc.h"
59 #include "gromacs/trajectory/energyframe.h"
60 #include "gromacs/trajectory/trajectoryframe.h"
61 #include "gromacs/utility/stringutil.h"
62
63 #include "testutils/testasserts.h"
64
65 #include "energycomparison.h"
66 #include "energyreader.h"
67 #include "mdruncomparison.h"
68 #include "moduletest.h"
69 #include "simulationdatabase.h"
70 #include "trajectorycomparison.h"
71 #include "trajectoryreader.h"
72
73 namespace gmx
74 {
75 namespace test
76 {
77 namespace
78 {
79
80 //! Functor for comparing reference and test frames on particular energies to match.
81 class EnergyComparator
82 {
83     public:
84         //! Constructor
85         EnergyComparator(const EnergyTolerances &energiesToMatch)
86             : energiesToMatch_(energiesToMatch) {}
87         //! The functor method.
88         void operator()(const EnergyFrame &reference, const EnergyFrame &test)
89         {
90             compareEnergyFrames(reference, test, energiesToMatch_);
91         }
92         //! Container of the energies to match and the tolerance required.
93         const EnergyTolerances &energiesToMatch_;
94 };
95
96 //! Run grompp for a normal mdrun, the run, and its rerun.
97 void executeRerunTest(TestFileManager        *fileManager,
98                       SimulationRunner       *runner,
99                       const std::string      &simulationName,
100                       int                     maxWarningsTolerated,
101                       const MdpFieldValues   &mdpFieldValues,
102                       const EnergyTolerances &energiesToMatch)
103 {
104     auto normalRunTrajectoryFileName = fileManager->getTemporaryFilePath("normal.trr");
105     auto normalRunEdrFileName        = fileManager->getTemporaryFilePath("normal.edr");
106     auto rerunTrajectoryFileName     = fileManager->getTemporaryFilePath("rerun.trr");
107     auto rerunEdrFileName            = fileManager->getTemporaryFilePath("rerun.edr");
108
109     // prepare the .tpr file
110     {
111         // TODO evolve grompp to report the number of warnings issued, so
112         // tests always expect the right number.
113         CommandLine caller;
114         caller.append("grompp");
115         caller.addOption("-maxwarn", maxWarningsTolerated);
116         runner->useTopGroAndNdxFromDatabase(simulationName);
117         runner->useStringAsMdpFile(prepareMdpFileContents(mdpFieldValues));
118         EXPECT_EQ(0, runner->callGrompp(caller));
119     }
120
121     // do the normal mdrun
122     {
123         runner->fullPrecisionTrajectoryFileName_ = normalRunTrajectoryFileName;
124         runner->edrFileName_                     = normalRunEdrFileName;
125         CommandLine normalRunCaller;
126         normalRunCaller.append("mdrun");
127         ASSERT_EQ(0, runner->callMdrun(normalRunCaller));
128     }
129
130     // do a rerun on the .trr just produced
131     {
132         runner->fullPrecisionTrajectoryFileName_ = rerunTrajectoryFileName;
133         runner->edrFileName_                     = rerunEdrFileName;
134         CommandLine rerunCaller;
135         rerunCaller.append("mdrun");
136         rerunCaller.addOption("-rerun", normalRunTrajectoryFileName);
137         ASSERT_EQ(0, runner->callMdrun(rerunCaller));
138     }
139
140     // Build the functor that will compare reference and test
141     // energy frames on the chosen energy fields.
142     //
143     // TODO It would be less code if we used a lambda for this, but either
144     // clang 3.4 or libstdc++ 5.2.1 have an issue with capturing a
145     // std::unordered_map
146     EnergyComparator energyComparator(energiesToMatch);
147     // Build the manager that will present matching pairs of frames to compare.
148     //
149     // TODO Here is an unnecessary copy of keys (ie. the energy field
150     // names), for convenience. In the future, use a range.
151     auto namesOfEnergiesToMatch = getKeys(energiesToMatch);
152     FramePairManager<EnergyFrameReader, EnergyFrame>
153          energyManager(openEnergyFileToReadFields(normalRunEdrFileName, namesOfEnergiesToMatch),
154                   openEnergyFileToReadFields(rerunEdrFileName, namesOfEnergiesToMatch));
155     // Compare the energy frames.
156     energyManager.compareAllFramePairs(energyComparator);
157
158     // Specify how trajectory frame matching must work.
159     TrajectoryFrameMatchSettings trajectoryMatchSettings {
160         true, true, true, true, true, true
161     };
162     /* Specify the default expected tolerances for trajectory
163      * components for all simulation systems. */
164     TrajectoryTolerances trajectoryTolerances {
165         defaultRealTolerance(),                                               // box
166         relativeToleranceAsFloatingPoint(1.0, 1.0e-3),                        // positions
167         defaultRealTolerance(),                                               // velocities - unused for rerun
168         relativeToleranceAsFloatingPoint(100.0, GMX_DOUBLE ? 1.0e-7 : 1.0e-5) // forces
169     };
170
171     // Build the functor that will compare reference and test
172     // trajectory frames in the chosen way.
173     auto trajectoryComparator = [&trajectoryMatchSettings, &trajectoryTolerances](const TrajectoryFrame &reference, const TrajectoryFrame &test)
174         {
175             compareTrajectoryFrames(reference, test, trajectoryMatchSettings, trajectoryTolerances);
176         };
177     // Build the manager that will present matching pairs of frames to compare
178     FramePairManager<TrajectoryFrameReader, TrajectoryFrame>
179     trajectoryManager(compat::make_unique<TrajectoryFrameReader>(normalRunTrajectoryFileName),
180                       compat::make_unique<TrajectoryFrameReader>(rerunTrajectoryFileName));
181     // Compare the trajectory frames.
182     trajectoryManager.compareAllFramePairs(trajectoryComparator);
183 }
184
185 /*! \brief Test fixture base for mdrun -rerun
186  *
187  * This test ensures mdrun can run a simulation, writing a trajectory
188  * and matching energies, and reproduce the same energies from a rerun
189  * to within a tight tolerance. It says nothing about whether a rerun
190  * can reproduce energies from a trajectory generated with older code,
191  * since that is not a useful property. Whether mdrun produced correct
192  * energies then and now needs different kinds of testing, but if
193  * true, this test ensures the rerun has the expected property.
194  *
195  * The limitations of mdrun and its output means that reproducing the
196  * same energies is currently only meaningful for integration without
197  * thermostats or barostats, however the present form of the test
198  * infrastructure has in-principle support for such, if that is ever
199  * needed/useful.
200  *
201  * We should also not compare pressure, because with constraints the
202  * non-search steps need a much larger tolerance, and per Redmine 1868
203  * we should stop computing pressure in reruns anyway.
204  *
205  * Similarly, per 1868, in the present implementation the kinetic
206  * energy quantities are not generally reproducible, either.
207  *
208  * The choices for tolerance are arbitrary but sufficient.  Rerun does
209  * pair search every frame, so it cannot in general exactly reproduce
210  * quantities from a normal run, because the accumulation order
211  * differs. (Nor does it reproduce pair-search frames exactly,
212  * either). */
213 class MdrunRerunTest : public MdrunTestFixture,
214                        public ::testing::WithParamInterface <
215                        std::tuple < std::string, std::string>>
216 {
217 };
218
219 TEST_P(MdrunRerunTest, WithinTolerances)
220 {
221     auto params         = GetParam();
222     auto simulationName = std::get<0>(params);
223     auto integrator     = std::get<1>(params);
224     SCOPED_TRACE(formatString("Comparing normal and rerun of simulation '%s' "
225                               "with integrator '%s'",
226                               simulationName.c_str(), integrator.c_str()));
227
228     auto mdpFieldValues = prepareMdpFieldValues(simulationName.c_str(),
229                                                 integrator.c_str(),
230                                                 "no", "no");
231
232     // bd is much less reproducible in a rerun than the other integrators
233     const int        toleranceScaleFactor = (integrator == "bd") ? 2 : 1;
234     EnergyTolerances energiesToMatch
235     {{
236          {
237              interaction_function[F_EPOT].longname,
238              relativeToleranceAsPrecisionDependentUlp(10.0, 24 * toleranceScaleFactor, 40 * toleranceScaleFactor)
239          },
240      }};
241
242     int numWarningsToTolerate = 0;
243     executeRerunTest(&fileManager_, &runner_,
244                      simulationName, numWarningsToTolerate, mdpFieldValues,
245                      energiesToMatch);
246 }
247
248 // TODO The time for OpenCL kernel compilation means these tests time
249 // out. Once that compilation is cached for the whole process, these
250 // tests can run in such configurations.
251 #if GMX_GPU != GMX_GPU_OPENCL
252 INSTANTIATE_TEST_CASE_P(NormalMdrunIsReproduced, MdrunRerunTest,
253                             ::testing::Combine(::testing::Values("argon12", "spc5", "alanine_vsite_vacuo"),
254                                                    ::testing::Values("md", "md-vv", "bd", "sd")));
255 #else
256 INSTANTIATE_TEST_CASE_P(DISABLED_NormalMdrunIsReproduced, MdrunRerunTest,
257                             ::testing::Combine(::testing::Values("argon12", "spc5", "alanine_vsite_vacuo"),
258                                                    ::testing::Values("md", "md-vv", "bd", "sd")));
259 #endif
260
261 class MdrunRerunFreeEnergyTest : public MdrunTestFixture,
262                                  public ::testing::WithParamInterface <
263                                  std::tuple < std::string, std::string, int>>
264 {
265 };
266
267 TEST_P(MdrunRerunFreeEnergyTest, WithinTolerances)
268 {
269     auto params          = GetParam();
270     auto simulationName  = std::get<0>(params);
271     auto integrator      = std::get<1>(params);
272     auto initLambdaState = std::get<2>(params);
273     SCOPED_TRACE(formatString("Comparing normal and rerun of simulation '%s' "
274                               "with integrator '%s' for initial lambda state %d",
275                               simulationName.c_str(), integrator.c_str(), initLambdaState));
276
277     auto mdpFieldValues = prepareMdpFieldValues(simulationName.c_str(),
278                                                 integrator.c_str(),
279                                                 "no", "no");
280     mdpFieldValues["other"] += formatString("\ninit-lambda-state %d", initLambdaState);
281
282     EnergyTolerances energiesToMatch
283     {{
284          {
285              interaction_function[F_EPOT].longname, relativeToleranceAsPrecisionDependentUlp(10.0, 24, 32)
286          },
287          {
288              interaction_function[F_DVDL_COUL].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8)
289          },
290          {
291              interaction_function[F_DVDL_VDW].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8)
292          },
293          {
294              interaction_function[F_DVDL_BONDED].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8)
295          },
296          {
297              interaction_function[F_DVDL_RESTRAINT].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8)
298          }
299      }};
300
301     // The md integrator triggers a warning for nearly decoupled
302     // states, which we need to suppress. TODO sometimes?
303     int numWarningsToTolerate = (integrator == "md") ? 1 : 0;
304     executeRerunTest(&fileManager_, &runner_,
305                      simulationName, numWarningsToTolerate, mdpFieldValues,
306                      energiesToMatch);
307 }
308
309 // TODO The time for OpenCL kernel compilation means these tests time
310 // out. Once that compilation is cached for the whole process, these
311 // tests can run in such configurations.
312 #if GMX_GPU != GMX_GPU_OPENCL
313 INSTANTIATE_TEST_CASE_P(MdrunIsReproduced, MdrunRerunFreeEnergyTest,
314                             ::testing::Combine(::testing::Values("nonanol_vacuo"),
315                                                    ::testing::Values("md", "md-vv", "sd"),
316                                                    ::testing::Range(0, 11)));
317 #else
318 INSTANTIATE_TEST_CASE_P(DISABLED_MdrunIsReproduced, MdrunRerunFreeEnergyTest,
319                             ::testing::Combine(::testing::Values("nonanol_vacuo"),
320                                                    ::testing::Values("md", "md-vv", "sd"),
321                                                    ::testing::Range(0, 11)));
322 #endif
323
324 } // namespace
325 } // namespace test
326 } // namespace gmx