Apply re-formatting to C++ in src/ tree.
[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/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"
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)
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     cr->mpiDefaultCommunicator    = communicator;
85     cr->sizeOfDefaultCommunicator = sizeOfCommunicator;
86     cr->rankInDefaultCommunicator = rankInCommunicator;
87
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)
90     cr->nnodes           = -1;
91     cr->nodeid           = -1;
92     cr->sim_nodeid       = -1;
93     cr->mpi_comm_mysim   = MPI_COMM_NULL;
94     cr->mpi_comm_mygroup = MPI_COMM_NULL;
95
96     // TODO cr->duty should not be initialized here
97     cr->duty = (DUTY_PP | DUTY_PME);
98
99     return handle;
100 }
101
102 void done_commrec(t_commrec* cr)
103 {
104     if (MASTER(cr))
105     {
106         if (nullptr != cr->dd)
107         {
108             // TODO: implement
109             // done_domdec(cr->dd);
110         }
111     }
112 #if GMX_MPI
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)
118     {
119         // TODO see above
120         // MPI_Comm_free(&cr->mpi_comm_mysim);
121     }
122     if (!mySimIsMyGroup && cr->mpi_comm_mygroup != MPI_COMM_NULL && cr->mpi_comm_mygroup != MPI_COMM_WORLD)
123     {
124         // TODO see above
125         // MPI_Comm_free(&cr->mpi_comm_mygroup);
126     }
127 #endif
128     sfree(cr);
129 }
130
131 void gmx_setup_nodecomm(FILE gmx_unused* fplog, t_commrec* cr)
132 {
133     gmx_nodecomm_t* nc;
134
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
144      */
145
146     nc = &cr->nc;
147
148     nc->bUse = FALSE;
149 #if !GMX_THREAD_MPI
150 #    if GMX_MPI
151     int n, rank;
152
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);
157
158     int nodehash = gmx_physicalnode_id_hash();
159
160     if (debug)
161     {
162         fprintf(debug, "In gmx_setup_nodecomm: splitting communicator of size %d\n", n);
163     }
164
165
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);
169     if (debug)
170     {
171         fprintf(debug, "In gmx_setup_nodecomm: node ID %d rank within node %d\n", rank, nc->rank_intra);
172     }
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.
176      */
177     MPI_Comm_split(cr->mpi_comm_mygroup, nc->rank_intra, rank, &nc->comm_inter);
178     /* Check if this really created two step communication */
179     int ng, ni;
180
181     MPI_Comm_size(nc->comm_inter, &ng);
182     MPI_Comm_size(nc->comm_intra, &ni);
183     if (debug)
184     {
185         fprintf(debug, "In gmx_setup_nodecomm: groups %d, my group size %d\n", ng, ni);
186     }
187
188     if (getenv("GMX_NO_NODECOMM") == nullptr && ((ng > 1 && ng < n) || (ni > 1 && ni < n)))
189     {
190         nc->bUse = TRUE;
191         if (fplog)
192         {
193             fprintf(fplog,
194                     "Using two step summing over %d groups of on average %.1f ranks\n\n",
195                     ng,
196                     (real)n / (real)ng);
197         }
198         if (nc->rank_intra > 0)
199         {
200             MPI_Comm_free(&nc->comm_inter);
201         }
202     }
203     else
204     {
205         /* One group or all processes in a separate group, use normal summing */
206         MPI_Comm_free(&nc->comm_inter);
207         MPI_Comm_free(&nc->comm_intra);
208         if (debug)
209         {
210             fprintf(debug,
211                     "In gmx_setup_nodecomm: not unsing separate inter- and intra-node "
212                     "communicators.\n");
213         }
214     }
215 #    endif
216 #else
217     /* tMPI runs only on a single node so just use the nodeid */
218     nc->rank_intra = cr->nodeid;
219 #endif
220 }
221
222 void gmx_barrier(MPI_Comm gmx_unused communicator)
223 {
224 #if !GMX_MPI
225     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_barrier");
226 #else
227     MPI_Barrier(communicator);
228 #endif
229 }
230
231 void gmx_bcast(int gmx_unused nbytes, void gmx_unused* b, MPI_Comm gmx_unused communicator)
232 {
233 #if !GMX_MPI
234     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_bcast");
235 #else
236     MPI_Bcast(b, nbytes, MPI_BYTE, 0, communicator);
237 #endif
238 }
239
240 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused* cr)
241 {
242 #if !GMX_MPI
243     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumd");
244 #else
245 #    if MPI_IN_PLACE_EXISTS
246     if (cr->nc.bUse)
247     {
248         if (cr->nc.rank_intra == 0)
249         {
250             /* Use two step summing. */
251             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
252             /* Sum the roots of the internal (intra) buffers. */
253             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
254         }
255         else
256         {
257             /* This is here because of the silly MPI specification
258                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
259             MPI_Reduce(r, nullptr, nr, MPI_DOUBLE, MPI_SUM, 0, cr->nc.comm_intra);
260         }
261         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
262     }
263     else
264     {
265         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
266     }
267 #    else
268     int i;
269
270     if (nr > cr->mpb->dbuf_alloc)
271     {
272         cr->mpb->dbuf_alloc = nr;
273         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
274     }
275     if (cr->nc.bUse)
276     {
277         /* Use two step summing */
278         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_intra);
279         if (cr->nc.rank_intra == 0)
280         {
281             /* Sum with the buffers reversed */
282             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
283         }
284         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
285     }
286     else
287     {
288         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
289         for (i = 0; i < nr; i++)
290         {
291             r[i] = cr->mpb->dbuf[i];
292         }
293     }
294 #    endif
295 #endif
296 }
297
298 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused* cr)
299 {
300 #if !GMX_MPI
301     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumf");
302 #else
303 #    if MPI_IN_PLACE_EXISTS
304     if (cr->nc.bUse)
305     {
306         /* Use two step summing.  */
307         if (cr->nc.rank_intra == 0)
308         {
309             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
310             /* Sum the roots of the internal (intra) buffers */
311             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
312         }
313         else
314         {
315             /* This is here because of the silly MPI specification
316                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
317             MPI_Reduce(r, nullptr, nr, MPI_FLOAT, MPI_SUM, 0, cr->nc.comm_intra);
318         }
319         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
320     }
321     else
322     {
323         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
324     }
325 #    else
326     int i;
327
328     if (nr > cr->mpb->fbuf_alloc)
329     {
330         cr->mpb->fbuf_alloc = nr;
331         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
332     }
333     if (cr->nc.bUse)
334     {
335         /* Use two step summing */
336         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_intra);
337         if (cr->nc.rank_intra == 0)
338         {
339             /* Sum with the buffers reversed */
340             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
341         }
342         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
343     }
344     else
345     {
346         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
347         for (i = 0; i < nr; i++)
348         {
349             r[i] = cr->mpb->fbuf[i];
350         }
351     }
352 #    endif
353 #endif
354 }
355
356 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused* cr)
357 {
358 #if !GMX_MPI
359     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumi");
360 #else
361 #    if MPI_IN_PLACE_EXISTS
362     if (cr->nc.bUse)
363     {
364         /* Use two step summing */
365         if (cr->nc.rank_intra == 0)
366         {
367             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
368             /* Sum with the buffers reversed */
369             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
370         }
371         else
372         {
373             /* This is here because of the silly MPI specification
374                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
375             MPI_Reduce(r, nullptr, nr, MPI_INT, MPI_SUM, 0, cr->nc.comm_intra);
376         }
377         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
378     }
379     else
380     {
381         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
382     }
383 #    else
384     int i;
385
386     if (nr > cr->mpb->ibuf_alloc)
387     {
388         cr->mpb->ibuf_alloc = nr;
389         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
390     }
391     if (cr->nc.bUse)
392     {
393         /* Use two step summing */
394         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->nc.comm_intra);
395         if (cr->nc.rank_intra == 0)
396         {
397             /* Sum with the buffers reversed */
398             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
399         }
400         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
401     }
402     else
403     {
404         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
405         for (i = 0; i < nr; i++)
406         {
407             r[i] = cr->mpb->ibuf[i];
408         }
409     }
410 #    endif
411 #endif
412 }
413
414 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused* cr)
415 {
416 #if !GMX_MPI
417     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumli");
418 #else
419 #    if MPI_IN_PLACE_EXISTS
420     if (cr->nc.bUse)
421     {
422         /* Use two step summing */
423         if (cr->nc.rank_intra == 0)
424         {
425             MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
426             /* Sum with the buffers reversed */
427             MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
428         }
429         else
430         {
431             /* This is here because of the silly MPI specification
432                 that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
433             MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
434         }
435         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
436     }
437     else
438     {
439         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
440     }
441 #    else
442     int i;
443
444     if (nr > cr->mpb->libuf_alloc)
445     {
446         cr->mpb->libuf_alloc = nr;
447         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
448     }
449     if (cr->nc.bUse)
450     {
451         /* Use two step summing */
452         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_intra);
453         if (cr->nc.rank_intra == 0)
454         {
455             /* Sum with the buffers reversed */
456             MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
457         }
458         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
459     }
460     else
461     {
462         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
463         for (i = 0; i < nr; i++)
464         {
465             r[i] = cr->mpb->libuf[i];
466         }
467     }
468 #    endif
469 #endif
470 }
471
472 const char* opt2fn_master(const char* opt, int nfile, const t_filenm fnm[], t_commrec* cr)
473 {
474     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
475 }
476
477 void gmx_fatal_collective(int                    f_errno,
478                           const char*            file,
479                           int                    line,
480                           MPI_Comm               comm,
481                           gmx_bool               bMaster,
482                           gmx_fmtstr const char* fmt,
483                           ...)
484 {
485     va_list  ap;
486     gmx_bool bFinalize;
487 #if GMX_MPI
488     int result;
489     /* Check if we are calling on all processes in MPI_COMM_WORLD */
490     MPI_Comm_compare(comm, MPI_COMM_WORLD, &result);
491     /* Any result except MPI_UNEQUAL allows us to call MPI_Finalize */
492     bFinalize = (result != MPI_UNEQUAL);
493 #else
494     GMX_UNUSED_VALUE(comm);
495     bFinalize = TRUE;
496 #endif
497
498     va_start(ap, fmt);
499     gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
500     va_end(ap);
501 }