Divide default communicator from DD communicators
[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, "Using two step summing over %d groups of on average %.1f ranks\n\n", ng,
194                     (real)n / (real)ng);
195         }
196         if (nc->rank_intra > 0)
197         {
198             MPI_Comm_free(&nc->comm_inter);
199         }
200     }
201     else
202     {
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);
206         if (debug)
207         {
208             fprintf(debug,
209                     "In gmx_setup_nodecomm: not unsing separate inter- and intra-node "
210                     "communicators.\n");
211         }
212     }
213 #    endif
214 #else
215     /* tMPI runs only on a single node so just use the nodeid */
216     nc->rank_intra = cr->nodeid;
217 #endif
218 }
219
220 void gmx_barrier(MPI_Comm gmx_unused communicator)
221 {
222 #if !GMX_MPI
223     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_barrier");
224 #else
225     MPI_Barrier(communicator);
226 #endif
227 }
228
229 void gmx_bcast(int gmx_unused nbytes, void gmx_unused* b, MPI_Comm gmx_unused communicator)
230 {
231 #if !GMX_MPI
232     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_bcast");
233 #else
234     MPI_Bcast(b, nbytes, MPI_BYTE, 0, communicator);
235 #endif
236 }
237
238 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused* cr)
239 {
240 #if !GMX_MPI
241     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumd");
242 #else
243 #    if MPI_IN_PLACE_EXISTS
244     if (cr->nc.bUse)
245     {
246         if (cr->nc.rank_intra == 0)
247         {
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);
252         }
253         else
254         {
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);
258         }
259         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
260     }
261     else
262     {
263         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
264     }
265 #    else
266     int i;
267
268     if (nr > cr->mpb->dbuf_alloc)
269     {
270         cr->mpb->dbuf_alloc = nr;
271         srenew(cr->mpb->dbuf, cr->mpb->dbuf_alloc);
272     }
273     if (cr->nc.bUse)
274     {
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)
278         {
279             /* Sum with the buffers reversed */
280             MPI_Allreduce(cr->mpb->dbuf, r, nr, MPI_DOUBLE, MPI_SUM, cr->nc.comm_inter);
281         }
282         MPI_Bcast(r, nr, MPI_DOUBLE, 0, cr->nc.comm_intra);
283     }
284     else
285     {
286         MPI_Allreduce(r, cr->mpb->dbuf, nr, MPI_DOUBLE, MPI_SUM, cr->mpi_comm_mygroup);
287         for (i = 0; i < nr; i++)
288         {
289             r[i] = cr->mpb->dbuf[i];
290         }
291     }
292 #    endif
293 #endif
294 }
295
296 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused* cr)
297 {
298 #if !GMX_MPI
299     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumf");
300 #else
301 #    if MPI_IN_PLACE_EXISTS
302     if (cr->nc.bUse)
303     {
304         /* Use two step summing.  */
305         if (cr->nc.rank_intra == 0)
306         {
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);
310         }
311         else
312         {
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);
316         }
317         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
318     }
319     else
320     {
321         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
322     }
323 #    else
324     int i;
325
326     if (nr > cr->mpb->fbuf_alloc)
327     {
328         cr->mpb->fbuf_alloc = nr;
329         srenew(cr->mpb->fbuf, cr->mpb->fbuf_alloc);
330     }
331     if (cr->nc.bUse)
332     {
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)
336         {
337             /* Sum with the buffers reversed */
338             MPI_Allreduce(cr->mpb->fbuf, r, nr, MPI_FLOAT, MPI_SUM, cr->nc.comm_inter);
339         }
340         MPI_Bcast(r, nr, MPI_FLOAT, 0, cr->nc.comm_intra);
341     }
342     else
343     {
344         MPI_Allreduce(r, cr->mpb->fbuf, nr, MPI_FLOAT, MPI_SUM, cr->mpi_comm_mygroup);
345         for (i = 0; i < nr; i++)
346         {
347             r[i] = cr->mpb->fbuf[i];
348         }
349     }
350 #    endif
351 #endif
352 }
353
354 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused* cr)
355 {
356 #if !GMX_MPI
357     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumi");
358 #else
359 #    if MPI_IN_PLACE_EXISTS
360     if (cr->nc.bUse)
361     {
362         /* Use two step summing */
363         if (cr->nc.rank_intra == 0)
364         {
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);
368         }
369         else
370         {
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);
374         }
375         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
376     }
377     else
378     {
379         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
380     }
381 #    else
382     int i;
383
384     if (nr > cr->mpb->ibuf_alloc)
385     {
386         cr->mpb->ibuf_alloc = nr;
387         srenew(cr->mpb->ibuf, cr->mpb->ibuf_alloc);
388     }
389     if (cr->nc.bUse)
390     {
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)
394         {
395             /* Sum with the buffers reversed */
396             MPI_Allreduce(cr->mpb->ibuf, r, nr, MPI_INT, MPI_SUM, cr->nc.comm_inter);
397         }
398         MPI_Bcast(r, nr, MPI_INT, 0, cr->nc.comm_intra);
399     }
400     else
401     {
402         MPI_Allreduce(r, cr->mpb->ibuf, nr, MPI_INT, MPI_SUM, cr->mpi_comm_mygroup);
403         for (i = 0; i < nr; i++)
404         {
405             r[i] = cr->mpb->ibuf[i];
406         }
407     }
408 #    endif
409 #endif
410 }
411
412 void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused* cr)
413 {
414 #if !GMX_MPI
415     GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumli");
416 #else
417 #    if MPI_IN_PLACE_EXISTS
418     if (cr->nc.bUse)
419     {
420         /* Use two step summing */
421         if (cr->nc.rank_intra == 0)
422         {
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);
426         }
427         else
428         {
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);
432         }
433         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
434     }
435     else
436     {
437         MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
438     }
439 #    else
440     int i;
441
442     if (nr > cr->mpb->libuf_alloc)
443     {
444         cr->mpb->libuf_alloc = nr;
445         srenew(cr->mpb->libuf, cr->mpb->libuf_alloc);
446     }
447     if (cr->nc.bUse)
448     {
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)
452         {
453             /* Sum with the buffers reversed */
454             MPI_Allreduce(cr->mpb->libuf, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
455         }
456         MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
457     }
458     else
459     {
460         MPI_Allreduce(r, cr->mpb->libuf, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
461         for (i = 0; i < nr; i++)
462         {
463             r[i] = cr->mpb->libuf[i];
464         }
465     }
466 #    endif
467 #endif
468 }
469
470 const char* opt2fn_master(const char* opt, int nfile, const t_filenm fnm[], t_commrec* cr)
471 {
472     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
473 }
474
475 void gmx_fatal_collective(int                    f_errno,
476                           const char*            file,
477                           int                    line,
478                           MPI_Comm               comm,
479                           gmx_bool               bMaster,
480                           gmx_fmtstr const char* fmt,
481                           ...)
482 {
483     va_list  ap;
484     gmx_bool bFinalize;
485 #if GMX_MPI
486     int result;
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);
491 #else
492     GMX_UNUSED_VALUE(comm);
493     bFinalize = TRUE;
494 #endif
495
496     va_start(ap, fmt);
497     gmx_fatal_mpi_va(f_errno, file, line, bMaster, bFinalize, fmt, ap);
498     va_end(ap);
499 }