Apply re-formatting to C++ in src/ tree.
[alexxy/gromacs.git] / src / programs / mdrun / tests / simple_mdrun.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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  * Simple tests for the mdrun functionality.
39  *
40  * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41  * \ingroup module_mdrun_integration_tests
42  */
43 #include "gmxpre.h"
44
45 #include <map>
46 #include <memory>
47 #include <string>
48 #include <tuple>
49 #include <unordered_map>
50 #include <vector>
51
52 #include <gtest/gtest.h>
53
54 #include "gromacs/options/filenameoption.h"
55 #include "gromacs/topology/idef.h"
56 #include "gromacs/topology/ifunc.h"
57 #include "gromacs/trajectory/trajectoryframe.h"
58 #include "gromacs/utility/basenetwork.h"
59 #include "gromacs/utility/filestream.h"
60 #include "gromacs/utility/strconvert.h"
61 #include "gromacs/utility/stringutil.h"
62
63 #include "testutils/mpitest.h"
64 #include "testutils/refdata.h"
65 #include "testutils/simulationdatabase.h"
66 #include "testutils/testasserts.h"
67 #include "testutils/xvgtest.h"
68
69 #include "energycomparison.h"
70 #include "moduletest.h"
71 #include "trajectoryreader.h"
72
73 namespace gmx
74 {
75 namespace test
76 {
77 namespace
78 {
79
80 /*! \brief Database of enerngy tolerances for MD integrator on the various systems. */
81 std::unordered_map<std::string, FloatingPointTolerance> energyToleranceForSystem_g = {
82     { { "angles1", relativeToleranceAsFloatingPoint(1, 1e-4) } }
83 };
84
85 /*! \brief Database of pressure
86    tolerances for MD integrator on the various systems. */
87 std::unordered_map<std::string, FloatingPointTolerance> pressureToleranceForSystem_g = {
88     { { "angles1", relativeToleranceAsFloatingPoint(1, 1e-4) } }
89 };
90
91 //! Helper type
92 using MdpField = MdpFieldValues::value_type;
93
94 /*! \brief Test fixture base for simple mdrun systems
95  *
96  * This test ensures mdrun can run a simulation, reaching
97  * reproducible energies.
98  *
99  * The choices for tolerance are arbitrary but sufficient. */
100 class SimpleMdrunTest :
101     public MdrunTestFixture,
102     public ::testing::WithParamInterface<std::tuple<std::string, std::string>>
103 {
104 };
105
106 TEST_P(SimpleMdrunTest, WithinTolerances)
107 {
108     auto params         = GetParam();
109     auto simulationName = std::get<0>(params);
110     auto integrator     = std::get<1>(params);
111     SCOPED_TRACE(formatString("Comparing simple mdrun for '%s'", simulationName.c_str()));
112
113     // TODO At some point we should also test PME-only ranks.
114     int numRanksAvailable = getNumberOfTestMpiRanks();
115     if (!isNumberOfPpRanksSupported(simulationName, numRanksAvailable))
116     {
117         fprintf(stdout,
118                 "Test system '%s' cannot run with %d ranks.\n"
119                 "The supported numbers are: %s\n",
120                 simulationName.c_str(),
121                 numRanksAvailable,
122                 reportNumbersOfPpRanksSupported(simulationName).c_str());
123         return;
124     }
125     auto mdpFieldValues =
126             prepareMdpFieldValues(simulationName.c_str(), integrator.c_str(), "no", "no");
127     mdpFieldValues["nsteps"]        = "50";
128     mdpFieldValues["nstfout"]       = "4";
129     mdpFieldValues["constraints"]   = "none";
130     mdpFieldValues["nstcalcenergy"] = "4";
131     mdpFieldValues["coulombtype"]   = "Cut-off";
132     mdpFieldValues["vdwtype"]       = "Cut-off";
133
134     // Prepare the .tpr file
135     {
136         CommandLine caller;
137         runner_.useTopGroAndNdxFromDatabase(simulationName);
138         runner_.useStringAsMdpFile(prepareMdpFileContents(mdpFieldValues));
139         EXPECT_EQ(0, runner_.callGrompp(caller));
140     }
141     // Do mdrun
142     {
143         CommandLine mdrunCaller;
144         ASSERT_EQ(0, runner_.callMdrun(mdrunCaller));
145         EnergyTermsToCompare energyTermsToCompare{ {
146                 { interaction_function[F_EPOT].longname, energyToleranceForSystem_g.at(simulationName) },
147                 { interaction_function[F_EKIN].longname, energyToleranceForSystem_g.at(simulationName) },
148                 { interaction_function[F_PRES].longname, pressureToleranceForSystem_g.at(simulationName) },
149         } };
150         TestReferenceData    refData;
151         auto                 checker = refData.rootChecker()
152                                .checkCompound("Simulation", simulationName)
153                                .checkCompound("Mdrun", integrator);
154         checkEnergiesAgainstReferenceData(runner_.edrFileName_, energyTermsToCompare, &checker);
155         // Now check the forces
156         TrajectoryFrameReader reader(runner_.fullPrecisionTrajectoryFileName_);
157         checker.setDefaultTolerance(relativeToleranceAsFloatingPoint(1, 1e-4));
158         do
159         {
160             auto frame = reader.frame();
161             auto force = frame.f();
162             int  atom  = 0;
163             for (auto& f : force)
164             {
165                 std::string forceName = frame.frameName() + " F[" + toString(atom) + "]";
166
167                 checker.checkVector(f, forceName.c_str());
168                 atom++;
169             }
170         } while (reader.readNextFrame());
171     }
172 }
173
174 // The time for OpenCL kernel compilation means these tests might time
175 // out. If that proves to be a problem, these can be disabled for
176 // OpenCL builds. However, once that compilation is cached for the
177 // lifetime of the whole test binary process, these tests should run in
178 // such configurations.
179 #if GMX_DOUBLE
180
181 //! Containers of systems to test.
182 //! \{
183 std::vector<std::string> systemsToTest_g = { "angles1" };
184 std::vector<std::string> md_g            = { "md", "md-vv" };
185 //! \}
186
187 INSTANTIATE_TEST_CASE_P(Angles1,
188                         SimpleMdrunTest,
189                         ::testing::Combine(::testing::ValuesIn(systemsToTest_g), ::testing::ValuesIn(md_g)));
190 #endif
191 } // namespace
192 } // namespace test
193 } // namespace gmx