Apply clang-format to source tree
[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,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.
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 "gromacs/topology/ifunc.h"
48 #include "gromacs/trajectory/energyframe.h"
49 #include "gromacs/utility/stringutil.h"
50
51 #include "simulatorcomparison.h"
52
53 namespace gmx
54 {
55 namespace test
56 {
57 namespace
58 {
59 /*! \brief Test fixture base for mdrun -rerun
60  *
61  * This test ensures mdrun can run a simulation, writing a trajectory
62  * and matching energies, and reproduce the same energies from a rerun
63  * to within a tight tolerance. It says nothing about whether a rerun
64  * can reproduce energies from a trajectory generated with older code,
65  * since that is not a useful property. Whether mdrun produced correct
66  * energies then and now needs different kinds of testing, but if
67  * true, this test ensures the rerun has the expected property.
68  *
69  * The limitations of mdrun and its output means that reproducing the
70  * same energies is currently only meaningful for integration without
71  * thermostats or barostats, however the present form of the test
72  * infrastructure has in-principle support for such, if that is ever
73  * needed/useful.
74  *
75  * We should also not compare pressure, because with constraints the
76  * non-search steps need a much larger tolerance, and per Redmine 1868
77  * we should stop computing pressure in reruns anyway.
78  *
79  * Similarly, per 1868, in the present implementation the kinetic
80  * energy quantities are not generally reproducible, either.
81  *
82  * The choices for tolerance are arbitrary but sufficient.  Rerun does
83  * pair search every frame, so it cannot in general exactly reproduce
84  * quantities from a normal run, because the accumulation order
85  * differs. (Nor does it reproduce pair-search frames exactly,
86  * either). */
87 class MdrunRerunTest :
88     public MdrunTestFixture,
89     public ::testing::WithParamInterface<std::tuple<std::string, std::string>>
90 {
91 public:
92     //! Trajectory components to compare
93     static const TrajectoryFrameMatchSettings trajectoryMatchSettings;
94 };
95
96 // Compare box, positions and forces, but not velocities
97 // (velocities are ignored in reruns)
98 const TrajectoryFrameMatchSettings MdrunRerunTest::trajectoryMatchSettings = {
99     true,
100     true,
101     true,
102     ComparisonConditions::MustCompare,
103     ComparisonConditions::NoComparison,
104     ComparisonConditions::MustCompare
105 };
106
107 TEST_P(MdrunRerunTest, WithinTolerances)
108 {
109     auto params         = GetParam();
110     auto simulationName = std::get<0>(params);
111     auto integrator     = std::get<1>(params);
112     SCOPED_TRACE(
113             formatString("Comparing normal and rerun of simulation '%s' "
114                          "with integrator '%s'",
115                          simulationName.c_str(), integrator.c_str()));
116
117     auto mdpFieldValues =
118             prepareMdpFieldValues(simulationName.c_str(), integrator.c_str(), "no", "no");
119
120     // bd is much less reproducible in a rerun than the other integrators
121     const int            toleranceScaleFactor = (integrator == "bd") ? 2 : 1;
122     EnergyTermsToCompare energyTermsToCompare{ {
123             { interaction_function[F_EPOT].longname,
124               relativeToleranceAsPrecisionDependentUlp(10.0, 24 * toleranceScaleFactor,
125                                                        40 * toleranceScaleFactor) },
126     } };
127
128     // Specify how trajectory frame matching must work
129     TrajectoryComparison trajectoryComparison{ trajectoryMatchSettings,
130                                                TrajectoryComparison::s_defaultTrajectoryTolerances };
131
132     int numWarningsToTolerate = 0;
133     executeRerunTest(&fileManager_, &runner_, simulationName, numWarningsToTolerate, mdpFieldValues,
134                      energyTermsToCompare, trajectoryComparison);
135 }
136
137 // TODO The time for OpenCL kernel compilation means these tests time
138 // out. Once that compilation is cached for the whole process, these
139 // tests can run in such configurations.
140 #if GMX_GPU != GMX_GPU_OPENCL
141 INSTANTIATE_TEST_CASE_P(
142         NormalMdrunIsReproduced,
143         MdrunRerunTest,
144         ::testing::Combine(::testing::Values("argon12", "tip3p5", "alanine_vsite_vacuo"),
145                            ::testing::Values("md", "md-vv", "bd", "sd")));
146 #else
147 INSTANTIATE_TEST_CASE_P(
148         DISABLED_NormalMdrunIsReproduced,
149         MdrunRerunTest,
150         ::testing::Combine(::testing::Values("argon12", "tip3p5", "alanine_vsite_vacuo"),
151                            ::testing::Values("md", "md-vv", "bd", "sd")));
152 #endif
153
154 class MdrunRerunFreeEnergyTest :
155     public MdrunTestFixture,
156     public ::testing::WithParamInterface<std::tuple<std::string, std::string, int>>
157 {
158 };
159
160 TEST_P(MdrunRerunFreeEnergyTest, WithinTolerances)
161 {
162     auto params          = GetParam();
163     auto simulationName  = std::get<0>(params);
164     auto integrator      = std::get<1>(params);
165     auto initLambdaState = std::get<2>(params);
166     SCOPED_TRACE(
167             formatString("Comparing normal and rerun of simulation '%s' "
168                          "with integrator '%s' for initial lambda state %d",
169                          simulationName.c_str(), integrator.c_str(), initLambdaState));
170
171     auto mdpFieldValues =
172             prepareMdpFieldValues(simulationName.c_str(), integrator.c_str(), "no", "no");
173     mdpFieldValues["other"] += formatString("\ninit-lambda-state = %d", initLambdaState);
174
175     EnergyTermsToCompare energyTermsToCompare{
176         { { interaction_function[F_EPOT].longname, relativeToleranceAsPrecisionDependentUlp(10.0, 24, 32) },
177           { interaction_function[F_DVDL_COUL].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) },
178           { interaction_function[F_DVDL_VDW].longname, relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) },
179           { interaction_function[F_DVDL_BONDED].longname,
180             relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) },
181           { interaction_function[F_DVDL_RESTRAINT].longname,
182             relativeToleranceAsPrecisionDependentUlp(1.0, 8, 8) } }
183     };
184
185     // Specify how trajectory frame matching must work
186     TrajectoryComparison trajectoryComparison{ MdrunRerunTest::trajectoryMatchSettings,
187                                                TrajectoryComparison::s_defaultTrajectoryTolerances };
188
189     // The md integrator triggers a warning for nearly decoupled
190     // states, which we need to suppress. TODO sometimes?
191     int numWarningsToTolerate = (integrator == "md") ? 1 : 0;
192     executeRerunTest(&fileManager_, &runner_, simulationName, numWarningsToTolerate, mdpFieldValues,
193                      energyTermsToCompare, trajectoryComparison);
194 }
195
196 // TODO The time for OpenCL kernel compilation means these tests time
197 // out. Once that compilation is cached for the whole process, these
198 // tests can run in such configurations.
199 #if GMX_GPU != GMX_GPU_OPENCL
200 INSTANTIATE_TEST_CASE_P(MdrunIsReproduced,
201                         MdrunRerunFreeEnergyTest,
202                         ::testing::Combine(::testing::Values("nonanol_vacuo"),
203                                            ::testing::Values("md", "md-vv", "sd"),
204                                            ::testing::Range(0, 11)));
205 #else
206 INSTANTIATE_TEST_CASE_P(DISABLED_MdrunIsReproduced,
207                         MdrunRerunFreeEnergyTest,
208                         ::testing::Combine(::testing::Values("nonanol_vacuo"),
209                                            ::testing::Values("md", "md-vv", "sd"),
210                                            ::testing::Range(0, 11)));
211 #endif
212
213 } // namespace
214 } // namespace test
215 } // namespace gmx