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