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