Split lines with many copyright years
[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(const t_commrec gmx_unused* cr)
243 {
244 #if !GMX_MPI
245     gmx_call("gmx_barrier");
246 #else
247     MPI_Barrier(cr->mpi_comm_mygroup);
248 #endif
249 }
250
251 void gmx_bcast(int gmx_unused nbytes, void gmx_unused* b, const t_commrec gmx_unused* cr)
252 {
253 #if !GMX_MPI
254     gmx_call("gmx_bast");
255 #else
256     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mygroup);
257 #endif
258 }
259
260 void gmx_bcast_sim(int gmx_unused nbytes, void gmx_unused* b, const t_commrec gmx_unused* cr)
261 {
262 #if !GMX_MPI
263     gmx_call("gmx_bast");
264 #else
265     MPI_Bcast(b, nbytes, MPI_BYTE, MASTERRANK(cr), cr->mpi_comm_mysim);
266 #endif
267 }
268
269 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused* cr)
270 {
271 #if !GMX_MPI
272     gmx_call("gmx_sumd");
273 #else
274 #    if MPI_IN_PLACE_EXISTS
275     if (cr->nc.bUse)
276     {
277         if (cr->nc.rank_intra == 0)
278         {
279             /* Use two step summing. */
280             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
281             /* Sum the roots of the internal (intra) buffers. */
282             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
283         }
284         else
285         {
286             /* This is here because of the silly MPI specification
287                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
288             MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
289         }
290         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
291     }
292     else
293     {
294         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
295     }
296 #    else
297     int i;
298
299     if (nr > cr->mpb->dbuf_alloc)
300     {
301         cr->mpb->dbuf_alloc = nr;
302         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
303     }
304     if (cr->nc.bUse)
305     {
306         /* Use two step summing */
307         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
308         if (cr->nc.rank_intra == 0)
309         {
310             /* Sum with the buffers reversed */
311             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
312         }
313         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
314     }
315     else
316     {
317         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
318         for (i = 0; i < nr; i++)
319         {
320             r[i] = cr->mpb->dbuf[i];
321         }
322     }
323 #    endif
324 #endif
325 }
326
327 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused* cr)
328 {
329 #if !GMX_MPI
330     gmx_call("gmx_sumf");
331 #else
332 #    if MPI_IN_PLACE_EXISTS
333     if (cr->nc.bUse)
334     {
335         /* Use two step summing.  */
336         if (cr->nc.rank_intra == 0)
337         {
338             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
339             /* Sum the roots of the internal (intra) buffers */
340             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
341         }
342         else
343         {
344             /* This is here because of the silly MPI specification
345                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
346             MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
347         }
348         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
349     }
350     else
351     {
352         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
353     }
354 #    else
355     int i;
356
357     if (nr > cr->mpb->fbuf_alloc)
358     {
359         cr->mpb->fbuf_alloc = nr;
360         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
361     }
362     if (cr->nc.bUse)
363     {
364         /* Use two step summing */
365         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
366         if (cr->nc.rank_intra == 0)
367         {
368             /* Sum with the buffers reversed */
369             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
370         }
371         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
372     }
373     else
374     {
375         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
376         for (i = 0; i < nr; i++)
377         {
378             r[i] = cr->mpb->fbuf[i];
379         }
380     }
381 #    endif
382 #endif
383 }
384
385 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused* cr)
386 {
387 #if !GMX_MPI
388     gmx_call("gmx_sumi");
389 #else
390 #    if MPI_IN_PLACE_EXISTS
391     if (cr->nc.bUse)
392     {
393         /* Use two step summing */
394         if (cr->nc.rank_intra == 0)
395         {
396             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
397             /* Sum with the buffers reversed */
398             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
399         }
400         else
401         {
402             /* This is here because of the silly MPI specification
403                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
404             MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
405         }
406         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
407     }
408     else
409     {
410         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
411     }
412 #    else
413     int i;
414
415     if (nr > cr->mpb->ibuf_alloc)
416     {
417         cr->mpb->ibuf_alloc = nr;
418         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
419     }
420     if (cr->nc.bUse)
421     {
422         /* Use two step summing */
423         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
424         if (cr->nc.rank_intra == 0)
425         {
426             /* Sum with the buffers reversed */
427             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
428         }
429         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
430     }
431     else
432     {
433         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
434         for (i = 0; i < nr; i++)
435         {
436             r[i] = cr->mpb->ibuf[i];
437         }
438     }
439 #    endif
440 #endif
441 }
442
443 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused* cr)
444 {
445 #if !GMX_MPI
446     gmx_call("gmx_sumli");
447 #else
448 #    if MPI_IN_PLACE_EXISTS
449     if (cr->nc.bUse)
450     {
451         /* Use two step summing */
452         if (cr->nc.rank_intra == 0)
453         {
454             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
455             /* Sum with the buffers reversed */
456             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
457         }
458         else
459         {
460             /* This is here because of the silly MPI specification
461                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
462             MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
463         }
464         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
465     }
466     else
467     {
468         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
469     }
470 #    else
471     int i;
472
473     if (nr > cr->mpb->libuf_alloc)
474     {
475         cr->mpb->libuf_alloc = nr;
476         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
477     }
478     if (cr->nc.bUse)
479     {
480         /* Use two step summing */
481         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_intra);
482         if (cr->nc.rank_intra == 0)
483         {
484             /* Sum with the buffers reversed */
485             MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
486         }
487         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
488     }
489     else
490     {
491         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
492         for (i = 0; i < nr; i++)
493         {
494             r[i] = cr->mpb->libuf[i];
495         }
496     }
497 #    endif
498 #endif
499 }
500
501 const char* opt2fn_master(const char* opt, int nfile, const t_filenm fnm[], t_commrec* cr)
502 {
503     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
504 }
505
506 void gmx_fatal_collective(int                    f_errno,
507                           const char*            file,
508                           int                    line,
509                           MPI_Comm               comm,
510                           gmx_bool               bMaster,
511                           gmx_fmtstr const char* fmt,
512                           ...)
513 {
514     va_list  ap;
515     gmx_bool bFinalize;
516 #if GMX_MPI
517     int result;
518     /* Check if we are calling on all processes in MPI_COMM_WORLD */
519     MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
520     /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
521     bFinalize = (result != MPI_UNEQUAL);
522 #else
523     GMX_UNUSED_VALUE(comm);
524     bFinalize = TRUE;
525 #endif
526
527     va_start(ap, fmt);
528     gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
529     va_end(ap);
530 }
531
532 void simulationBarrier(const t_commrec* cr)
533 {
534     if (PAR(cr))
535     {
536 #if GMX_MPI
537         MPI_Barrier(cr->mpi_comm_mysim);
538 #endif
539     }
540 }