2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
49 #include "gromacs/commandline/filenm.h"
50 #include "gromacs/mdtypes/commrec.h"
51 #include "gromacs/utility/basenetwork.h"
52 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/gmxmpi.h"
57 #include "gromacs/utility/real.h"
58 #include "gromacs/utility/smalloc.h"
60 /* The source code in this file should be thread-safe.
61 Please keep it that way. */
63 CommrecHandle init_commrec(MPI_Comm communicator)
71 int rankInCommunicator, sizeOfCommunicator;
74 GMX_RELEASE_ASSERT(gmx_mpi_initialized(), "Must have initialized MPI before building commrec");
76 MPI_Comm_rank(communicator, &rankInCommunicator);
77 MPI_Comm_size(communicator, &sizeOfCommunicator);
79 GMX_UNUSED_VALUE(communicator);
80 rankInCommunicator = 0;
81 sizeOfCommunicator = 1;
84 cr->mpiDefaultCommunicator = communicator;
85 cr->sizeOfDefaultCommunicator = sizeOfCommunicator;
86 cr->rankInDefaultCommunicator = rankInCommunicator;
88 // For now, we want things to go horribly wrong if this is used too early...
89 // TODO: Remove when communicators are removed from commrec (#2395)
91 cr->sizeOfMyGroupCommunicator = -1;
94 cr->mpi_comm_mysim = MPI_COMM_NULL;
95 cr->mpi_comm_mygroup = MPI_COMM_NULL;
97 // TODO cr->duty should not be initialized here
98 cr->duty = (DUTY_PP | DUTY_PME);
103 void done_commrec(t_commrec* cr)
107 if (nullptr != cr->dd)
110 // done_domdec(cr->dd);
114 // TODO We need to be able to free communicators, but the
115 // structure of the commrec and domdec initialization code makes
116 // it hard to avoid both leaks and double frees.
117 bool mySimIsMyGroup = (cr->mpi_comm_mysim == cr->mpi_comm_mygroup);
118 if (cr->mpi_comm_mysim != MPI_COMM_NULL && cr->mpi_comm_mysim != MPI_COMM_WORLD)
121 // MPI_Comm_free(&cr->mpi_comm_mysim);
123 if (!mySimIsMyGroup && cr->mpi_comm_mygroup != MPI_COMM_NULL && cr->mpi_comm_mygroup != MPI_COMM_WORLD)
126 // MPI_Comm_free(&cr->mpi_comm_mygroup);
132 void gmx_setup_nodecomm(FILE gmx_unused* fplog, t_commrec* cr)
136 /* Many MPI implementations do not optimize MPI_Allreduce
137 * (and probably also other global communication calls)
138 * for multi-core nodes connected by a network.
139 * We can optimize such communication by using one MPI call
140 * within each node and one between the nodes.
141 * For MVAPICH2 and Intel MPI this reduces the time for
142 * the global_stat communication by 25%
143 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
144 * B. Hess, November 2007
154 // TODO PhysicalNodeCommunicator could be extended/used to handle
155 // the need for per-node per-group communicators.
156 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
157 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
159 int nodehash = gmx_physicalnode_id_hash();
163 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
167 /* The intra-node communicator, split on node number */
168 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
169 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
172 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n", rank, nc->rank_intra);
174 /* The inter-node communicator, split on rank_intra.
175 * We actually only need the one for rank=0,
176 * but it is easier to create them all.
178 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
179 /* Check if this really created two step communication */
182 MPI_Comm_size(nc->comm_inter, &ng);
183 MPI_Comm_size(nc->comm_intra, &ni);
186 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n", ng, ni);
189 if (getenv("GMX_NO_NODECOMM") == nullptr && ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
195 "Using two step summing over %d groups of on average %.1f ranks\n\n",
199 if (nc->rank_intra > 0)
201 MPI_Comm_free(&nc->comm_inter);
206 /* One group or all processes in a separate group, use normal summing */
207 MPI_Comm_free(&nc->comm_inter);
208 MPI_Comm_free(&nc->comm_intra);
212 "In gmx_setup_nodecomm: not unsing separate inter- and intra-node "
218 /* tMPI runs only on a single node so just use the nodeid */
219 nc->rank_intra = cr->nodeid;
223 void gmx_barrier(MPI_Comm gmx_unused communicator)
225 if (communicator == MPI_COMM_NULL)
230 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_barrier");
232 MPI_Barrier(communicator);
236 void gmx_bcast(int gmx_unused nbytes, void gmx_unused* b, MPI_Comm gmx_unused communicator)
238 // Without MPI we have a single rank, so bcast is a no-op
240 MPI_Bcast(b, nbytes, MPI_BYTE, 0, communicator);
244 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused* cr)
246 // Without MPI we have a single rank, so sum is a no-op
248 if (cr->sizeOfMyGroupCommunicator == 1)
255 if (cr->nc.rank_intra == 0)
257 /* Use two step summing. */
258 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
259 /* Sum the roots of the internal (intra) buffers. */
260 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
264 /* This is here because of the silly MPI specification
265 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
266 MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
268 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
272 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
277 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused* cr)
279 // Without MPI we have a single rank, so sum is a no-op
281 if (cr->sizeOfMyGroupCommunicator == 1)
288 /* Use two step summing. */
289 if (cr->nc.rank_intra == 0)
291 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
292 /* Sum the roots of the internal (intra) buffers */
293 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
297 /* This is here because of the silly MPI specification
298 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
299 MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
301 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
305 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
310 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused* cr)
312 // Without MPI we have a single rank, so sum is a no-op
314 if (cr->sizeOfMyGroupCommunicator == 1)
321 /* Use two step summing */
322 if (cr->nc.rank_intra == 0)
324 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
325 /* Sum with the buffers reversed */
326 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
330 /* This is here because of the silly MPI specification
331 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
332 MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
334 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
338 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
343 const char* opt2fn_master(const char* opt, int nfile, const t_filenm fnm[], t_commrec* cr)
345 return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
348 void gmx_fatal_collective(int f_errno,
353 gmx_fmtstr const char* fmt,
360 /* Check if we are calling on all processes in MPI_COMM_WORLD */
361 MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
362 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
363 bFinalize = (result != MPI_UNEQUAL);
365 GMX_UNUSED_VALUE(comm);
370 gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);