2 * This file is part of the GROMACS molecular simulation package.
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.
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/gmxlib/network.h"
49 #include "gromacs/mdtypes/commrec.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/utility/smalloc.h"
54 gmx_multisim_t *init_multisystem(MPI_Comm comm,
55 gmx::ArrayRef<const std::string> multidirs)
59 MPI_Group mpi_group_world;
63 if (multidirs.empty())
68 if (!GMX_LIB_MPI && !multidirs.empty())
70 gmx_fatal(FARGS, "mdrun -multidir is only supported when GROMACS has been "
71 "configured with a proper external MPI library.");
74 if (multidirs.size() == 1)
76 /* NOTE: It would be nice if this special case worked, but this requires checks/tests. */
77 gmx_fatal(FARGS, "To run mdrun in multiple simulation mode, more then one "
78 "actual simulation is required. The single simulation case is not supported.");
83 MPI_Comm_size(comm, &numRanks);
84 if (numRanks % multidirs.size() != 0)
86 gmx_fatal(FARGS, "The number of ranks (%d) is not a multiple of the number of simulations (%td)", 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", multidirs.ssize(), numRanksPerSim, rankWithinComm/numRanksPerSim);
98 ms = new gmx_multisim_t;
99 ms->nsim = multidirs.size();
100 ms->sim = rankWithinComm/numRanksPerSim;
101 /* Create a communicator for the master nodes */
102 snew(rank, ms->nsim);
103 for (int i = 0; i < ms->nsim; i++)
105 rank[i] = i*numRanksPerSim;
107 MPI_Comm_group(comm, &mpi_group_world);
108 MPI_Group_incl(mpi_group_world, ms->nsim, rank, &ms->mpi_group_masters);
110 MPI_Comm_create(MPI_COMM_WORLD, ms->mpi_group_masters,
111 &ms->mpi_comm_masters);
113 #if !MPI_IN_PLACE_EXISTS
114 /* initialize the MPI_IN_PLACE replacement buffers */
116 ms->mpb->ibuf = NULL;
117 ms->mpb->libuf = NULL;
118 ms->mpb->fbuf = NULL;
119 ms->mpb->dbuf = NULL;
120 ms->mpb->ibuf_alloc = 0;
121 ms->mpb->libuf_alloc = 0;
122 ms->mpb->fbuf_alloc = 0;
123 ms->mpb->dbuf_alloc = 0;
126 // TODO This should throw upon error
127 gmx_chdir(multidirs[ms->sim].c_str());
129 GMX_UNUSED_VALUE(comm);
136 void done_multisim(gmx_multisim_t *ms)
140 done_mpi_in_place_buf(ms->mpb);
145 bool isMasterSim(const gmx_multisim_t *ms)
147 return !isMultiSim(ms) || ms->sim == 0;
150 bool isMasterSimMasterRank(const gmx_multisim_t *ms,
153 return (isMaster && isMasterSim(ms));