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