Port basic pull test to gtest framework
[alexxy/gromacs.git] / src / programs / mdrun / tests / moduletest.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief
38  * Implements classes in moduletest.h.
39  *
40  * \author Mark Abraham <mark.j.abraham@gmail.com>
41  * \ingroup module_mdrun_integration_tests
42  */
43 #include "gmxpre.h"
44
45 #include "moduletest.h"
46
47 #include "config.h"
48
49 #include <cstdio>
50
51 #include <utility>
52
53 #include "gromacs/gmxana/gmx_ana.h"
54 #include "gromacs/gmxpreprocess/grompp.h"
55 #include "gromacs/hardware/detecthardware.h"
56 #include "gromacs/hardware/hw_info.h"
57 #include "gromacs/options/basicoptions.h"
58 #include "gromacs/options/ioptionscontainer.h"
59 #include "gromacs/tools/convert_tpr.h"
60 #include "gromacs/utility/basedefinitions.h"
61 #include "gromacs/utility/basenetwork.h"
62 #include "gromacs/utility/gmxmpi.h"
63 #include "gromacs/utility/physicalnodecommunicator.h"
64 #include "gromacs/utility/textwriter.h"
65 #include "programs/mdrun/mdrun_main.h"
66
67 #include "testutils/cmdlinetest.h"
68 #include "testutils/mpitest.h"
69 #include "testutils/testfilemanager.h"
70 #include "testutils/testoptions.h"
71
72 namespace gmx
73 {
74 namespace test
75 {
76
77 /********************************************************************
78  * MdrunTestFixture
79  */
80
81 namespace
82 {
83
84 #if GMX_OPENMP || defined(DOXYGEN)
85 //! Number of OpenMP threads for child mdrun call.
86 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
87 int g_numOpenMPThreads = 1;
88 #endif
89 //! \cond
90 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
91 GMX_TEST_OPTIONS(MdrunTestOptions, options)
92 {
93     GMX_UNUSED_VALUE(options);
94 #if GMX_OPENMP
95     options->addOption(
96             IntegerOption("ntomp").store(&g_numOpenMPThreads).description("Number of OpenMP threads for child mdrun calls"));
97 #endif
98 }
99 //! \endcond
100
101 } // namespace
102
103 SimulationRunner::SimulationRunner(TestFileManager* fileManager) :
104     fullPrecisionTrajectoryFileName_(fileManager->getTemporaryFilePath(".trr")),
105     mdpOutputFileName_(fileManager->getTemporaryFilePath("output.mdp")),
106     tprFileName_(fileManager->getTemporaryFilePath(".tpr")),
107     logFileName_(fileManager->getTemporaryFilePath(".log")),
108     edrFileName_(fileManager->getTemporaryFilePath(".edr")),
109     mtxFileName_(fileManager->getTemporaryFilePath(".mtx")),
110
111     nsteps_(-2),
112     mdpSource_(SimulationRunnerMdpSource::Undefined),
113     fileManager_(*fileManager)
114 {
115 #if GMX_LIB_MPI
116     GMX_RELEASE_ASSERT(gmx_mpi_initialized(), "MPI system not initialized for mdrun tests");
117
118     // It would be better to also detect this in a thread-MPI build,
119     // but there is no way to do that currently, and it is also not a
120     // problem for such a build. Any code based on such an invalid
121     // test fixture will be found in CI testing, however.
122     GMX_RELEASE_ASSERT(MdrunTestFixtureBase::communicator_ != MPI_COMM_NULL,
123                        "SimulationRunner may only be used from a test fixture that inherits from "
124                        "MdrunTestFixtureBase");
125 #endif
126 }
127
128 // TODO The combination of defaulting to Verlet cut-off scheme, NVE,
129 // and verlet-buffer-tolerance = -1 gives a grompp error. If we keep
130 // things that way, this function should be renamed. For now,
131 // we use the Verlet scheme and hard-code a tolerance.
132 // TODO There is possible outstanding unexplained behaviour of mdp
133 // input parsing e.g. Issue #2074, so this particular set of mdp
134 // contents is also tested with GetIrTest in gmxpreprocess-test.
135 void SimulationRunner::useEmptyMdpFile()
136 {
137     useStringAsMdpFile("");
138 }
139
140 void SimulationRunner::useStringAsMdpFile(const char* mdpString)
141 {
142     useStringAsMdpFile(std::string(mdpString));
143 }
144
145 void SimulationRunner::useStringAsMdpFile(const std::string& mdpString)
146 {
147     GMX_RELEASE_ASSERT(mdpSource_ != SimulationRunnerMdpSource::File,
148                        "Cannot mix .mdp file from database with options set via string.");
149     mdpSource_        = SimulationRunnerMdpSource::String;
150     mdpInputContents_ = mdpString;
151 }
152
153 void SimulationRunner::useStringAsNdxFile(const char* ndxString) const
154 {
155     gmx::TextWriter::writeFileFromString(ndxFileName_, ndxString);
156 }
157
158 void SimulationRunner::useTopG96AndNdxFromDatabase(const std::string& name)
159 {
160     topFileName_ = gmx::test::TestFileManager::getInputFilePath(name + ".top");
161     groFileName_ = gmx::test::TestFileManager::getInputFilePath(name + ".g96");
162     ndxFileName_ = gmx::test::TestFileManager::getInputFilePath(name + ".ndx");
163 }
164
165 void SimulationRunner::useTopGroAndNdxFromDatabase(const std::string& name)
166 {
167     topFileName_ = gmx::test::TestFileManager::getInputFilePath(name + ".top");
168     groFileName_ = gmx::test::TestFileManager::getInputFilePath(name + ".gro");
169     ndxFileName_ = gmx::test::TestFileManager::getInputFilePath(name + ".ndx");
170 }
171
172 void SimulationRunner::useGroFromDatabase(const char* name)
173 {
174     groFileName_ = gmx::test::TestFileManager::getInputFilePath((std::string(name) + ".gro").c_str());
175 }
176
177 void SimulationRunner::useNdxFromDatabase(const std::string& name)
178 {
179     ndxFileName_ = gmx::test::TestFileManager::getInputFilePath(name + ".ndx");
180 }
181
182 void SimulationRunner::useTopGroAndMdpFromFepTestDatabase(const std::string& name)
183 {
184     GMX_RELEASE_ASSERT(mdpSource_ != SimulationRunnerMdpSource::String,
185                        "Cannot mix .mdp file from database with options set via string.");
186     mdpSource_   = SimulationRunnerMdpSource::File;
187     topFileName_ = gmx::test::TestFileManager::getInputFilePath("freeenergy/" + name + "/topol.top");
188     groFileName_ = gmx::test::TestFileManager::getInputFilePath("freeenergy/" + name + "/conf.gro");
189     mdpFileName_ =
190             gmx::test::TestFileManager::getInputFilePath("freeenergy/" + name + "/grompp.mdp");
191 }
192
193 int SimulationRunner::callGromppOnThisRank(const CommandLine& callerRef)
194 {
195     std::string mdpInputFileName;
196     if (mdpSource_ == SimulationRunnerMdpSource::File)
197     {
198         mdpInputFileName = mdpFileName_;
199     }
200     else
201     {
202         mdpInputFileName = fileManager_.getTemporaryFilePath("input.mdp");
203         gmx::TextWriter::writeFileFromString(mdpInputFileName, mdpInputContents_);
204     }
205
206     CommandLine caller;
207     caller.append("grompp");
208     caller.merge(callerRef);
209     caller.addOption("-f", mdpInputFileName);
210     if (!ndxFileName_.empty())
211     {
212         caller.addOption("-n", ndxFileName_);
213     }
214     caller.addOption("-p", topFileName_);
215     caller.addOption("-c", groFileName_);
216     caller.addOption("-r", groFileName_);
217
218     caller.addOption("-po", mdpOutputFileName_);
219     caller.addOption("-o", tprFileName_);
220
221     return gmx_grompp(caller.argc(), caller.argv());
222 }
223
224 int SimulationRunner::callGromppOnThisRank()
225 {
226     return callGromppOnThisRank(CommandLine());
227 }
228
229 int SimulationRunner::callGrompp(const CommandLine& callerRef)
230 {
231     int returnValue = 0;
232 #if GMX_LIB_MPI
233     // When compiled with external MPI, we're trying to run mdrun with
234     // MPI, but we need to make sure that we only do grompp on one
235     // rank
236     if (0 == gmx_node_rank())
237 #endif
238     {
239         returnValue = callGromppOnThisRank(callerRef);
240     }
241 #if GMX_LIB_MPI
242     // Make sure rank zero has written the .tpr file before other
243     // ranks try to read it. Thread-MPI and serial do this just fine
244     // on their own.
245     MPI_Barrier(MdrunTestFixtureBase::communicator_);
246 #endif
247     return returnValue;
248 }
249
250 int SimulationRunner::callGrompp()
251 {
252     return callGrompp(CommandLine());
253 }
254
255 int SimulationRunner::changeTprNsteps(int nsteps) const
256 {
257     CommandLine caller;
258     caller.append("convert-tpr");
259     caller.addOption("-nsteps", nsteps);
260     // Because the operation is to change the .tpr, we replace the
261     // file. TODO Do we need to delete an automatic backup?
262     caller.addOption("-s", tprFileName_);
263     caller.addOption("-o", tprFileName_);
264
265     return gmx::test::CommandLineTestHelper::runModuleFactory(&gmx::ConvertTprInfo::create, &caller);
266 }
267
268 int SimulationRunner::callNmeig() const
269 {
270     /* Conforming to style guide by not passing a non-const reference
271        to this function. Passing a non-const reference might make it
272        easier to write code that incorrectly re-uses callerRef after
273        the call to this function. */
274
275     CommandLine caller;
276     caller.append("nmeig");
277     caller.addOption("-s", tprFileName_);
278     caller.addOption("-f", mtxFileName_);
279     // Ignore the overall translation and rotation in the
280     // first six eigenvectors.
281     caller.addOption("-first", "7");
282     // No need to check more than a number of output values.
283     caller.addOption("-last", "50");
284     caller.addOption("-xvg", "none");
285
286     return gmx_nmeig(caller.argc(), caller.argv());
287 }
288
289 int SimulationRunner::callMdrun(const CommandLine& callerRef)
290 {
291     /* Conforming to style guide by not passing a non-const reference
292        to this function. Passing a non-const reference might make it
293        easier to write code that incorrectly re-uses callerRef after
294        the call to this function. */
295
296     CommandLine caller;
297     caller.append("mdrun");
298     caller.merge(callerRef);
299     caller.addOption("-s", tprFileName_);
300
301     caller.addOption("-g", logFileName_);
302     caller.addOption("-e", edrFileName_);
303     caller.addOption("-mtx", mtxFileName_);
304     caller.addOption("-o", fullPrecisionTrajectoryFileName_);
305     caller.addOption("-x", reducedPrecisionTrajectoryFileName_);
306     if (!dhdlFileName_.empty())
307     {
308         caller.addOption("-dhdl", dhdlFileName_);
309     }
310
311     caller.addOption("-deffnm", fileManager_.getTemporaryFilePath("state"));
312
313     if (nsteps_ > -2)
314     {
315         caller.addOption("-nsteps", nsteps_);
316     }
317
318 #if GMX_THREAD_MPI
319     caller.addOption("-ntmpi", getNumberOfTestMpiRanks());
320 #endif
321
322 #if GMX_OPENMP
323     caller.addOption("-ntomp", g_numOpenMPThreads);
324 #endif
325
326     return gmx_mdrun(MdrunTestFixtureBase::communicator_,
327                      *MdrunTestFixtureBase::hwinfo_,
328                      caller.argc(),
329                      caller.argv());
330 }
331
332 int SimulationRunner::callMdrun()
333 {
334     return callMdrun(CommandLine());
335 }
336
337 // ====
338
339 // static
340 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
341 MPI_Comm MdrunTestFixtureBase::communicator_ = MPI_COMM_NULL;
342 // static
343 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
344 std::unique_ptr<gmx_hw_info_t> MdrunTestFixtureBase::hwinfo_;
345
346 // static
347 void MdrunTestFixtureBase::SetUpTestCase()
348 {
349     communicator_ = MPI_COMM_WORLD;
350     auto newHwinfo =
351             gmx_detect_hardware(PhysicalNodeCommunicator{ communicator_, gmx_physicalnode_id_hash() });
352     std::swap(hwinfo_, newHwinfo);
353 }
354
355 // static
356 void MdrunTestFixtureBase::TearDownTestCase()
357 {
358     hwinfo_.reset(nullptr);
359 }
360
361 MdrunTestFixtureBase::MdrunTestFixtureBase()
362 {
363 #if GMX_LIB_MPI
364     GMX_RELEASE_ASSERT(gmx_mpi_initialized(), "MPI system not initialized for mdrun tests");
365 #endif
366 }
367
368 MdrunTestFixtureBase::~MdrunTestFixtureBase() {}
369
370 // ====
371
372 MdrunTestFixture::MdrunTestFixture() : runner_(&fileManager_) {}
373
374 MdrunTestFixture::~MdrunTestFixture()
375 {
376 #if GMX_LIB_MPI
377     // fileManager_ should only clean up after all the ranks are done.
378     MPI_Barrier(MdrunTestFixtureBase::communicator_);
379 #endif
380 }
381
382 } // namespace test
383 } // namespace gmx