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, 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)
93 cr->mpi_comm_mysim = MPI_COMM_NULL;
94 cr->mpi_comm_mygroup = MPI_COMM_NULL;
96 // TODO cr->duty should not be initialized here
97 cr->duty = (DUTY_PP | DUTY_PME);
102 void done_commrec(t_commrec* cr)
106 if (nullptr != cr->dd)
109 // done_domdec(cr->dd);
113 // TODO We need to be able to free communicators, but the
114 // structure of the commrec and domdec initialization code makes
115 // it hard to avoid both leaks and double frees.
116 bool mySimIsMyGroup = (cr->mpi_comm_mysim == cr->mpi_comm_mygroup);
117 if (cr->mpi_comm_mysim != MPI_COMM_NULL && cr->mpi_comm_mysim != MPI_COMM_WORLD)
120 // MPI_Comm_free(&cr->mpi_comm_mysim);
122 if (!mySimIsMyGroup && cr->mpi_comm_mygroup != MPI_COMM_NULL && cr->mpi_comm_mygroup != MPI_COMM_WORLD)
125 // MPI_Comm_free(&cr->mpi_comm_mygroup);
131 void gmx_setup_nodecomm(FILE gmx_unused* fplog, t_commrec* cr)
135 /* Many MPI implementations do not optimize MPI_Allreduce
136 * (and probably also other global communication calls)
137 * for multi-core nodes connected by a network.
138 * We can optimize such communication by using one MPI call
139 * within each node and one between the nodes.
140 * For MVAPICH2 and Intel MPI this reduces the time for
141 * the global_stat communication by 25%
142 * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
143 * B. Hess, November 2007
153 // TODO PhysicalNodeCommunicator could be extended/used to handle
154 // the need for per-node per-group communicators.
155 MPI_Comm_size(cr->mpi_comm_mygroup, &n);
156 MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
158 int nodehash = gmx_physicalnode_id_hash();
162 fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
166 /* The intra-node communicator, split on node number */
167 MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
168 MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
171 fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n", rank, nc->rank_intra);
173 /* The inter-node communicator, split on rank_intra.
174 * We actually only need the one for rank=0,
175 * but it is easier to create them all.
177 MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
178 /* Check if this really created two step communication */
181 MPI_Comm_size(nc->comm_inter, &ng);
182 MPI_Comm_size(nc->comm_intra, &ni);
185 fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n", ng, ni);
188 if (getenv("GMX_NO_NODECOMM") == nullptr && ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
193 fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n", ng,
196 if (nc->rank_intra > 0)
198 MPI_Comm_free(&nc->comm_inter);
203 /* One group or all processes in a separate group, use normal summing */
204 MPI_Comm_free(&nc->comm_inter);
205 MPI_Comm_free(&nc->comm_intra);
209 "In gmx_setup_nodecomm: not unsing separate inter- and intra-node "
215 /* tMPI runs only on a single node so just use the nodeid */
216 nc->rank_intra = cr->nodeid;
220 void gmx_barrier(MPI_Comm gmx_unused communicator)
223 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_barrier");
225 MPI_Barrier(communicator);
229 void gmx_bcast(int gmx_unused nbytes, void gmx_unused* b, MPI_Comm gmx_unused communicator)
232 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_bcast");
234 MPI_Bcast(b, nbytes, MPI_BYTE, 0, communicator);
238 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused* cr)
241 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumd");
243 # if MPI_IN_PLACE_EXISTS
246 if (cr->nc.rank_intra == 0)
248 /* Use two step summing. */
249 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
250 /* Sum the roots of the internal (intra) buffers. */
251 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
255 /* This is here because of the silly MPI specification
256 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
257 MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
259 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
263 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
268 if (nr > cr->mpb->dbuf_alloc)
270 cr->mpb->dbuf_alloc = nr;
271 srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
275 /* Use two step summing */
276 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
277 if (cr->nc.rank_intra == 0)
279 /* Sum with the buffers reversed */
280 MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
282 MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
286 MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
287 for (i = 0; i < nr; i++)
289 r[i] = cr->mpb->dbuf[i];
296 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused* cr)
299 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumf");
301 # if MPI_IN_PLACE_EXISTS
304 /* Use two step summing. */
305 if (cr->nc.rank_intra == 0)
307 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
308 /* Sum the roots of the internal (intra) buffers */
309 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
313 /* This is here because of the silly MPI specification
314 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
315 MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
317 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
321 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
326 if (nr > cr->mpb->fbuf_alloc)
328 cr->mpb->fbuf_alloc = nr;
329 srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
333 /* Use two step summing */
334 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
335 if (cr->nc.rank_intra == 0)
337 /* Sum with the buffers reversed */
338 MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
340 MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
344 MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
345 for (i = 0; i < nr; i++)
347 r[i] = cr->mpb->fbuf[i];
354 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused* cr)
357 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumi");
359 # if MPI_IN_PLACE_EXISTS
362 /* Use two step summing */
363 if (cr->nc.rank_intra == 0)
365 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
366 /* Sum with the buffers reversed */
367 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
371 /* This is here because of the silly MPI specification
372 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
373 MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
375 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
379 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
384 if (nr > cr->mpb->ibuf_alloc)
386 cr->mpb->ibuf_alloc = nr;
387 srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
391 /* Use two step summing */
392 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
393 if (cr->nc.rank_intra == 0)
395 /* Sum with the buffers reversed */
396 MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
398 MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
402 MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
403 for (i = 0; i < nr; i++)
405 r[i] = cr->mpb->ibuf[i];
412 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused* cr)
415 GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumli");
417 # if MPI_IN_PLACE_EXISTS
420 /* Use two step summing */
421 if (cr->nc.rank_intra == 0)
423 MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
424 /* Sum with the buffers reversed */
425 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
429 /* This is here because of the silly MPI specification
430 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
431 MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
433 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
437 MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
442 if (nr > cr->mpb->libuf_alloc)
444 cr->mpb->libuf_alloc = nr;
445 srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
449 /* Use two step summing */
450 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_intra);
451 if (cr->nc.rank_intra == 0)
453 /* Sum with the buffers reversed */
454 MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
456 MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
460 MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
461 for (i = 0; i < nr; i++)
463 r[i] = cr->mpb->libuf[i];
470 const char* opt2fn_master(const char* opt, int nfile, const t_filenm fnm[], t_commrec* cr)
472 return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
475 void gmx_fatal_collective(int f_errno,
480 gmx_fmtstr const char* fmt,
487 /* Check if we are calling on all processes in MPI_COMM_WORLD */
488 MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
489 /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
490 bFinalize = (result != MPI_UNEQUAL);
492 GMX_UNUSED_VALUE(comm);
497 gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);