2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019,2020, 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.
37 * \brief Implements the multi-simulation support routines.
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_mdrunutility
48 #include "gromacs/mdtypes/commrec.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/futil.h"
51 #include "gromacs/utility/gmxassert.h"
52 #include "gromacs/utility/mpiinplacebuffers.h"
53 #include "gromacs/utility/smalloc.h"
55 gmx_multisim_t::gmx_multisim_t() = default;
57 gmx_multisim_t::gmx_multisim_t(MPI_Comm comm, gmx::ArrayRef<const std::string> multidirs)
59 if (multidirs.empty())
64 if (!GMX_LIB_MPI && !multidirs.empty())
67 "mdrun -multidir is only supported when GROMACS has been "
68 "configured with a proper external MPI library.");
71 if (multidirs.size() == 1)
73 /* NOTE: It would be nice if this special case worked, but this requires checks/tests. */
75 "To run mdrun in multiple simulation mode, more then one "
76 "actual simulation is required. The single simulation case is not supported.");
81 MPI_Comm_size(comm, &numRanks);
82 if (numRanks % multidirs.size() != 0)
85 "The number of ranks (%d) is not a multiple of the number of simulations (%td)",
86 numRanks, multidirs.ssize());
89 int numRanksPerSim = numRanks / multidirs.size();
91 MPI_Comm_rank(comm, &rankWithinComm);
95 fprintf(debug, "We have %td simulations, %d ranks per simulation, local simulation is %d\n",
96 multidirs.ssize(), numRanksPerSim, rankWithinComm / numRanksPerSim);
99 nsim = multidirs.size();
100 sim = rankWithinComm / numRanksPerSim;
101 /* Create a communicator for the master nodes */
102 std::vector<int> rank(nsim);
103 for (int i = 0; i < nsim; i++)
105 rank[i] = i * numRanksPerSim;
107 MPI_Group mpi_group_world;
108 MPI_Comm_group(comm, &mpi_group_world);
109 MPI_Group_incl(mpi_group_world, nsim, rank.data(), &mpi_group_masters);
110 MPI_Comm_create(comm, mpi_group_masters, &mpi_comm_masters);
112 # if !MPI_IN_PLACE_EXISTS
113 /* initialize the MPI_IN_PLACE replacement buffers */
116 mpb->libuf = nullptr;
120 mpb->libuf_alloc = 0;
125 // TODO This should throw upon error
126 gmx_chdir(multidirs[sim].c_str());
128 GMX_UNUSED_VALUE(comm);
132 gmx_multisim_t::~gmx_multisim_t()
134 done_mpi_in_place_buf(mpb);
137 // TODO This would work better if the result of MPI_Comm_split was
138 // put into an RAII-style guard, such as gmx::unique_cptr.
139 if (mpi_comm_masters != MPI_COMM_NULL && mpi_comm_masters != MPI_COMM_WORLD)
141 MPI_Comm_free(&mpi_comm_masters);
143 if (mpi_group_masters != MPI_GROUP_NULL)
145 MPI_Group_free(&mpi_group_masters);
151 static void gmx_sumd_comm(int nr, double r[], MPI_Comm mpi_comm)
153 # if MPI_IN_PLACE_EXISTS
154 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
156 /* this function is only used in code that is not performance critical,
157 (during setup, when comm_rec is not the appropriate communication
158 structure), so this isn't as bad as it looks. */
163 MPI_Allreduce(r, buf, nr, MPI_DOUBLE, MPI_SUM, mpi_comm);
164 for (i = 0; i < nr; i++)
174 static void gmx_sumf_comm(int nr, float r[], MPI_Comm mpi_comm)
176 # if MPI_IN_PLACE_EXISTS
177 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
179 /* this function is only used in code that is not performance critical,
180 (during setup, when comm_rec is not the appropriate communication
181 structure), so this isn't as bad as it looks. */
186 MPI_Allreduce(r, buf, nr, MPI_FLOAT, MPI_SUM, mpi_comm);
187 for (i = 0; i < nr; i++)
196 void gmx_sumd_sim(int gmx_unused nr, double gmx_unused r[], const gmx_multisim_t gmx_unused* ms)
199 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumd_sim");
201 gmx_sumd_comm(nr, r, ms->mpi_comm_masters);
205 void gmx_sumf_sim(int gmx_unused nr, float gmx_unused r[], const gmx_multisim_t gmx_unused* ms)
208 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumf_sim");
210 gmx_sumf_comm(nr, r, ms->mpi_comm_masters);
214 void gmx_sumi_sim(int gmx_unused nr, int gmx_unused r[], const gmx_multisim_t gmx_unused* ms)
217 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumi_sim");
219 # if MPI_IN_PLACE_EXISTS
220 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
222 /* this is thread-unsafe, but it will do for now: */
225 if (nr > ms->mpb->ibuf_alloc)
227 ms->mpb->ibuf_alloc = nr;
228 srenew(ms->mpb->ibuf, ms->mpb->ibuf_alloc);
230 MPI_Allreduce(r, ms->mpb->ibuf, nr, MPI_INT, MPI_SUM, ms->mpi_comm_masters);
231 for (i = 0; i < nr; i++)
233 r[i] = ms->mpb->ibuf[i];
239 void gmx_sumli_sim(int gmx_unused nr, int64_t gmx_unused r[], const gmx_multisim_t gmx_unused* ms)
242 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumli_sim");
244 # if MPI_IN_PLACE_EXISTS
245 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, ms->mpi_comm_masters);
247 /* this is thread-unsafe, but it will do for now: */
250 if (nr > ms->mpb->libuf_alloc)
252 ms->mpb->libuf_alloc = nr;
253 srenew(ms->mpb->libuf, ms->mpb->libuf_alloc);
255 MPI_Allreduce(r, ms->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, ms->mpi_comm_masters);
256 for (i = 0; i < nr; i++)
258 r[i] = ms->mpb->libuf[i];
264 std::vector<int> gatherIntFromMultiSimulation(const gmx_multisim_t* ms, const int localValue)
266 std::vector<int> valuesFromAllRanks;
270 valuesFromAllRanks.resize(ms->nsim);
271 valuesFromAllRanks[ms->sim] = localValue;
272 gmx_sumi_sim(ms->nsim, valuesFromAllRanks.data(), ms);
276 valuesFromAllRanks.emplace_back(localValue);
279 GMX_UNUSED_VALUE(ms);
280 valuesFromAllRanks.emplace_back(localValue);
282 return valuesFromAllRanks;
285 void check_multi_int(FILE* log, const gmx_multisim_t* ms, int val, const char* name, gmx_bool bQuiet)
288 gmx_bool bCompatible;
290 if (nullptr != log && !bQuiet)
292 fprintf(log, "Multi-checking %s ... ", name);
297 gmx_fatal(FARGS, "check_multi_int called with a NULL communication pointer");
300 snew(ibuf, ms->nsim);
302 gmx_sumi_sim(ms->nsim, ibuf, ms);
305 for (p = 1; p < ms->nsim; p++)
307 bCompatible = bCompatible && (ibuf[p - 1] == ibuf[p]);
312 if (nullptr != log && !bQuiet)
314 fprintf(log, "OK\n");
321 fprintf(log, "\n%s is not equal for all subsystems\n", name);
322 for (p = 0; p < ms->nsim; p++)
324 fprintf(log, " subsystem %d: %d\n", p, ibuf[p]);
327 gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
333 void check_multi_int64(FILE* log, const gmx_multisim_t* ms, int64_t val, const char* name, gmx_bool bQuiet)
337 gmx_bool bCompatible;
339 if (nullptr != log && !bQuiet)
341 fprintf(log, "Multi-checking %s ... ", name);
346 gmx_fatal(FARGS, "check_multi_int called with a NULL communication pointer");
349 snew(ibuf, ms->nsim);
351 gmx_sumli_sim(ms->nsim, ibuf, ms);
354 for (p = 1; p < ms->nsim; p++)
356 bCompatible = bCompatible && (ibuf[p - 1] == ibuf[p]);
361 if (nullptr != log && !bQuiet)
363 fprintf(log, "OK\n");
368 // TODO Part of this error message would also be good to go to
369 // stderr (from one rank of one sim only)
372 fprintf(log, "\n%s is not equal for all subsystems\n", name);
373 for (p = 0; p < ms->nsim; p++)
376 /* first make the format string */
377 snprintf(strbuf, 255, " subsystem %%d: %s\n", "%" PRId64);
378 fprintf(log, strbuf, p, ibuf[p]);
381 gmx_fatal(FARGS, "The %d subsystems are not compatible\n", ms->nsim);
387 bool findIsSimulationMasterRank(const gmx_multisim_t* ms, MPI_Comm communicator)
391 // Ranks of multi-simulations know whether they are a master
392 // rank. Ranks of non-multi simulation do not know until a
393 // t_commrec is available.
394 if ((ms != nullptr) && (ms->nsim > 1))
396 return ms->mpi_comm_masters != MPI_COMM_NULL;
402 MPI_Comm_rank(communicator, &rank);
407 else if (GMX_THREAD_MPI)
409 GMX_RELEASE_ASSERT(communicator == MPI_COMM_NULL || communicator == MPI_COMM_WORLD,
410 "Invalid communicator");
411 // Spawned threads have MPI_COMM_WORLD upon creation, so if
412 // the communicator is MPI_COMM_NULL this is not a spawned thread,
413 // ie is the master thread
414 return (communicator == MPI_COMM_NULL);
418 // No MPI means it must be the master (and only) rank.
423 bool isMasterSim(const gmx_multisim_t* ms)
425 return !isMultiSim(ms) || ms->sim == 0;
428 bool isMasterSimMasterRank(const gmx_multisim_t* ms, const bool isMaster)
430 return (isMaster && isMasterSim(ms));