Added new testing code for normal mode analysis.
[gromacs.git] / src / programs / mdrun / tests / normalmodes.cpp
blob4535c62334e6310aa572412afe59cfa6ace9fa58
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 for the normal modes functionality.
40 * \author David van der Spoel <david.vanderspoel@icm.uu.se>
41 * \ingroup module_mdrun_integration_tests
43 #include "gmxpre.h"
45 #include <map>
46 #include <memory>
47 #include <string>
48 #include <tuple>
49 #include <vector>
51 #include <gtest/gtest.h>
53 #include "gromacs/compat/make_unique.h"
54 #include "gromacs/options/filenameoption.h"
55 #include "gromacs/topology/idef.h"
56 #include "gromacs/topology/ifunc.h"
57 #include "gromacs/trajectory/energyframe.h"
58 #include "gromacs/trajectory/trajectoryframe.h"
59 #include "gromacs/utility/basenetwork.h"
60 #include "gromacs/utility/filestream.h"
61 #include "gromacs/utility/stringutil.h"
63 #include "testutils/mpitest.h"
64 #include "testutils/refdata.h"
65 #include "testutils/testasserts.h"
66 #include "testutils/xvgtest.h"
68 #include "energycomparison.h"
69 #include "energyreader.h"
70 #include "moduletest.h"
71 #include "simulationdatabase.h"
73 namespace gmx
75 namespace test
77 namespace
80 //! Helper type
81 using MdpField = MdpFieldValues::value_type;
83 /*! \brief Test fixture base for normal mode analysis
85 * This test ensures mdrun can run a normal mode analys, reaching
86 * a reproducible eigenvalues following diagonalization.
88 * The choices for tolerance are arbitrary but sufficient. */
89 class NormalModesTest : public MdrunTestFixture,
90 public ::testing::WithParamInterface <
91 std::tuple < std::string, std::string>>
95 TEST_P(NormalModesTest, WithinTolerances)
97 auto params = GetParam();
98 auto simulationName = std::get<0>(params);
99 auto integrator = std::get<1>(params);
100 SCOPED_TRACE(formatString("Comparing normal modes for '%s'",
101 simulationName.c_str()));
103 // TODO At some point we should also test PME-only ranks.
104 int numRanksAvailable = getNumberOfTestMpiRanks();
105 if (!isNumberOfPpRanksSupported(simulationName, numRanksAvailable))
107 fprintf(stdout, "Test system '%s' cannot run with %d ranks.\n"
108 "The supported numbers are: %s\n",
109 simulationName.c_str(), numRanksAvailable,
110 reportNumbersOfPpRanksSupported(simulationName).c_str());
111 return;
113 auto mdpFieldValues = prepareMdpFieldValues(simulationName.c_str(),
114 integrator.c_str(),
115 "no", "no");
116 mdpFieldValues["nsteps"] = "1";
117 mdpFieldValues["rcoulomb"] = "5.6";
118 mdpFieldValues["rlist"] = "5.6";
119 mdpFieldValues["rvdw"] = "5.6";
120 mdpFieldValues["constraints"] = "none";
121 mdpFieldValues.insert(MdpField("coulombtype", "Cut-off"));
122 mdpFieldValues.insert(MdpField("vdwtype", "Cut-off"));
124 // prepare the .tpr file
126 CommandLine caller;
127 caller.append("grompp");
128 runner_.useTopG96AndNdxFromDatabase(simulationName);
129 runner_.useStringAsMdpFile(prepareMdpFileContents(mdpFieldValues));
130 EXPECT_EQ(0, runner_.callGrompp(caller));
132 // Do mdrun, preparing to check the normal modes later
134 CommandLine mdrunCaller;
135 mdrunCaller.append("mdrun");
136 ASSERT_EQ(0, runner_.callMdrun(mdrunCaller));
138 // Now run gmx nmeig and check the output
140 ASSERT_EQ(0, runner_.callNmeig());
141 TestReferenceData refData;
142 auto checker = refData.rootChecker()
143 .checkCompound("System", simulationName)
144 .checkCompound("Integrator", integrator);
145 auto settings = XvgMatchSettings();
146 TextInputFile input("eigenval.xvg");
147 checkXvgFile(&input, &checker, settings);
151 //! Containers of systems and integrators to test.
152 //! \{
153 std::vector<std::string> systemsToTest_g = { "scaled-water" };
154 std::vector<std::string> integratorsToTest_g = { "nm" };
156 //! \}
158 // The time for OpenCL kernel compilation means these tests might time
159 // out. If that proves to be a problem, these can be disabled for
160 // OpenCL builds. However, once that compilation is cached for the
161 // lifetime of the whole test binary process, these tests should run in
162 // such configurations.
163 #if GMX_DOUBLE
164 INSTANTIATE_TEST_CASE_P(NormalModesWorks, NormalModesTest,
165 ::testing::Combine(::testing::ValuesIn(systemsToTest_g),
166 ::testing::ValuesIn(integratorsToTest_g)));
167 #endif
168 } // namespace
169 } // namespace test
170 } // namespace gmx