Generalize constraints on MPI rank counts for tests
[alexxy/gromacs.git] / src / programs / mdrun / tests / multisimtest.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,2014,2015,2016,2018 by the GROMACS development team.
5  * Copyright (c) 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
37 /*! \internal \file
38  * \brief
39  * Tests for the mdrun multi-simulation functionality
40  *
41  * \author Mark Abraham <mark.j.abraham@gmail.com>
42  * \ingroup module_mdrun_integration_tests
43  */
44 #include "gmxpre.h"
45
46 #include "multisimtest.h"
47
48 #include <cmath>
49
50 #include <algorithm>
51 #include <string>
52
53 #include <gtest/gtest.h>
54
55 #include "gromacs/mdtypes/md_enums.h"
56 #include "gromacs/utility/basenetwork.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/path.h"
59 #include "gromacs/utility/real.h"
60 #include "gromacs/utility/stringutil.h"
61
62 #include "testutils/cmdlinetest.h"
63 #include "testutils/mpitest.h"
64
65 #include "moduletest.h"
66 #include "terminationhelper.h"
67
68 namespace gmx
69 {
70 namespace test
71 {
72
73 MultiSimTest::MultiSimTest() :
74     size_(gmx_node_num()),
75     rank_(gmx_node_rank()),
76     numRanksPerSimulation_(std::get<0>(GetParam())),
77     simulationNumber_(rank_ / numRanksPerSimulation_),
78     mdrunCaller_(new CommandLine)
79
80 {
81     // Zero or less ranks doesn't make sense
82     GMX_RELEASE_ASSERT(numRanksPerSimulation_ > 0, "Invalid number of ranks per simulation.");
83
84     const char* directoryNameFormat = "sim_%d";
85
86     // Modify the file manager to have a temporary directory unique to
87     // each simulation. No need to have a mutex on this, nobody else
88     // can access the fileManager_ yet because we only just
89     // constructed it.
90     std::string originalTempDirectory = fileManager_.getOutputTempDirectory();
91     std::string newTempDirectory =
92             Path::join(originalTempDirectory, formatString(directoryNameFormat, simulationNumber_));
93     if (rank_ % numRanksPerSimulation_ == 0)
94     {
95         // Only one rank per simulation creates directory
96         Directory::create(newTempDirectory);
97     }
98 #if GMX_LIB_MPI
99     // Make sure directories got created.
100     MPI_Barrier(MdrunTestFixtureBase::communicator_);
101 #endif
102     fileManager_.setOutputTempDirectory(newTempDirectory);
103
104     mdrunCaller_->append("mdrun");
105     mdrunCaller_->addOption("-multidir");
106     for (int i = 0; i < size_ / numRanksPerSimulation_; ++i)
107     {
108         mdrunCaller_->append(Path::join(originalTempDirectory, formatString(directoryNameFormat, i)));
109     }
110 }
111
112 bool MultiSimTest::mpiSetupValid() const
113 {
114     // Single simulation case is not implemented in multi-sim
115     const bool haveAtLeastTwoSimulations = ((size_ / numRanksPerSimulation_) >= 2);
116     // Mdrun will throw error if simulations don't have identical number of ranks
117     const bool simulationsHaveIdenticalRankNumber = ((size_ % numRanksPerSimulation_) == 0);
118
119     return (haveAtLeastTwoSimulations && simulationsHaveIdenticalRankNumber);
120 }
121
122 void MultiSimTest::organizeMdpFile(SimulationRunner*    runner,
123                                    IntegrationAlgorithm integrator,
124                                    TemperatureCoupling  tcoupl,
125                                    PressureCoupling     pcoupl,
126                                    int                  numSteps,
127                                    bool                 doRegression) const
128 {
129     GMX_RELEASE_ASSERT(mpiSetupValid(), "Creating the mdp file without valid MPI setup is useless.");
130     const real  baseTemperature = 298;
131     const real  basePressure    = 1;
132     std::string mdpFileContents = formatString(
133             "integrator = %s\n"
134             "tcoupl = %s\n"
135             "pcoupl = %s\n"
136             "nsteps = %d\n"
137             "nstlog = 1\n"
138             "nstcalcenergy = 1\n"
139             "tc-grps = System\n"
140             "tau-t = 1\n"
141             "ref-t = %f\n"
142             // pressure coupling (if active)
143             "tau-p = 2\n"
144             "ref-p = %f\n"
145             "compressibility = 4.5e-5\n"
146             // velocity generation
147             "gen-vel = yes\n"
148             "gen-temp = %f\n"
149             "gen-seed = %d\n"
150             // v-rescale and c-rescale use ld-seed
151             "ld-seed = %d\n"
152             // Two systems are used: spc2    non-interacting also at cutoff 1nm
153             //                       tip3p5  has box length of 1.86
154             "rcoulomb = 0.7\n"
155             "rvdw = 0.7\n"
156             // Trajectory output if required
157             "nstxout = %d\n"
158             "nstvout = %d\n"
159             "nstfout = %d\n"
160             "nstenergy = %d\n",
161             enumValueToString(integrator),
162             enumValueToString(tcoupl),
163             enumValueToString(pcoupl),
164             numSteps,
165             baseTemperature + 0.0001 * rank_,
166             basePressure * std::pow(1.01, rank_),
167             /* Set things up so that the initial KE decreases with
168                increasing replica number, so that the (identical)
169                starting PE decreases on the first step more for the
170                replicas with higher number, which will tend to force
171                replica exchange to occur. */
172             std::max(baseTemperature - 10 * rank_, real(0)),
173             // If we do regression, we need reproducible velocity
174             // generation, which can be different per simulation
175             (doRegression ? 671324 + simulationNumber_ : -1),
176             // If we do regression, we need reproducible temperature and
177             // pressure coupling, which can be different per simulation
178             (doRegression ? 51203 + simulationNumber_ : -1),
179             // If we do regression, write one intermediate point
180             (doRegression ? int(numSteps / 2) : 0),
181             (doRegression ? int(numSteps / 2) : 0),
182             (doRegression ? int(numSteps / 2) : 0),
183             // If we do regression, print energies every step so
184             // we're sure to catch the replica exchange steps
185             (doRegression ? 1 : 1000));
186     runner->useStringAsMdpFile(mdpFileContents);
187 }
188
189 void MultiSimTest::runGrompp(SimulationRunner* runner, int numSteps, bool doRegression, int maxWarnings) const
190 {
191     // Call grompp once per simulation
192     if (rank_ % numRanksPerSimulation_ == 0)
193     {
194         const auto& simulator = std::get<1>(GetParam());
195         const auto& tcoupl    = std::get<2>(GetParam());
196         const auto& pcoupl    = std::get<3>(GetParam());
197         organizeMdpFile(runner, simulator, tcoupl, pcoupl, numSteps, doRegression);
198         CommandLine caller;
199         caller.addOption("-maxwarn", maxWarnings);
200         EXPECT_EQ(0, runner->callGromppOnThisRank(caller));
201     }
202
203 #if GMX_LIB_MPI
204     // Make sure simulation masters have written the .tpr file before other ranks try to read it.
205     MPI_Barrier(MdrunTestFixtureBase::communicator_);
206 #endif
207 }
208
209 void MultiSimTest::runExitsNormallyTest()
210 {
211     if (!mpiSetupValid())
212     {
213         // Can't test multi-sim without multiple simulations
214         return;
215     }
216
217     SimulationRunner runner(&fileManager_);
218     runner.useTopGroAndNdxFromDatabase("spc2");
219
220     runGrompp(&runner);
221
222     ASSERT_EQ(0, runner.callMdrun(*mdrunCaller_));
223 }
224
225 void MultiSimTest::runMaxhTest()
226 {
227     if (!mpiSetupValid())
228     {
229         // Can't test multi-sim without multiple simulations
230         return;
231     }
232
233     SimulationRunner runner(&fileManager_);
234     runner.useTopGroAndNdxFromDatabase("spc2");
235
236     TerminationHelper helper(&fileManager_, mdrunCaller_.get(), &runner);
237     // Make sure -maxh has a chance to propagate
238     int numSteps = 100;
239     runGrompp(&runner, numSteps);
240
241     helper.runFirstMdrun(runner.cptFileName_);
242     helper.runSecondMdrun();
243 }
244
245 } // namespace test
246 } // namespace gmx