Require explicit MPI_COMM for gmx_bcast and gmx_barrier
[alexxy/gromacs.git] / src / gromacs / gmxlib / network.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #include "gmxpre.h"
39
40 #include "network.h"
41
42 #include "config.h"
43
44 #include <cctype>
45 #include <cstdarg>
46 #include <cstdlib>
47 #include <cstring>
48
49 #include "gromacs/commandline/filenm.h"
50 #include "gromacs/mdrunutility/multisim.h"
51 #include "gromacs/mdtypes/commrec.h"
52 #include "gromacs/utility/basenetwork.h"
53 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/gmxmpi.h"
57 #include "gromacs/utility/mpiinplacebuffers.h"
58 #include "gromacs/utility/real.h"
59 #include "gromacs/utility/smalloc.h"
60
61 /* The source code in this file should be thread-safe.
62       Please keep it that way. */
63
64 CommrecHandle init_commrec(MPI_Comm communicator, const gmx_multisim_t* ms)
65 {
66     CommrecHandle handle;
67     t_commrec*    cr;
68
69     snew(cr, 1);
70     handle.reset(cr);
71
72     int rankInCommunicator, sizeOfCommunicator;
73 #if GMX_MPI
74 #    if GMX_LIB_MPI
75     GMX_RELEASE_ASSERT(gmx_mpi_initialized(), "Must have initialized MPI before building commrec");
76 #    endif
77     MPI_Comm_rank(communicator, &rankInCommunicator);
78     MPI_Comm_size(communicator, &sizeOfCommunicator);
79 #else
80     GMX_UNUSED_VALUE(communicator);
81     rankInCommunicator = 0;
82     sizeOfCommunicator = 1;
83 #endif
84
85     if (ms != nullptr)
86     {
87 #if GMX_MPI
88         cr->nnodes = sizeOfCommunicator / ms->nsim;
89         MPI_Comm_split(communicator, ms->sim, rankInCommunicator, &cr->mpi_comm_mysim);
90         cr->mpi_comm_mygroup = cr->mpi_comm_mysim;
91         MPI_Comm_rank(cr->mpi_comm_mysim, &cr->sim_nodeid);
92         MPI_Comm_rank(cr->mpi_comm_mygroup, &cr->nodeid);
93 #endif
94     }
95     else
96     {
97         cr->nnodes           = sizeOfCommunicator;
98         cr->nodeid           = rankInCommunicator;
99         cr->sim_nodeid       = cr->nodeid;
100         cr->mpi_comm_mysim   = communicator;
101         cr->mpi_comm_mygroup = communicator;
102     }
103
104     // TODO cr->duty should not be initialized here
105     cr->duty = (DUTY_PP | DUTY_PME);
106
107 #if GMX_MPI && !MPI_IN_PLACE_EXISTS
108     /* initialize the MPI_IN_PLACE replacement buffers */
109     snew(cr->mpb, 1);
110     cr->mpb->ibuf        = nullptr;
111     cr->mpb->libuf       = nullptr;
112     cr->mpb->fbuf        = nullptr;
113     cr->mpb->dbuf        = nullptr;
114     cr->mpb->ibuf_alloc  = 0;
115     cr->mpb->libuf_alloc = 0;
116     cr->mpb->fbuf_alloc  = 0;
117     cr->mpb->dbuf_alloc  = 0;
118 #endif
119
120     return handle;
121 }
122
123 void done_commrec(t_commrec* cr)
124 {
125     if (MASTER(cr))
126     {
127         if (nullptr != cr->dd)
128         {
129             // TODO: implement
130             // done_domdec(cr->dd);
131         }
132         done_mpi_in_place_buf(cr->mpb);
133     }
134 #if GMX_MPI
135     // TODO We need to be able to free communicators, but the
136     // structure of the commrec and domdec initialization code makes
137     // it hard to avoid both leaks and double frees.
138     bool mySimIsMyGroup = (cr->mpi_comm_mysim == cr->mpi_comm_mygroup);
139     if (cr->mpi_comm_mysim != MPI_COMM_NULL && cr->mpi_comm_mysim != MPI_COMM_WORLD)
140     {
141         // TODO see above
142         // MPI_Comm_free(&cr->mpi_comm_mysim);
143     }
144     if (!mySimIsMyGroup && cr->mpi_comm_mygroup != MPI_COMM_NULL && cr->mpi_comm_mygroup != MPI_COMM_WORLD)
145     {
146         // TODO see above
147         // MPI_Comm_free(&cr->mpi_comm_mygroup);
148     }
149 #endif
150     sfree(cr);
151 }
152
153 void gmx_setup_nodecomm(FILE gmx_unused* fplog, t_commrec* cr)
154 {
155     gmx_nodecomm_t* nc;
156
157     /* Many MPI implementations do not optimize MPI_Allreduce
158      * (and probably also other global communication calls)
159      * for multi-core nodes connected by a network.
160      * We can optimize such communication by using one MPI call
161      * within each node and one between the nodes.
162      * For MVAPICH2 and Intel MPI this reduces the time for
163      * the global_stat communication by 25%
164      * for 2x2-core 3 GHz Woodcrest connected by mixed DDR/SDR Infiniband.
165      * B. Hess, November 2007
166      */
167
168     nc = &cr->nc;
169
170     nc->bUse = FALSE;
171 #if !GMX_THREAD_MPI
172 #    if GMX_MPI
173     int n, rank;
174
175     // TODO PhysicalNodeCommunicator could be extended/used to handle
176     // the need for per-node per-group communicators.
177     MPI_Comm_size(cr->mpi_comm_mygroup, &n);
178     MPI_Comm_rank(cr->mpi_comm_mygroup, &rank);
179
180     int nodehash = gmx_physicalnode_id_hash();
181
182     if (debug)
183     {
184         fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
185     }
186
187
188     /* The intra-node communicator, split on node number */
189     MPI_Comm_split(cr->mpi_comm_mygroup, nodehash, rank, &nc->comm_intra);
190     MPI_Comm_rank(nc->comm_intra, &nc->rank_intra);
191     if (debug)
192     {
193         fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n", rank, nc->rank_intra);
194     }
195     /* The inter-node communicator, split on rank_intra.
196      * We actually only need the one for rank=0,
197      * but it is easier to create them all.
198      */
199     MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
200     /* Check if this really created two step communication */
201     int ng, ni;
202
203     MPI_Comm_size(nc->comm_inter, &ng);
204     MPI_Comm_size(nc->comm_intra, &ni);
205     if (debug)
206     {
207         fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n", ng, ni);
208     }
209
210     if (getenv("GMX_NO_NODECOMM") == nullptr && ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
211     {
212         nc->bUse = TRUE;
213         if (fplog)
214         {
215             fprintf(fplog, "Using two step summing over %d groups of on average %.1f ranks\n\n", ng,
216                     (real)n / (real)ng);
217         }
218         if (nc->rank_intra > 0)
219         {
220             MPI_Comm_free(&nc->comm_inter);
221         }
222     }
223     else
224     {
225         /* One group or all processes in a separate group, use normal summing */
226         MPI_Comm_free(&nc->comm_inter);
227         MPI_Comm_free(&nc->comm_intra);
228         if (debug)
229         {
230             fprintf(debug,
231                     "In gmx_setup_nodecomm: not unsing separate inter- and intra-node "
232                     "communicators.\n");
233         }
234     }
235 #    endif
236 #else
237     /* tMPI runs only on a single node so just use the nodeid */
238     nc->rank_intra = cr->nodeid;
239 #endif
240 }
241
242 void gmx_barrier(MPI_Comm gmx_unused communicator)
243 {
244 #if !GMX_MPI
245     gmx_call("gmx_barrier");
246 #else
247     MPI_Barrier(communicator);
248 #endif
249 }
250
251 void gmx_bcast(int gmx_unused nbytes, void gmx_unused* b, MPI_Comm gmx_unused communicator)
252 {
253 #if !GMX_MPI
254     gmx_call("gmx_bcast");
255 #else
256     MPI_Bcast(b, nbytes, MPI_BYTE, 0, communicator);
257 #endif
258 }
259
260 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused* cr)
261 {
262 #if !GMX_MPI
263     gmx_call("gmx_sumd");
264 #else
265 #    if MPI_IN_PLACE_EXISTS
266     if (cr->nc.bUse)
267     {
268         if (cr->nc.rank_intra == 0)
269         {
270             /* Use two step summing. */
271             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
272             /* Sum the roots of the internal (intra) buffers. */
273             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
274         }
275         else
276         {
277             /* This is here because of the silly MPI specification
278                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
279             MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
280         }
281         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
282     }
283     else
284     {
285         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
286     }
287 #    else
288     int i;
289
290     if (nr > cr->mpb->dbuf_alloc)
291     {
292         cr->mpb->dbuf_alloc = nr;
293         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
294     }
295     if (cr->nc.bUse)
296     {
297         /* Use two step summing */
298         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
299         if (cr->nc.rank_intra == 0)
300         {
301             /* Sum with the buffers reversed */
302             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
303         }
304         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
305     }
306     else
307     {
308         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
309         for (i = 0; i < nr; i++)
310         {
311             r[i] = cr->mpb->dbuf[i];
312         }
313     }
314 #    endif
315 #endif
316 }
317
318 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused* cr)
319 {
320 #if !GMX_MPI
321     gmx_call("gmx_sumf");
322 #else
323 #    if MPI_IN_PLACE_EXISTS
324     if (cr->nc.bUse)
325     {
326         /* Use two step summing.  */
327         if (cr->nc.rank_intra == 0)
328         {
329             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
330             /* Sum the roots of the internal (intra) buffers */
331             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
332         }
333         else
334         {
335             /* This is here because of the silly MPI specification
336                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
337             MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
338         }
339         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
340     }
341     else
342     {
343         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
344     }
345 #    else
346     int i;
347
348     if (nr > cr->mpb->fbuf_alloc)
349     {
350         cr->mpb->fbuf_alloc = nr;
351         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
352     }
353     if (cr->nc.bUse)
354     {
355         /* Use two step summing */
356         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
357         if (cr->nc.rank_intra == 0)
358         {
359             /* Sum with the buffers reversed */
360             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
361         }
362         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
363     }
364     else
365     {
366         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
367         for (i = 0; i < nr; i++)
368         {
369             r[i] = cr->mpb->fbuf[i];
370         }
371     }
372 #    endif
373 #endif
374 }
375
376 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused* cr)
377 {
378 #if !GMX_MPI
379     gmx_call("gmx_sumi");
380 #else
381 #    if MPI_IN_PLACE_EXISTS
382     if (cr->nc.bUse)
383     {
384         /* Use two step summing */
385         if (cr->nc.rank_intra == 0)
386         {
387             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
388             /* Sum with the buffers reversed */
389             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
390         }
391         else
392         {
393             /* This is here because of the silly MPI specification
394                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
395             MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
396         }
397         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
398     }
399     else
400     {
401         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
402     }
403 #    else
404     int i;
405
406     if (nr > cr->mpb->ibuf_alloc)
407     {
408         cr->mpb->ibuf_alloc = nr;
409         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
410     }
411     if (cr->nc.bUse)
412     {
413         /* Use two step summing */
414         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
415         if (cr->nc.rank_intra == 0)
416         {
417             /* Sum with the buffers reversed */
418             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
419         }
420         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
421     }
422     else
423     {
424         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
425         for (i = 0; i < nr; i++)
426         {
427             r[i] = cr->mpb->ibuf[i];
428         }
429     }
430 #    endif
431 #endif
432 }
433
434 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused* cr)
435 {
436 #if !GMX_MPI
437     gmx_call("gmx_sumli");
438 #else
439 #    if MPI_IN_PLACE_EXISTS
440     if (cr->nc.bUse)
441     {
442         /* Use two step summing */
443         if (cr->nc.rank_intra == 0)
444         {
445             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
446             /* Sum with the buffers reversed */
447             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
448         }
449         else
450         {
451             /* This is here because of the silly MPI specification
452                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
453             MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
454         }
455         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
456     }
457     else
458     {
459         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
460     }
461 #    else
462     int i;
463
464     if (nr > cr->mpb->libuf_alloc)
465     {
466         cr->mpb->libuf_alloc = nr;
467         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
468     }
469     if (cr->nc.bUse)
470     {
471         /* Use two step summing */
472         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_intra);
473         if (cr->nc.rank_intra == 0)
474         {
475             /* Sum with the buffers reversed */
476             MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
477         }
478         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
479     }
480     else
481     {
482         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
483         for (i = 0; i < nr; i++)
484         {
485             r[i] = cr->mpb->libuf[i];
486         }
487     }
488 #    endif
489 #endif
490 }
491
492 const char* opt2fn_master(const char* opt, int nfile, const t_filenm fnm[], t_commrec* cr)
493 {
494     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
495 }
496
497 void gmx_fatal_collective(int                    f_errno,
498                           const char*            file,
499                           int                    line,
500                           MPI_Comm               comm,
501                           gmx_bool               bMaster,
502                           gmx_fmtstr const char* fmt,
503                           ...)
504 {
505     va_list  ap;
506     gmx_bool bFinalize;
507 #if GMX_MPI
508     int result;
509     /* Check if we are calling on all processes in MPI_COMM_WORLD */
510     MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
511     /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
512     bFinalize = (result != MPI_UNEQUAL);
513 #else
514     GMX_UNUSED_VALUE(comm);
515     bFinalize = TRUE;
516 #endif
517
518     va_start(ap, fmt);
519     gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
520     va_end(ap);
521 }