2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
39 * Tests for the mdrun multi-simulation functionality
41 * \author Mark Abraham <mark.j.abraham@gmail.com>
42 * \ingroup module_mdrun_integration_tests
46 #include "multisimtest.h"
53 #include <gtest/gtest.h>
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"
62 #include "testutils/cmdlinetest.h"
64 #include "moduletest.h"
65 #include "terminationhelper.h"
72 MultiSimTest::MultiSimTest() :
73 size_(gmx_node_num()),
74 rank_(gmx_node_rank()),
75 numRanksPerSimulation_(std::get<0>(GetParam())),
76 simulationNumber_(rank_ / numRanksPerSimulation_),
77 mdrunCaller_(new CommandLine)
80 // Zero or less ranks doesn't make sense
81 GMX_RELEASE_ASSERT(numRanksPerSimulation_ > 0, "Invalid number of ranks per simulation.");
83 const char* directoryNameFormat = "sim_%d";
85 // Modify the file manager to have a temporary directory unique to
86 // each simulation. No need to have a mutex on this, nobody else
87 // can access the fileManager_ yet because we only just
89 std::string originalTempDirectory = fileManager_.getOutputTempDirectory();
90 std::string newTempDirectory =
91 Path::join(originalTempDirectory, formatString(directoryNameFormat, simulationNumber_));
92 if (rank_ % numRanksPerSimulation_ == 0)
94 // Only one rank per simulation creates directory
95 Directory::create(newTempDirectory);
98 // Make sure directories got created.
99 MPI_Barrier(MdrunTestFixtureBase::communicator_);
101 fileManager_.setOutputTempDirectory(newTempDirectory);
103 mdrunCaller_->append("mdrun");
104 mdrunCaller_->addOption("-multidir");
105 for (int i = 0; i < size_ / numRanksPerSimulation_; ++i)
107 mdrunCaller_->append(Path::join(originalTempDirectory, formatString(directoryNameFormat, i)));
111 bool MultiSimTest::mpiSetupValid() const
113 // Single simulation case is not implemented in multi-sim
114 const bool haveAtLeastTwoSimulations = ((size_ / numRanksPerSimulation_) >= 2);
115 // Mdrun will throw error if simulations don't have identical number of ranks
116 const bool simulationsHaveIdenticalRankNumber = ((size_ % numRanksPerSimulation_) == 0);
118 return (haveAtLeastTwoSimulations && simulationsHaveIdenticalRankNumber);
121 void MultiSimTest::organizeMdpFile(SimulationRunner* runner,
122 IntegrationAlgorithm integrator,
123 TemperatureCoupling tcoupl,
124 PressureCoupling pcoupl,
126 bool doRegression) const
128 GMX_RELEASE_ASSERT(mpiSetupValid(), "Creating the mdp file without valid MPI setup is useless.");
129 const real baseTemperature = 298;
130 const real basePressure = 1;
131 std::string mdpFileContents = formatString(
137 "nstcalcenergy = 1\n"
141 // pressure coupling (if active)
144 "compressibility = 4.5e-5\n"
145 // velocity generation
149 // v-rescale and c-rescale use ld-seed
151 // Two systems are used: spc2 non-interacting also at cutoff 1nm
152 // tip3p5 has box length of 1.86
155 // Trajectory output if required
160 enumValueToString(integrator),
161 enumValueToString(tcoupl),
162 enumValueToString(pcoupl),
164 baseTemperature + 0.0001 * rank_,
165 basePressure * std::pow(1.01, rank_),
166 /* Set things up so that the initial KE decreases with
167 increasing replica number, so that the (identical)
168 starting PE decreases on the first step more for the
169 replicas with higher number, which will tend to force
170 replica exchange to occur. */
171 std::max(baseTemperature - 10 * rank_, real(0)),
172 // If we do regression, we need reproducible velocity
173 // generation, which can be different per simulation
174 (doRegression ? 671324 + simulationNumber_ : -1),
175 // If we do regression, we need reproducible temperature and
176 // pressure coupling, which can be different per simulation
177 (doRegression ? 51203 + simulationNumber_ : -1),
178 // If we do regression, write one intermediate point
179 (doRegression ? int(numSteps / 2) : 0),
180 (doRegression ? int(numSteps / 2) : 0),
181 (doRegression ? int(numSteps / 2) : 0),
182 // If we do regression, print energies every step so
183 // we're sure to catch the replica exchange steps
184 (doRegression ? 1 : 1000));
185 runner->useStringAsMdpFile(mdpFileContents);
188 void MultiSimTest::runGrompp(SimulationRunner* runner, int numSteps, bool doRegression, int maxWarnings) const
190 // Call grompp once per simulation
191 if (rank_ % numRanksPerSimulation_ == 0)
193 const auto& simulator = std::get<1>(GetParam());
194 const auto& tcoupl = std::get<2>(GetParam());
195 const auto& pcoupl = std::get<3>(GetParam());
196 organizeMdpFile(runner, simulator, tcoupl, pcoupl, numSteps, doRegression);
198 caller.addOption("-maxwarn", maxWarnings);
199 EXPECT_EQ(0, runner->callGromppOnThisRank(caller));
203 // Make sure simulation masters have written the .tpr file before other ranks try to read it.
204 MPI_Barrier(MdrunTestFixtureBase::communicator_);
208 void MultiSimTest::runExitsNormallyTest()
210 if (!mpiSetupValid())
212 // Can't test multi-sim without multiple simulations
216 SimulationRunner runner(&fileManager_);
217 runner.useTopGroAndNdxFromDatabase("spc2");
221 ASSERT_EQ(0, runner.callMdrun(*mdrunCaller_));
224 void MultiSimTest::runMaxhTest()
226 if (!mpiSetupValid())
228 // Can't test multi-sim without multiple simulations
232 SimulationRunner runner(&fileManager_);
233 runner.useTopGroAndNdxFromDatabase("spc2");
235 TerminationHelper helper(&fileManager_, mdrunCaller_.get(), &runner);
236 // Make sure -maxh has a chance to propagate
238 runGrompp(&runner, numSteps);
240 helper.runFirstMdrun(runner.cptFileName_);
241 helper.runSecondMdrun();