xref: /petsc/src/sys/objects/subcomm.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
1 
2 /*
3      Provides utility routines for split MPI communicator.
4 */
5 #include <petscsys.h> /*I   "petscsys.h"    I*/
6 #include <petscviewer.h>
7 
8 const char *const PetscSubcommTypes[] = {"GENERAL", "CONTIGUOUS", "INTERLACED", "PetscSubcommType", "PETSC_SUBCOMM_", NULL};
9 
10 static PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm);
11 static PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm);
12 
13 /*@
14    PetscSubcommSetFromOptions - Allows setting options for a `PetscSubcomm`
15 
16    Collective on psubcomm
17 
18    Input Parameter:
19 .  psubcomm - `PetscSubcomm` context
20 
21    Level: beginner
22 
23 .seealso: `PetscSubcomm`, `PetscSubcommCreate()`
24 @*/
25 PetscErrorCode PetscSubcommSetFromOptions(PetscSubcomm psubcomm)
26 {
27   PetscSubcommType type;
28   PetscBool        flg;
29 
30   PetscFunctionBegin;
31   PetscCheck(psubcomm, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "Must call PetscSubcommCreate firt");
32 
33   PetscOptionsBegin(psubcomm->parent, psubcomm->subcommprefix, "Options for PetscSubcomm", NULL);
34   PetscCall(PetscOptionsEnum("-psubcomm_type", NULL, NULL, PetscSubcommTypes, (PetscEnum)psubcomm->type, (PetscEnum *)&type, &flg));
35   if (flg && psubcomm->type != type) {
36     /* free old structures */
37     PetscCall(PetscCommDestroy(&(psubcomm)->dupparent));
38     PetscCall(PetscCommDestroy(&(psubcomm)->child));
39     PetscCall(PetscFree((psubcomm)->subsize));
40     switch (type) {
41     case PETSC_SUBCOMM_GENERAL:
42       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Runtime option PETSC_SUBCOMM_GENERAL is not supported, use PetscSubcommSetTypeGeneral()");
43     case PETSC_SUBCOMM_CONTIGUOUS:
44       PetscCall(PetscSubcommCreate_contiguous(psubcomm));
45       break;
46     case PETSC_SUBCOMM_INTERLACED:
47       PetscCall(PetscSubcommCreate_interlaced(psubcomm));
48       break;
49     default:
50       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "PetscSubcommType %s is not supported yet", PetscSubcommTypes[type]);
51     }
52   }
53 
54   PetscCall(PetscOptionsName("-psubcomm_view", "Triggers display of PetscSubcomm context", "PetscSubcommView", &flg));
55   if (flg) PetscCall(PetscSubcommView(psubcomm, PETSC_VIEWER_STDOUT_(psubcomm->parent)));
56   PetscOptionsEnd();
57   PetscFunctionReturn(0);
58 }
59 
60 /*@C
61   PetscSubcommSetOptionsPrefix - Sets the prefix used for searching for options in the options database for this object
62 
63   Logically collective on psubcomm
64 
65   Level: Intermediate
66 
67   Input Parameters:
68 +   psubcomm - `PetscSubcomm` context
69 -   prefix - the prefix to prepend all PetscSubcomm item names with.
70 
71 .seealso: `PetscSubcomm`, `PetscSubcommCreate()`
72 @*/
73 PetscErrorCode PetscSubcommSetOptionsPrefix(PetscSubcomm psubcomm, const char pre[])
74 {
75   PetscFunctionBegin;
76   if (!pre) {
77     PetscCall(PetscFree(psubcomm->subcommprefix));
78   } else {
79     PetscCheck(pre[0] != '-', PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Options prefix should not begin with a hyphen");
80     PetscCall(PetscFree(psubcomm->subcommprefix));
81     PetscCall(PetscStrallocpy(pre, &(psubcomm->subcommprefix)));
82   }
83   PetscFunctionReturn(0);
84 }
85 
86 /*@C
87    PetscSubcommView - Views a `PetscSubcomm`
88 
89    Collective on psubcomm
90 
91    Input Parameters:
92 +  psubcomm - `PetscSubcomm` context
93 -  viewer - `PetscViewer` to display the information
94 
95    Level: beginner
96 
97 .seealso: `PetscSubcomm`, `PetscSubcommCreate()`, `PetscViewer`
98 @*/
99 PetscErrorCode PetscSubcommView(PetscSubcomm psubcomm, PetscViewer viewer)
100 {
101   PetscBool         iascii;
102   PetscViewerFormat format;
103 
104   PetscFunctionBegin;
105   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
106   if (iascii) {
107     PetscCall(PetscViewerGetFormat(viewer, &format));
108     if (format == PETSC_VIEWER_DEFAULT) {
109       MPI_Comm    comm = psubcomm->parent;
110       PetscMPIInt rank, size, subsize, subrank, duprank;
111 
112       PetscCallMPI(MPI_Comm_size(comm, &size));
113       PetscCall(PetscViewerASCIIPrintf(viewer, "PetscSubcomm type %s with total %d MPI processes:\n", PetscSubcommTypes[psubcomm->type], size));
114       PetscCallMPI(MPI_Comm_rank(comm, &rank));
115       PetscCallMPI(MPI_Comm_size(psubcomm->child, &subsize));
116       PetscCallMPI(MPI_Comm_rank(psubcomm->child, &subrank));
117       PetscCallMPI(MPI_Comm_rank(psubcomm->dupparent, &duprank));
118       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
119       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d], color %d, sub-size %d, sub-rank %d, duprank %d\n", rank, psubcomm->color, subsize, subrank, duprank));
120       PetscCall(PetscViewerFlush(viewer));
121       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
122     }
123   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not supported yet");
124   PetscFunctionReturn(0);
125 }
126 
127 /*@
128   PetscSubcommSetNumber - Set total number of subcommunicators desired in the given `PetscSubcomm`
129 
130    Collective
131 
132    Input Parameters:
133 +  psubcomm - `PetscSubcomm` context
134 -  nsubcomm - the total number of subcommunicators in psubcomm
135 
136    Level: advanced
137 
138 .seealso: `PetscSubcomm`, `PetscSubcommCreate()`, `PetscSubcommDestroy()`, `PetscSubcommSetType()`, `PetscSubcommSetTypeGeneral()`
139 @*/
140 PetscErrorCode PetscSubcommSetNumber(PetscSubcomm psubcomm, PetscInt nsubcomm)
141 {
142   MPI_Comm    comm = psubcomm->parent;
143   PetscMPIInt msub, size;
144 
145   PetscFunctionBegin;
146   PetscCheck(psubcomm, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "PetscSubcomm is not created. Call PetscSubcommCreate() first");
147   PetscCallMPI(MPI_Comm_size(comm, &size));
148   PetscCall(PetscMPIIntCast(nsubcomm, &msub));
149   PetscCheck(msub >= 1 && msub <= size, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Num of subcommunicators %d cannot be < 1 or > input comm size %d", msub, size);
150 
151   psubcomm->n = msub;
152   PetscFunctionReturn(0);
153 }
154 
155 /*@
156   PetscSubcommSetType - Set the way the original MPI communicator is divided up in the `PetscSubcomm`
157 
158    Collective
159 
160    Input Parameters:
161 +  psubcomm - `PetscSubcomm` context
162 -  subcommtype - `PetscSubcommType` `PETSC_SUBCOMM_CONTIGUOUS` or `PETSC_SUBCOMM_INTERLACED`
163 
164    Level: advanced
165 
166 .seealso: `PetscSubcommType`, `PETSC_SUBCOMM_CONTIGUOUS`, `PETSC_SUBCOMM_INTERLACED`,
167           `PetscSubcommCreate()`, `PetscSubcommDestroy()`, `PetscSubcommSetNumber()`, `PetscSubcommSetTypeGeneral()`, `PetscSubcommType`
168 @*/
169 PetscErrorCode PetscSubcommSetType(PetscSubcomm psubcomm, PetscSubcommType subcommtype)
170 {
171   PetscFunctionBegin;
172   PetscCheck(psubcomm, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "PetscSubcomm is not created. Call PetscSubcommCreate()");
173   PetscCheck(psubcomm->n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()", psubcomm->n);
174 
175   if (subcommtype == PETSC_SUBCOMM_CONTIGUOUS) {
176     PetscCall(PetscSubcommCreate_contiguous(psubcomm));
177   } else if (subcommtype == PETSC_SUBCOMM_INTERLACED) {
178     PetscCall(PetscSubcommCreate_interlaced(psubcomm));
179   } else SETERRQ(psubcomm->parent, PETSC_ERR_SUP, "PetscSubcommType %s is not supported yet", PetscSubcommTypes[subcommtype]);
180   PetscFunctionReturn(0);
181 }
182 
183 /*@
184   PetscSubcommSetTypeGeneral - Divides up a communicator based on a specific user's specification
185 
186    Collective
187 
188    Input Parameters:
189 +  psubcomm - `PetscSubcomm` context
190 .  color   - control of subset assignment (nonnegative integer). Processes with the same color are in the same subcommunicator.
191 -  subrank - rank in the subcommunicator
192 
193    Level: advanced
194 
195 .seealso: `PetscSubcommType`, `PETSC_SUBCOMM_CONTIGUOUS`, `PETSC_SUBCOMM_INTERLACED`, `PetscSubcommCreate()`, `PetscSubcommDestroy()`, `PetscSubcommSetNumber()`, `PetscSubcommSetType()`
196 @*/
197 PetscErrorCode PetscSubcommSetTypeGeneral(PetscSubcomm psubcomm, PetscMPIInt color, PetscMPIInt subrank)
198 {
199   MPI_Comm    subcomm = 0, dupcomm = 0, comm = psubcomm->parent;
200   PetscMPIInt size, icolor, duprank, *recvbuf, sendbuf[3], mysubsize, rank, *subsize;
201   PetscMPIInt i, nsubcomm = psubcomm->n;
202 
203   PetscFunctionBegin;
204   PetscCheck(psubcomm, PETSC_COMM_SELF, PETSC_ERR_ARG_NULL, "PetscSubcomm is not created. Call PetscSubcommCreate()");
205   PetscCheck(nsubcomm >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of subcommunicators %d is incorrect. Call PetscSubcommSetNumber()", nsubcomm);
206 
207   PetscCallMPI(MPI_Comm_split(comm, color, subrank, &subcomm));
208 
209   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
210   /* TODO: this can be done in an ostensibly scalale way (i.e., without allocating an array of size 'size') as is done in PetscObjectsCreateGlobalOrdering(). */
211   PetscCallMPI(MPI_Comm_size(comm, &size));
212   PetscCall(PetscMalloc1(2 * size, &recvbuf));
213 
214   PetscCallMPI(MPI_Comm_rank(comm, &rank));
215   PetscCallMPI(MPI_Comm_size(subcomm, &mysubsize));
216 
217   sendbuf[0] = color;
218   sendbuf[1] = mysubsize;
219   PetscCallMPI(MPI_Allgather(sendbuf, 2, MPI_INT, recvbuf, 2, MPI_INT, comm));
220 
221   PetscCall(PetscCalloc1(nsubcomm, &subsize));
222   for (i = 0; i < 2 * size; i += 2) subsize[recvbuf[i]] = recvbuf[i + 1];
223   PetscCall(PetscFree(recvbuf));
224 
225   duprank = 0;
226   for (icolor = 0; icolor < nsubcomm; icolor++) {
227     if (icolor != color) { /* not color of this process */
228       duprank += subsize[icolor];
229     } else {
230       duprank += subrank;
231       break;
232     }
233   }
234   PetscCallMPI(MPI_Comm_split(comm, 0, duprank, &dupcomm));
235 
236   PetscCall(PetscCommDuplicate(dupcomm, &psubcomm->dupparent, NULL));
237   PetscCall(PetscCommDuplicate(subcomm, &psubcomm->child, NULL));
238   PetscCallMPI(MPI_Comm_free(&dupcomm));
239   PetscCallMPI(MPI_Comm_free(&subcomm));
240 
241   psubcomm->color   = color;
242   psubcomm->subsize = subsize;
243   psubcomm->type    = PETSC_SUBCOMM_GENERAL;
244   PetscFunctionReturn(0);
245 }
246 
247 /*@
248   PetscSubcommDestroy - Destroys a `PetscSubcomm` object
249 
250    Collective on psubcomm
251 
252    Input Parameter:
253    .  psubcomm - the `PetscSubcomm` context
254 
255    Level: advanced
256 
257 .seealso: `PetscSubcommCreate()`, `PetscSubcommSetType()`
258 @*/
259 PetscErrorCode PetscSubcommDestroy(PetscSubcomm *psubcomm)
260 {
261   PetscFunctionBegin;
262   if (!*psubcomm) PetscFunctionReturn(0);
263   PetscCall(PetscCommDestroy(&(*psubcomm)->dupparent));
264   PetscCall(PetscCommDestroy(&(*psubcomm)->child));
265   PetscCall(PetscFree((*psubcomm)->subsize));
266   if ((*psubcomm)->subcommprefix) PetscCall(PetscFree((*psubcomm)->subcommprefix));
267   PetscCall(PetscFree((*psubcomm)));
268   PetscFunctionReturn(0);
269 }
270 
271 /*@
272   PetscSubcommCreate - Create a `PetscSubcomm` context. This object is used to manage the division of a `MPI_Comm` into subcommunicators
273 
274    Collective
275 
276    Input Parameter:
277 .  comm - MPI communicator
278 
279    Output Parameter:
280 .  psubcomm - location to store the `PetscSubcomm` context
281 
282    Level: advanced
283 
284 .seealso: `PetscSubcomm`, `PetscSubcommDestroy()`, `PetscSubcommSetTypeGeneral()`, `PetscSubcommSetFromOptions()`, `PetscSubcommSetType()`,
285           `PetscSubcommSetNumber()`
286 @*/
287 PetscErrorCode PetscSubcommCreate(MPI_Comm comm, PetscSubcomm *psubcomm)
288 {
289   PetscMPIInt rank, size;
290 
291   PetscFunctionBegin;
292   PetscCall(PetscNew(psubcomm));
293 
294   /* set defaults */
295   PetscCallMPI(MPI_Comm_rank(comm, &rank));
296   PetscCallMPI(MPI_Comm_size(comm, &size));
297 
298   (*psubcomm)->parent    = comm;
299   (*psubcomm)->dupparent = comm;
300   (*psubcomm)->child     = PETSC_COMM_SELF;
301   (*psubcomm)->n         = size;
302   (*psubcomm)->color     = rank;
303   (*psubcomm)->subsize   = NULL;
304   (*psubcomm)->type      = PETSC_SUBCOMM_INTERLACED;
305   PetscFunctionReturn(0);
306 }
307 
308 /*@C
309   PetscSubcommGetParent - Gets the communicator that was used to create the `PetscSubcomm`
310 
311    Collective
312 
313    Input Parameter:
314 .  scomm - the `PetscSubcomm`
315 
316    Output Parameter:
317 .  pcomm - location to store the parent communicator
318 
319    Level: intermediate
320 
321 .seealso: `PetscSubcommDestroy()`, `PetscSubcommSetTypeGeneral()`, `PetscSubcommSetFromOptions()`, `PetscSubcommSetType()`,
322           `PetscSubcommSetNumber()`, `PetscSubcommGetChild()`, `PetscSubcommContiguousParent()`
323 @*/
324 PetscErrorCode PetscSubcommGetParent(PetscSubcomm scomm, MPI_Comm *pcomm)
325 {
326   *pcomm = PetscSubcommParent(scomm);
327   return 0;
328 }
329 
330 /*@C
331   PetscSubcommGetContiguousParent - Gets a communicator that that is a duplicate of the parent but has the ranks
332                                     reordered by the order they are in the children
333 
334    Collective
335 
336    Input Parameter:
337 .  scomm - the `PetscSubcomm`
338 
339    Output Parameter:
340 .  pcomm - location to store the parent communicator
341 
342    Level: intermediate
343 
344 .seealso: `PetscSubcommDestroy()`, `PetscSubcommSetTypeGeneral()`, `PetscSubcommSetFromOptions()`, `PetscSubcommSetType()`,
345           `PetscSubcommSetNumber()`, `PetscSubcommGetChild()`, `PetscSubcommContiguousParent()`
346 @*/
347 PetscErrorCode PetscSubcommGetContiguousParent(PetscSubcomm scomm, MPI_Comm *pcomm)
348 {
349   *pcomm = PetscSubcommContiguousParent(scomm);
350   return 0;
351 }
352 
353 /*@C
354   PetscSubcommGetChild - Gets the communicator created by the `PetscSubcomm`. This is part of one of the subcommunicators created by the `PetscSubcomm`
355 
356    Collective
357 
358    Input Parameter:
359 .  scomm - the `PetscSubcomm`
360 
361    Output Parameter:
362 .  ccomm - location to store the child communicator
363 
364    Level: intermediate
365 
366 .seealso: `PetscSubcommDestroy()`, `PetscSubcommSetTypeGeneral()`, `PetscSubcommSetFromOptions()`, `PetscSubcommSetType()`,
367           `PetscSubcommSetNumber()`, `PetscSubcommGetParent()`, `PetscSubcommContiguousParent()`
368 @*/
369 PetscErrorCode PetscSubcommGetChild(PetscSubcomm scomm, MPI_Comm *ccomm)
370 {
371   *ccomm = PetscSubcommChild(scomm);
372   return 0;
373 }
374 
375 static PetscErrorCode PetscSubcommCreate_contiguous(PetscSubcomm psubcomm)
376 {
377   PetscMPIInt rank, size, *subsize, duprank = -1, subrank = -1;
378   PetscMPIInt np_subcomm, nleftover, i, color = -1, rankstart, nsubcomm = psubcomm->n;
379   MPI_Comm    subcomm = 0, dupcomm = 0, comm = psubcomm->parent;
380 
381   PetscFunctionBegin;
382   PetscCallMPI(MPI_Comm_rank(comm, &rank));
383   PetscCallMPI(MPI_Comm_size(comm, &size));
384 
385   /* get size of each subcommunicator */
386   PetscCall(PetscMalloc1(1 + nsubcomm, &subsize));
387 
388   np_subcomm = size / nsubcomm;
389   nleftover  = size - nsubcomm * np_subcomm;
390   for (i = 0; i < nsubcomm; i++) {
391     subsize[i] = np_subcomm;
392     if (i < nleftover) subsize[i]++;
393   }
394 
395   /* get color and subrank of this proc */
396   rankstart = 0;
397   for (i = 0; i < nsubcomm; i++) {
398     if (rank >= rankstart && rank < rankstart + subsize[i]) {
399       color   = i;
400       subrank = rank - rankstart;
401       duprank = rank;
402       break;
403     } else rankstart += subsize[i];
404   }
405 
406   PetscCallMPI(MPI_Comm_split(comm, color, subrank, &subcomm));
407 
408   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
409   PetscCallMPI(MPI_Comm_split(comm, 0, duprank, &dupcomm));
410   PetscCall(PetscCommDuplicate(dupcomm, &psubcomm->dupparent, NULL));
411   PetscCall(PetscCommDuplicate(subcomm, &psubcomm->child, NULL));
412   PetscCallMPI(MPI_Comm_free(&dupcomm));
413   PetscCallMPI(MPI_Comm_free(&subcomm));
414 
415   psubcomm->color   = color;
416   psubcomm->subsize = subsize;
417   psubcomm->type    = PETSC_SUBCOMM_CONTIGUOUS;
418   PetscFunctionReturn(0);
419 }
420 
421 /*
422    Note:
423    In PCREDUNDANT, to avoid data scattering from subcomm back to original comm, we create subcommunicators
424    by iteratively taking a process into a subcommunicator.
425    Example: size=4, nsubcomm=(*psubcomm)->n=3
426      comm=(*psubcomm)->parent:
427       rank:     [0]  [1]  [2]  [3]
428       color:     0    1    2    0
429 
430      subcomm=(*psubcomm)->comm:
431       subrank:  [0]  [0]  [0]  [1]
432 
433      dupcomm=(*psubcomm)->dupparent:
434       duprank:  [0]  [2]  [3]  [1]
435 
436      Here, subcomm[color = 0] has subsize=2, owns process [0] and [3]
437            subcomm[color = 1] has subsize=1, owns process [1]
438            subcomm[color = 2] has subsize=1, owns process [2]
439            dupcomm has same number of processes as comm, and its duprank maps
440            processes in subcomm contiguously into a 1d array:
441             duprank: [0] [1]      [2]         [3]
442             rank:    [0] [3]      [1]         [2]
443                     subcomm[0] subcomm[1]  subcomm[2]
444 */
445 
446 static PetscErrorCode PetscSubcommCreate_interlaced(PetscSubcomm psubcomm)
447 {
448   PetscMPIInt rank, size, *subsize, duprank, subrank;
449   PetscMPIInt np_subcomm, nleftover, i, j, color, nsubcomm = psubcomm->n;
450   MPI_Comm    subcomm = 0, dupcomm = 0, comm = psubcomm->parent;
451 
452   PetscFunctionBegin;
453   PetscCallMPI(MPI_Comm_rank(comm, &rank));
454   PetscCallMPI(MPI_Comm_size(comm, &size));
455 
456   /* get size of each subcommunicator */
457   PetscCall(PetscMalloc1(1 + nsubcomm, &subsize));
458 
459   np_subcomm = size / nsubcomm;
460   nleftover  = size - nsubcomm * np_subcomm;
461   for (i = 0; i < nsubcomm; i++) {
462     subsize[i] = np_subcomm;
463     if (i < nleftover) subsize[i]++;
464   }
465 
466   /* find color for this proc */
467   color   = rank % nsubcomm;
468   subrank = rank / nsubcomm;
469 
470   PetscCallMPI(MPI_Comm_split(comm, color, subrank, &subcomm));
471 
472   j       = 0;
473   duprank = 0;
474   for (i = 0; i < nsubcomm; i++) {
475     if (j == color) {
476       duprank += subrank;
477       break;
478     }
479     duprank += subsize[i];
480     j++;
481   }
482 
483   /* create dupcomm with same size as comm, but its rank, duprank, maps subcomm's contiguously into dupcomm */
484   PetscCallMPI(MPI_Comm_split(comm, 0, duprank, &dupcomm));
485   PetscCall(PetscCommDuplicate(dupcomm, &psubcomm->dupparent, NULL));
486   PetscCall(PetscCommDuplicate(subcomm, &psubcomm->child, NULL));
487   PetscCallMPI(MPI_Comm_free(&dupcomm));
488   PetscCallMPI(MPI_Comm_free(&subcomm));
489 
490   psubcomm->color   = color;
491   psubcomm->subsize = subsize;
492   psubcomm->type    = PETSC_SUBCOMM_INTERLACED;
493   PetscFunctionReturn(0);
494 }
495