Actual source code: iscoloring.c
1: #include <petsc/private/isimpl.h>
2: #include <petscviewer.h>
3: #include <petscsf.h>
5: const char *const ISColoringTypes[] = {"global", "ghosted", "ISColoringType", "IS_COLORING_", NULL};
7: PetscErrorCode ISColoringReference(ISColoring coloring)
8: {
9: PetscFunctionBegin;
10: coloring->refct++;
11: PetscFunctionReturn(PETSC_SUCCESS);
12: }
14: /*@
15: ISColoringSetType - indicates if the coloring is for the local representation (including ghost points) or the global representation of a `Mat`
17: Collective
19: Input Parameters:
20: + coloring - the coloring object
21: - type - either `IS_COLORING_LOCAL` or `IS_COLORING_GLOBAL`
23: Level: intermediate
25: Notes:
26: `IS_COLORING_LOCAL` can lead to faster computations since parallel ghost point updates are not needed for each color
28: With `IS_COLORING_LOCAL` the coloring is in the numbering of the local vector, for `IS_COLORING_GLOBAL` it is in the numbering of the global vector
30: .seealso: `MatFDColoringCreate()`, `ISColoring`, `ISColoringType`, `ISColoringCreate()`, `IS_COLORING_LOCAL`, `IS_COLORING_GLOBAL`, `ISColoringGetType()`
31: @*/
32: PetscErrorCode ISColoringSetType(ISColoring coloring, ISColoringType type)
33: {
34: PetscFunctionBegin;
35: coloring->ctype = type;
36: PetscFunctionReturn(PETSC_SUCCESS);
37: }
39: /*@
40: ISColoringGetType - gets if the coloring is for the local representation (including ghost points) or the global representation
42: Collective
44: Input Parameter:
45: . coloring - the coloring object
47: Output Parameter:
48: . type - either `IS_COLORING_LOCAL` or `IS_COLORING_GLOBAL`
50: Level: intermediate
52: .seealso: `MatFDColoringCreate()`, `ISColoring`, `ISColoringType`, `ISColoringCreate()`, `IS_COLORING_LOCAL`, `IS_COLORING_GLOBAL`, `ISColoringSetType()`
53: @*/
54: PetscErrorCode ISColoringGetType(ISColoring coloring, ISColoringType *type)
55: {
56: PetscFunctionBegin;
57: *type = coloring->ctype;
58: PetscFunctionReturn(PETSC_SUCCESS);
59: }
61: /*@
62: ISColoringDestroy - Destroys an `ISColoring` coloring context.
64: Collective
66: Input Parameter:
67: . iscoloring - the coloring context
69: Level: advanced
71: .seealso: `ISColoring`, `ISColoringView()`, `MatColoring`
72: @*/
73: PetscErrorCode ISColoringDestroy(ISColoring *iscoloring)
74: {
75: PetscInt i;
77: PetscFunctionBegin;
78: if (!*iscoloring) PetscFunctionReturn(PETSC_SUCCESS);
79: PetscAssertPointer(*iscoloring, 1);
80: if (--(*iscoloring)->refct > 0) {
81: *iscoloring = NULL;
82: PetscFunctionReturn(PETSC_SUCCESS);
83: }
85: if ((*iscoloring)->is) {
86: for (i = 0; i < (*iscoloring)->n; i++) PetscCall(ISDestroy(&(*iscoloring)->is[i]));
87: PetscCall(PetscFree((*iscoloring)->is));
88: }
89: if ((*iscoloring)->allocated) PetscCall(PetscFree((*iscoloring)->colors));
90: PetscCall(PetscCommDestroy(&(*iscoloring)->comm));
91: PetscCall(PetscFree(*iscoloring));
92: PetscFunctionReturn(PETSC_SUCCESS);
93: }
95: /*@C
96: ISColoringViewFromOptions - Processes command line options to determine if/how an `ISColoring` object is to be viewed.
98: Collective
100: Input Parameters:
101: + obj - the `ISColoring` object
102: . bobj - prefix to use for viewing, or `NULL` to use prefix of `mat`
103: - optionname - option to activate viewing
105: Level: intermediate
107: Developer Notes:
108: This cannot use `PetscObjectViewFromOptions()` because `ISColoring` is not a `PetscObject`
110: .seealso: `ISColoring`, `ISColoringView()`
111: @*/
112: PetscErrorCode ISColoringViewFromOptions(ISColoring obj, PetscObject bobj, const char optionname[])
113: {
114: PetscViewer viewer;
115: PetscBool flg;
116: PetscViewerFormat format;
117: char *prefix;
119: PetscFunctionBegin;
120: prefix = bobj ? bobj->prefix : NULL;
121: PetscCall(PetscOptionsGetViewer(obj->comm, NULL, prefix, optionname, &viewer, &format, &flg));
122: if (flg) {
123: PetscCall(PetscViewerPushFormat(viewer, format));
124: PetscCall(ISColoringView(obj, viewer));
125: PetscCall(PetscViewerPopFormat(viewer));
126: PetscCall(PetscOptionsRestoreViewer(&viewer));
127: }
128: PetscFunctionReturn(PETSC_SUCCESS);
129: }
131: /*@C
132: ISColoringView - Views an `ISColoring` coloring context.
134: Collective
136: Input Parameters:
137: + iscoloring - the coloring context
138: - viewer - the viewer
140: Level: advanced
142: .seealso: `ISColoring()`, `ISColoringViewFromOptions()`, `ISColoringDestroy()`, `ISColoringGetIS()`, `MatColoring`
143: @*/
144: PetscErrorCode ISColoringView(ISColoring iscoloring, PetscViewer viewer)
145: {
146: PetscInt i;
147: PetscBool iascii;
148: IS *is;
150: PetscFunctionBegin;
151: PetscAssertPointer(iscoloring, 1);
152: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(iscoloring->comm, &viewer));
155: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
156: if (iascii) {
157: MPI_Comm comm;
158: PetscMPIInt size, rank;
160: PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm));
161: PetscCallMPI(MPI_Comm_size(comm, &size));
162: PetscCallMPI(MPI_Comm_rank(comm, &rank));
163: PetscCall(PetscViewerASCIIPrintf(viewer, "ISColoring Object: %d MPI processes\n", size));
164: PetscCall(PetscViewerASCIIPrintf(viewer, "ISColoringType: %s\n", ISColoringTypes[iscoloring->ctype]));
165: PetscCall(PetscViewerASCIIPushSynchronized(viewer));
166: PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] Number of colors %" PetscInt_FMT "\n", rank, iscoloring->n));
167: PetscCall(PetscViewerFlush(viewer));
168: PetscCall(PetscViewerASCIIPopSynchronized(viewer));
169: }
171: PetscCall(ISColoringGetIS(iscoloring, PETSC_USE_POINTER, PETSC_IGNORE, &is));
172: for (i = 0; i < iscoloring->n; i++) PetscCall(ISView(iscoloring->is[i], viewer));
173: PetscCall(ISColoringRestoreIS(iscoloring, PETSC_USE_POINTER, &is));
174: PetscFunctionReturn(PETSC_SUCCESS);
175: }
177: /*@C
178: ISColoringGetColors - Returns an array with the color for each local node
180: Not Collective
182: Input Parameter:
183: . iscoloring - the coloring context
185: Output Parameters:
186: + n - number of nodes
187: . nc - number of colors
188: - colors - color for each node
190: Level: advanced
192: Notes:
193: Do not free the `colors` array.
195: The `colors` array will only be valid for the lifetime of the `ISColoring`
197: .seealso: `ISColoring`, `ISColoringValue`, `ISColoringRestoreIS()`, `ISColoringView()`, `ISColoringGetIS()`
198: @*/
199: PetscErrorCode ISColoringGetColors(ISColoring iscoloring, PetscInt *n, PetscInt *nc, const ISColoringValue **colors)
200: {
201: PetscFunctionBegin;
202: PetscAssertPointer(iscoloring, 1);
204: if (n) *n = iscoloring->N;
205: if (nc) *nc = iscoloring->n;
206: if (colors) *colors = iscoloring->colors;
207: PetscFunctionReturn(PETSC_SUCCESS);
208: }
210: /*@C
211: ISColoringGetIS - Extracts index sets from the coloring context. Each is contains the nodes of one color
213: Collective
215: Input Parameters:
216: + iscoloring - the coloring context
217: - mode - if this value is `PETSC_OWN_POINTER` then the caller owns the pointer and must free the array of `IS` and each `IS` in the array
219: Output Parameters:
220: + nn - number of index sets in the coloring context
221: - isis - array of index sets
223: Level: advanced
225: Note:
226: If mode is `PETSC_USE_POINTER` then `ISColoringRestoreIS()` must be called when the `IS` are no longer needed
228: .seealso: `ISColoring`, `IS`, `ISColoringRestoreIS()`, `ISColoringView()`, `ISColoringGetColoring()`, `ISColoringGetColors()`
229: @*/
230: PetscErrorCode ISColoringGetIS(ISColoring iscoloring, PetscCopyMode mode, PetscInt *nn, IS *isis[])
231: {
232: PetscFunctionBegin;
233: PetscAssertPointer(iscoloring, 1);
235: if (nn) *nn = iscoloring->n;
236: if (isis) {
237: if (!iscoloring->is) {
238: PetscInt *mcolors, **ii, nc = iscoloring->n, i, base, n = iscoloring->N;
239: ISColoringValue *colors = iscoloring->colors;
240: IS *is;
242: if (PetscDefined(USE_DEBUG)) {
243: for (i = 0; i < n; i++) PetscCheck(((PetscInt)colors[i]) < nc, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Coloring is our of range index %d value %d number colors %d", (int)i, (int)colors[i], (int)nc);
244: }
246: /* generate the lists of nodes for each color */
247: PetscCall(PetscCalloc1(nc, &mcolors));
248: for (i = 0; i < n; i++) mcolors[colors[i]]++;
250: PetscCall(PetscMalloc1(nc, &ii));
251: PetscCall(PetscMalloc1(n, &ii[0]));
252: for (i = 1; i < nc; i++) ii[i] = ii[i - 1] + mcolors[i - 1];
253: PetscCall(PetscArrayzero(mcolors, nc));
255: if (iscoloring->ctype == IS_COLORING_GLOBAL) {
256: PetscCallMPI(MPI_Scan(&iscoloring->N, &base, 1, MPIU_INT, MPI_SUM, iscoloring->comm));
257: base -= iscoloring->N;
258: for (i = 0; i < n; i++) ii[colors[i]][mcolors[colors[i]]++] = i + base; /* global idx */
259: } else if (iscoloring->ctype == IS_COLORING_LOCAL) {
260: for (i = 0; i < n; i++) ii[colors[i]][mcolors[colors[i]]++] = i; /* local idx */
261: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not provided for this ISColoringType type");
263: PetscCall(PetscMalloc1(nc, &is));
264: for (i = 0; i < nc; i++) PetscCall(ISCreateGeneral(iscoloring->comm, mcolors[i], ii[i], PETSC_COPY_VALUES, is + i));
266: if (mode != PETSC_OWN_POINTER) iscoloring->is = is;
267: *isis = is;
268: PetscCall(PetscFree(ii[0]));
269: PetscCall(PetscFree(ii));
270: PetscCall(PetscFree(mcolors));
271: } else {
272: *isis = iscoloring->is;
273: }
274: }
275: PetscFunctionReturn(PETSC_SUCCESS);
276: }
278: /*@C
279: ISColoringRestoreIS - Restores the index sets extracted from the coloring context with `ISColoringGetIS()` using `PETSC_USE_POINTER`
281: Collective
283: Input Parameters:
284: + iscoloring - the coloring context
285: . mode - who retains ownership of the is
286: - is - array of index sets
288: Level: advanced
290: .seealso: `ISColoring()`, `IS`, `ISColoringGetIS()`, `ISColoringView()`, `PetscCopyMode`
291: @*/
292: PetscErrorCode ISColoringRestoreIS(ISColoring iscoloring, PetscCopyMode mode, IS *is[])
293: {
294: PetscFunctionBegin;
295: PetscAssertPointer(iscoloring, 1);
297: /* currently nothing is done here */
298: PetscFunctionReturn(PETSC_SUCCESS);
299: }
301: /*@
302: ISColoringCreate - Generates an `ISColoring` context from lists (provided by each MPI process) of colors for each node.
304: Collective
306: Input Parameters:
307: + comm - communicator for the processors creating the coloring
308: . ncolors - max color value
309: . n - number of nodes on this processor
310: . colors - array containing the colors for this MPI rank, color numbers begin at 0, for each local node
311: - mode - see `PetscCopyMode` for meaning of this flag.
313: Output Parameter:
314: . iscoloring - the resulting coloring data structure
316: Options Database Key:
317: . -is_coloring_view - Activates `ISColoringView()`
319: Level: advanced
321: Notes:
322: By default sets coloring type to `IS_COLORING_GLOBAL`
324: .seealso: `ISColoring`, `ISColoringValue`, `MatColoringCreate()`, `ISColoringView()`, `ISColoringDestroy()`, `ISColoringSetType()`
325: @*/
326: PetscErrorCode ISColoringCreate(MPI_Comm comm, PetscInt ncolors, PetscInt n, const ISColoringValue colors[], PetscCopyMode mode, ISColoring *iscoloring)
327: {
328: PetscMPIInt size, rank, tag;
329: PetscInt base, top, i;
330: PetscInt nc, ncwork;
331: MPI_Status status;
333: PetscFunctionBegin;
334: if (ncolors != PETSC_DECIDE && ncolors > IS_COLORING_MAX) {
335: PetscCheck(ncolors <= PETSC_MAX_UINT16, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Max color value exceeds %d limit. This number is unrealistic. Perhaps a bug in code? Current max: %d user requested: %" PetscInt_FMT, PETSC_MAX_UINT16, PETSC_IS_COLORING_MAX, ncolors);
336: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Max color value exceeds limit. Perhaps reconfigure PETSc with --with-is-color-value-type=short? Current max: %d user requested: %" PetscInt_FMT, PETSC_IS_COLORING_MAX, ncolors);
337: }
338: PetscCall(PetscNew(iscoloring));
339: PetscCall(PetscCommDuplicate(comm, &(*iscoloring)->comm, &tag));
340: comm = (*iscoloring)->comm;
342: /* compute the number of the first node on my processor */
343: PetscCallMPI(MPI_Comm_size(comm, &size));
345: /* should use MPI_Scan() */
346: PetscCallMPI(MPI_Comm_rank(comm, &rank));
347: if (rank == 0) {
348: base = 0;
349: top = n;
350: } else {
351: PetscCallMPI(MPI_Recv(&base, 1, MPIU_INT, rank - 1, tag, comm, &status));
352: top = base + n;
353: }
354: if (rank < size - 1) PetscCallMPI(MPI_Send(&top, 1, MPIU_INT, rank + 1, tag, comm));
356: /* compute the total number of colors */
357: ncwork = 0;
358: for (i = 0; i < n; i++) {
359: if (ncwork < colors[i]) ncwork = colors[i];
360: }
361: ncwork++;
362: PetscCall(MPIU_Allreduce(&ncwork, &nc, 1, MPIU_INT, MPI_MAX, comm));
363: PetscCheck(nc <= ncolors, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Number of colors passed in %" PetscInt_FMT " is less than the actual number of colors in array %" PetscInt_FMT, ncolors, nc);
364: (*iscoloring)->n = nc;
365: (*iscoloring)->is = NULL;
366: (*iscoloring)->N = n;
367: (*iscoloring)->refct = 1;
368: (*iscoloring)->ctype = IS_COLORING_GLOBAL;
369: if (mode == PETSC_COPY_VALUES) {
370: PetscCall(PetscMalloc1(n, &(*iscoloring)->colors));
371: PetscCall(PetscArraycpy((*iscoloring)->colors, colors, n));
372: (*iscoloring)->allocated = PETSC_TRUE;
373: } else if (mode == PETSC_OWN_POINTER) {
374: (*iscoloring)->colors = (ISColoringValue *)colors;
375: (*iscoloring)->allocated = PETSC_TRUE;
376: } else {
377: (*iscoloring)->colors = (ISColoringValue *)colors;
378: (*iscoloring)->allocated = PETSC_FALSE;
379: }
380: PetscCall(ISColoringViewFromOptions(*iscoloring, NULL, "-is_coloring_view"));
381: PetscCall(PetscInfo(0, "Number of colors %" PetscInt_FMT "\n", nc));
382: PetscFunctionReturn(PETSC_SUCCESS);
383: }
385: /*@
386: ISBuildTwoSided - Takes an `IS` that describes where each element will be mapped globally over all ranks.
387: Generates an `IS` that contains new numbers from remote or local on the `IS`.
389: Collective
391: Input Parameters:
392: + ito - an `IS` describes where each entry will be mapped. Negative target rank will be ignored
393: - toindx - an `IS` describes what indices should send. `NULL` means sending natural numbering
395: Output Parameter:
396: . rows - contains new numbers from remote or local
398: Level: advanced
400: Developer Notes:
401: This manual page is incomprehensible and still needs to be fixed
403: .seealso: [](sec_scatter), `IS`, `MatPartitioningCreate()`, `ISPartitioningToNumbering()`, `ISPartitioningCount()`
404: @*/
405: PetscErrorCode ISBuildTwoSided(IS ito, IS toindx, IS *rows)
406: {
407: const PetscInt *ito_indices, *toindx_indices;
408: PetscInt *send_indices, rstart, *recv_indices, nrecvs, nsends;
409: PetscInt *tosizes, *fromsizes, i, j, *tosizes_tmp, *tooffsets_tmp, ito_ln;
410: PetscMPIInt *toranks, *fromranks, size, target_rank, *fromperm_newtoold, nto, nfrom;
411: PetscLayout isrmap;
412: MPI_Comm comm;
413: PetscSF sf;
414: PetscSFNode *iremote;
416: PetscFunctionBegin;
417: PetscCall(PetscObjectGetComm((PetscObject)ito, &comm));
418: PetscCallMPI(MPI_Comm_size(comm, &size));
419: PetscCall(ISGetLocalSize(ito, &ito_ln));
420: PetscCall(ISGetLayout(ito, &isrmap));
421: PetscCall(PetscLayoutGetRange(isrmap, &rstart, NULL));
422: PetscCall(ISGetIndices(ito, &ito_indices));
423: PetscCall(PetscCalloc2(size, &tosizes_tmp, size + 1, &tooffsets_tmp));
424: for (i = 0; i < ito_ln; i++) {
425: if (ito_indices[i] < 0) continue;
426: else PetscCheck(ito_indices[i] < size, comm, PETSC_ERR_ARG_OUTOFRANGE, "target rank %" PetscInt_FMT " is larger than communicator size %d ", ito_indices[i], size);
427: tosizes_tmp[ito_indices[i]]++;
428: }
429: nto = 0;
430: for (i = 0; i < size; i++) {
431: tooffsets_tmp[i + 1] = tooffsets_tmp[i] + tosizes_tmp[i];
432: if (tosizes_tmp[i] > 0) nto++;
433: }
434: PetscCall(PetscCalloc2(nto, &toranks, 2 * nto, &tosizes));
435: nto = 0;
436: for (i = 0; i < size; i++) {
437: if (tosizes_tmp[i] > 0) {
438: toranks[nto] = i;
439: tosizes[2 * nto] = tosizes_tmp[i]; /* size */
440: tosizes[2 * nto + 1] = tooffsets_tmp[i]; /* offset */
441: nto++;
442: }
443: }
444: nsends = tooffsets_tmp[size];
445: PetscCall(PetscCalloc1(nsends, &send_indices));
446: if (toindx) PetscCall(ISGetIndices(toindx, &toindx_indices));
447: for (i = 0; i < ito_ln; i++) {
448: if (ito_indices[i] < 0) continue;
449: target_rank = ito_indices[i];
450: send_indices[tooffsets_tmp[target_rank]] = toindx ? toindx_indices[i] : (i + rstart);
451: tooffsets_tmp[target_rank]++;
452: }
453: if (toindx) PetscCall(ISRestoreIndices(toindx, &toindx_indices));
454: PetscCall(ISRestoreIndices(ito, &ito_indices));
455: PetscCall(PetscFree2(tosizes_tmp, tooffsets_tmp));
456: PetscCall(PetscCommBuildTwoSided(comm, 2, MPIU_INT, nto, toranks, tosizes, &nfrom, &fromranks, &fromsizes));
457: PetscCall(PetscFree2(toranks, tosizes));
458: PetscCall(PetscMalloc1(nfrom, &fromperm_newtoold));
459: for (i = 0; i < nfrom; i++) fromperm_newtoold[i] = i;
460: PetscCall(PetscSortMPIIntWithArray(nfrom, fromranks, fromperm_newtoold));
461: nrecvs = 0;
462: for (i = 0; i < nfrom; i++) nrecvs += fromsizes[i * 2];
463: PetscCall(PetscCalloc1(nrecvs, &recv_indices));
464: PetscCall(PetscMalloc1(nrecvs, &iremote));
465: nrecvs = 0;
466: for (i = 0; i < nfrom; i++) {
467: for (j = 0; j < fromsizes[2 * fromperm_newtoold[i]]; j++) {
468: iremote[nrecvs].rank = fromranks[i];
469: iremote[nrecvs++].index = fromsizes[2 * fromperm_newtoold[i] + 1] + j;
470: }
471: }
472: PetscCall(PetscSFCreate(comm, &sf));
473: PetscCall(PetscSFSetGraph(sf, nsends, nrecvs, NULL, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
474: PetscCall(PetscSFSetType(sf, PETSCSFBASIC));
475: /* how to put a prefix ? */
476: PetscCall(PetscSFSetFromOptions(sf));
477: PetscCall(PetscSFBcastBegin(sf, MPIU_INT, send_indices, recv_indices, MPI_REPLACE));
478: PetscCall(PetscSFBcastEnd(sf, MPIU_INT, send_indices, recv_indices, MPI_REPLACE));
479: PetscCall(PetscSFDestroy(&sf));
480: PetscCall(PetscFree(fromranks));
481: PetscCall(PetscFree(fromsizes));
482: PetscCall(PetscFree(fromperm_newtoold));
483: PetscCall(PetscFree(send_indices));
484: if (rows) {
485: PetscCall(PetscSortInt(nrecvs, recv_indices));
486: PetscCall(ISCreateGeneral(comm, nrecvs, recv_indices, PETSC_OWN_POINTER, rows));
487: } else {
488: PetscCall(PetscFree(recv_indices));
489: }
490: PetscFunctionReturn(PETSC_SUCCESS);
491: }
493: /*@
494: ISPartitioningToNumbering - Takes an `IS' that represents a partitioning (the MPI rank that each local entry belongs to) and on each MPI process
495: generates an `IS` that contains a new global node number in the new ordering for each entry
497: Collective
499: Input Parameter:
500: . part - a partitioning as generated by `MatPartitioningApply()` or `MatPartitioningApplyND()`
502: Output Parameter:
503: . is - on each processor the index set that defines the global numbers
504: (in the new numbering) for all the nodes currently (before the partitioning)
505: on that processor
507: Level: advanced
509: Note:
510: The resulting `IS` tells where each local entry is mapped to in a new global ordering
512: .seealso: [](sec_scatter), `IS`, `MatPartitioningCreate()`, `AOCreateBasic()`, `ISPartitioningCount()`
513: @*/
514: PetscErrorCode ISPartitioningToNumbering(IS part, IS *is)
515: {
516: MPI_Comm comm;
517: IS ndorder;
518: PetscInt i, np, npt, n, *starts = NULL, *sums = NULL, *lsizes = NULL, *newi = NULL;
519: const PetscInt *indices = NULL;
521: PetscFunctionBegin;
523: PetscAssertPointer(is, 2);
524: /* see if the partitioning comes from nested dissection */
525: PetscCall(PetscObjectQuery((PetscObject)part, "_petsc_matpartitioning_ndorder", (PetscObject *)&ndorder));
526: if (ndorder) {
527: PetscCall(PetscObjectReference((PetscObject)ndorder));
528: *is = ndorder;
529: PetscFunctionReturn(PETSC_SUCCESS);
530: }
532: PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
533: /* count the number of partitions, i.e., virtual processors */
534: PetscCall(ISGetLocalSize(part, &n));
535: PetscCall(ISGetIndices(part, &indices));
536: np = 0;
537: for (i = 0; i < n; i++) np = PetscMax(np, indices[i]);
538: PetscCall(MPIU_Allreduce(&np, &npt, 1, MPIU_INT, MPI_MAX, comm));
539: np = npt + 1; /* so that it looks like a MPI_Comm_size output */
541: /*
542: lsizes - number of elements of each partition on this particular processor
543: sums - total number of "previous" nodes for any particular partition
544: starts - global number of first element in each partition on this processor
545: */
546: PetscCall(PetscMalloc3(np, &lsizes, np, &starts, np, &sums));
547: PetscCall(PetscArrayzero(lsizes, np));
548: for (i = 0; i < n; i++) lsizes[indices[i]]++;
549: PetscCall(MPIU_Allreduce(lsizes, sums, np, MPIU_INT, MPI_SUM, comm));
550: PetscCallMPI(MPI_Scan(lsizes, starts, np, MPIU_INT, MPI_SUM, comm));
551: for (i = 0; i < np; i++) starts[i] -= lsizes[i];
552: for (i = 1; i < np; i++) {
553: sums[i] += sums[i - 1];
554: starts[i] += sums[i - 1];
555: }
557: /*
558: For each local index give it the new global number
559: */
560: PetscCall(PetscMalloc1(n, &newi));
561: for (i = 0; i < n; i++) newi[i] = starts[indices[i]]++;
562: PetscCall(PetscFree3(lsizes, starts, sums));
564: PetscCall(ISRestoreIndices(part, &indices));
565: PetscCall(ISCreateGeneral(comm, n, newi, PETSC_OWN_POINTER, is));
566: PetscCall(ISSetPermutation(*is));
567: PetscFunctionReturn(PETSC_SUCCESS);
568: }
570: /*@
571: ISPartitioningCount - Takes a `IS` that represents a partitioning (the MPI rank that each local entry belongs to) and determines the number of
572: resulting elements on each (partition) rank
574: Collective
576: Input Parameters:
577: + part - a partitioning as generated by `MatPartitioningApply()` or `MatPartitioningApplyND()`
578: - len - length of the array count, this is the total number of partitions
580: Output Parameter:
581: . count - array of length size, to contain the number of elements assigned
582: to each partition, where size is the number of partitions generated
583: (see notes below).
585: Level: advanced
587: Notes:
588: By default the number of partitions generated (and thus the length
589: of count) is the size of the communicator associated with `IS`,
590: but it can be set by `MatPartitioningSetNParts()`.
592: The resulting array of lengths can for instance serve as input of `PCBJacobiSetTotalBlocks()`.
594: If the partitioning has been obtained by `MatPartitioningApplyND()`, the returned count does not include the separators.
596: .seealso: [](sec_scatter), `IS`, `MatPartitioningCreate()`, `AOCreateBasic()`, `ISPartitioningToNumbering()`,
597: `MatPartitioningSetNParts()`, `MatPartitioningApply()`, `MatPartitioningApplyND()`
598: @*/
599: PetscErrorCode ISPartitioningCount(IS part, PetscInt len, PetscInt count[])
600: {
601: MPI_Comm comm;
602: PetscInt i, n, *lsizes;
603: const PetscInt *indices;
604: PetscMPIInt npp;
606: PetscFunctionBegin;
607: PetscCall(PetscObjectGetComm((PetscObject)part, &comm));
608: if (len == PETSC_DEFAULT) {
609: PetscMPIInt size;
610: PetscCallMPI(MPI_Comm_size(comm, &size));
611: len = (PetscInt)size;
612: }
614: /* count the number of partitions */
615: PetscCall(ISGetLocalSize(part, &n));
616: PetscCall(ISGetIndices(part, &indices));
617: if (PetscDefined(USE_DEBUG)) {
618: PetscInt np = 0, npt;
619: for (i = 0; i < n; i++) np = PetscMax(np, indices[i]);
620: PetscCall(MPIU_Allreduce(&np, &npt, 1, MPIU_INT, MPI_MAX, comm));
621: np = npt + 1; /* so that it looks like a MPI_Comm_size output */
622: PetscCheck(np <= len, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Length of count array %" PetscInt_FMT " is less than number of partitions %" PetscInt_FMT, len, np);
623: }
625: /*
626: lsizes - number of elements of each partition on this particular processor
627: sums - total number of "previous" nodes for any particular partition
628: starts - global number of first element in each partition on this processor
629: */
630: PetscCall(PetscCalloc1(len, &lsizes));
631: for (i = 0; i < n; i++) {
632: if (indices[i] > -1) lsizes[indices[i]]++;
633: }
634: PetscCall(ISRestoreIndices(part, &indices));
635: PetscCall(PetscMPIIntCast(len, &npp));
636: PetscCall(MPIU_Allreduce(lsizes, count, npp, MPIU_INT, MPI_SUM, comm));
637: PetscCall(PetscFree(lsizes));
638: PetscFunctionReturn(PETSC_SUCCESS);
639: }
641: /*@
642: ISAllGather - Given an index set `IS` on each processor, generates a large
643: index set (same on each processor) by concatenating together each
644: processors index set.
646: Collective
648: Input Parameter:
649: . is - the distributed index set
651: Output Parameter:
652: . isout - the concatenated index set (same on all processors)
654: Level: intermediate
656: Notes:
657: `ISAllGather()` is clearly not scalable for large index sets.
659: The `IS` created on each processor must be created with a common
660: communicator (e.g., `PETSC_COMM_WORLD`). If the index sets were created
661: with `PETSC_COMM_SELF`, this routine will not work as expected, since
662: each process will generate its own new `IS` that consists only of
663: itself.
665: The communicator for this new `IS` is `PETSC_COMM_SELF`
667: .seealso: [](sec_scatter), `IS`, `ISCreateGeneral()`, `ISCreateStride()`, `ISCreateBlock()`
668: @*/
669: PetscErrorCode ISAllGather(IS is, IS *isout)
670: {
671: PetscInt *indices, n, i, N, step, first;
672: const PetscInt *lindices;
673: MPI_Comm comm;
674: PetscMPIInt size, *sizes = NULL, *offsets = NULL, nn;
675: PetscBool stride;
677: PetscFunctionBegin;
679: PetscAssertPointer(isout, 2);
681: PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
682: PetscCallMPI(MPI_Comm_size(comm, &size));
683: PetscCall(ISGetLocalSize(is, &n));
684: PetscCall(PetscObjectTypeCompare((PetscObject)is, ISSTRIDE, &stride));
685: if (size == 1 && stride) { /* should handle parallel ISStride also */
686: PetscCall(ISStrideGetInfo(is, &first, &step));
687: PetscCall(ISCreateStride(PETSC_COMM_SELF, n, first, step, isout));
688: } else {
689: PetscCall(PetscMalloc2(size, &sizes, size, &offsets));
691: PetscCall(PetscMPIIntCast(n, &nn));
692: PetscCallMPI(MPI_Allgather(&nn, 1, MPI_INT, sizes, 1, MPI_INT, comm));
693: offsets[0] = 0;
694: for (i = 1; i < size; i++) {
695: PetscInt s = offsets[i - 1] + sizes[i - 1];
696: PetscCall(PetscMPIIntCast(s, &offsets[i]));
697: }
698: N = offsets[size - 1] + sizes[size - 1];
700: PetscCall(PetscMalloc1(N, &indices));
701: PetscCall(ISGetIndices(is, &lindices));
702: PetscCallMPI(MPI_Allgatherv((void *)lindices, nn, MPIU_INT, indices, sizes, offsets, MPIU_INT, comm));
703: PetscCall(ISRestoreIndices(is, &lindices));
704: PetscCall(PetscFree2(sizes, offsets));
706: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, N, indices, PETSC_OWN_POINTER, isout));
707: }
708: PetscFunctionReturn(PETSC_SUCCESS);
709: }
711: /*@C
712: ISAllGatherColors - Given a set of colors on each processor, generates a large
713: set (same on each processor) by concatenating together each processors colors
715: Collective
717: Input Parameters:
718: + comm - communicator to share the indices
719: . n - local size of set
720: - lindices - local colors
722: Output Parameters:
723: + outN - total number of indices
724: - outindices - all of the colors
726: Level: intermediate
728: Note:
729: `ISAllGatherColors()` is clearly not scalable for large index sets.
731: .seealso: `ISCOloringValue`, `ISColoring()`, `ISCreateGeneral()`, `ISCreateStride()`, `ISCreateBlock()`, `ISAllGather()`
732: @*/
733: PetscErrorCode ISAllGatherColors(MPI_Comm comm, PetscInt n, ISColoringValue lindices[], PetscInt *outN, ISColoringValue *outindices[])
734: {
735: ISColoringValue *indices;
736: PetscInt i, N;
737: PetscMPIInt size, *offsets = NULL, *sizes = NULL, nn = n;
739: PetscFunctionBegin;
740: PetscCallMPI(MPI_Comm_size(comm, &size));
741: PetscCall(PetscMalloc2(size, &sizes, size, &offsets));
743: PetscCallMPI(MPI_Allgather(&nn, 1, MPI_INT, sizes, 1, MPI_INT, comm));
744: offsets[0] = 0;
745: for (i = 1; i < size; i++) offsets[i] = offsets[i - 1] + sizes[i - 1];
746: N = offsets[size - 1] + sizes[size - 1];
747: PetscCall(PetscFree2(sizes, offsets));
749: PetscCall(PetscMalloc1(N + 1, &indices));
750: PetscCallMPI(MPI_Allgatherv(lindices, (PetscMPIInt)n, MPIU_COLORING_VALUE, indices, sizes, offsets, MPIU_COLORING_VALUE, comm));
752: *outindices = indices;
753: if (outN) *outN = N;
754: PetscFunctionReturn(PETSC_SUCCESS);
755: }
757: /*@
758: ISComplement - Given an index set `IS` generates the complement index set. That is
759: all indices that are NOT in the given set.
761: Collective
763: Input Parameters:
764: + is - the index set
765: . nmin - the first index desired in the local part of the complement
766: - nmax - the largest index desired in the local part of the complement (note that all indices in `is` must be greater or equal to `nmin` and less than `nmax`)
768: Output Parameter:
769: . isout - the complement
771: Level: intermediate
773: Notes:
774: The communicator for `isout` is the same as for the input `is`
776: For a parallel `is`, this will generate the local part of the complement on each process
778: To generate the entire complement (on each process) of a parallel `is`, first call `ISAllGather()` and then
779: call this routine.
781: .seealso: [](sec_scatter), `IS`, `ISCreateGeneral()`, `ISCreateStride()`, `ISCreateBlock()`, `ISAllGather()`
782: @*/
783: PetscErrorCode ISComplement(IS is, PetscInt nmin, PetscInt nmax, IS *isout)
784: {
785: const PetscInt *indices;
786: PetscInt n, i, j, unique, cnt, *nindices;
787: PetscBool sorted;
789: PetscFunctionBegin;
791: PetscAssertPointer(isout, 4);
792: PetscCheck(nmin >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nmin %" PetscInt_FMT " cannot be negative", nmin);
793: PetscCheck(nmin <= nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "nmin %" PetscInt_FMT " cannot be greater than nmax %" PetscInt_FMT, nmin, nmax);
794: PetscCall(ISSorted(is, &sorted));
795: PetscCheck(sorted, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index set must be sorted");
797: PetscCall(ISGetLocalSize(is, &n));
798: PetscCall(ISGetIndices(is, &indices));
799: if (PetscDefined(USE_DEBUG)) {
800: for (i = 0; i < n; i++) {
801: PetscCheck(indices[i] >= nmin, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT "'s value %" PetscInt_FMT " is smaller than minimum given %" PetscInt_FMT, i, indices[i], nmin);
802: PetscCheck(indices[i] < nmax, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Index %" PetscInt_FMT "'s value %" PetscInt_FMT " is larger than maximum given %" PetscInt_FMT, i, indices[i], nmax);
803: }
804: }
805: /* Count number of unique entries */
806: unique = (n > 0);
807: for (i = 0; i < n - 1; i++) {
808: if (indices[i + 1] != indices[i]) unique++;
809: }
810: PetscCall(PetscMalloc1(nmax - nmin - unique, &nindices));
811: cnt = 0;
812: for (i = nmin, j = 0; i < nmax; i++) {
813: if (j < n && i == indices[j]) do {
814: j++;
815: } while (j < n && i == indices[j]);
816: else nindices[cnt++] = i;
817: }
818: PetscCheck(cnt == nmax - nmin - unique, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of entries found in complement %" PetscInt_FMT " does not match expected %" PetscInt_FMT, cnt, nmax - nmin - unique);
819: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)is), cnt, nindices, PETSC_OWN_POINTER, isout));
820: PetscCall(ISRestoreIndices(is, &indices));
821: PetscFunctionReturn(PETSC_SUCCESS);
822: }