Actual source code: index.c
1: /*
2: Defines the abstract operations on index sets, i.e. the public interface.
3: */
4: #include <petsc/private/isimpl.h>
5: #include <petscviewer.h>
6: #include <petscsf.h>
8: /* Logging support */
9: PetscClassId IS_CLASSID;
10: /* TODO: Much more events are missing! */
11: PetscLogEvent IS_View;
12: PetscLogEvent IS_Load;
14: /*@
15: ISRenumber - Renumbers the non-negative entries of an index set in a contiguous way, starting from 0.
17: Collective
19: Input Parameters:
20: + subset - the index set
21: - subset_mult - the multiplicity of each entry in subset (optional, can be NULL)
23: Output Parameters:
24: + N - one past the largest entry of the new IS
25: - subset_n - the new IS
27: Notes: All negative entries are mapped to -1. Indices with non positive multiplicities are skipped.
29: Level: intermediate
31: .seealso:
32: @*/
33: PetscErrorCode ISRenumber(IS subset, IS subset_mult, PetscInt *N, IS *subset_n)
34: {
35: PetscSF sf;
36: PetscLayout map;
37: const PetscInt *idxs, *idxs_mult = NULL;
38: PetscInt *leaf_data, *root_data, *gidxs, *ilocal, *ilocalneg;
39: PetscInt N_n, n, i, lbounds[2], gbounds[2], Nl, ibs;
40: PetscInt n_n, nlocals, start, first_index, npos, nneg;
41: PetscMPIInt commsize;
42: PetscBool first_found, isblock;
44: PetscFunctionBegin;
48: else if (!subset_n) PetscFunctionReturn(PETSC_SUCCESS);
49: PetscCall(ISGetLocalSize(subset, &n));
50: if (subset_mult) {
51: PetscCall(ISGetLocalSize(subset_mult, &i));
52: PetscCheck(i == n, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Local subset and multiplicity sizes don't match! %" PetscInt_FMT " != %" PetscInt_FMT, n, i);
53: }
54: /* create workspace layout for computing global indices of subset */
55: PetscCall(PetscMalloc1(n, &ilocal));
56: PetscCall(PetscMalloc1(n, &ilocalneg));
57: PetscCall(ISGetIndices(subset, &idxs));
58: PetscCall(ISGetBlockSize(subset, &ibs));
59: PetscCall(PetscObjectTypeCompare((PetscObject)subset, ISBLOCK, &isblock));
60: if (subset_mult) PetscCall(ISGetIndices(subset_mult, &idxs_mult));
61: lbounds[0] = PETSC_MAX_INT;
62: lbounds[1] = PETSC_MIN_INT;
63: for (i = 0, npos = 0, nneg = 0; i < n; i++) {
64: if (idxs[i] < 0) {
65: ilocalneg[nneg++] = i;
66: continue;
67: }
68: if (idxs[i] < lbounds[0]) lbounds[0] = idxs[i];
69: if (idxs[i] > lbounds[1]) lbounds[1] = idxs[i];
70: ilocal[npos++] = i;
71: }
72: if (npos == n) {
73: PetscCall(PetscFree(ilocal));
74: PetscCall(PetscFree(ilocalneg));
75: }
77: /* create sf : leaf_data == multiplicity of indexes, root data == global index in layout */
78: PetscCall(PetscMalloc1(n, &leaf_data));
79: for (i = 0; i < n; i++) leaf_data[i] = idxs_mult ? PetscMax(idxs_mult[i], 0) : 1;
81: /* local size of new subset */
82: n_n = 0;
83: for (i = 0; i < n; i++) n_n += leaf_data[i];
84: if (ilocalneg)
85: for (i = 0; i < nneg; i++) leaf_data[ilocalneg[i]] = 0;
86: PetscCall(PetscFree(ilocalneg));
87: PetscCall(PetscMalloc1(PetscMax(n_n, n), &gidxs)); /* allocating extra space to reuse gidxs */
88: /* check for early termination (all negative) */
89: PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)subset), lbounds, gbounds));
90: if (gbounds[1] < gbounds[0]) {
91: if (N) *N = 0;
92: if (subset_n) { /* all negative */
93: for (i = 0; i < n_n; i++) gidxs[i] = -1;
94: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)subset), n_n, gidxs, PETSC_COPY_VALUES, subset_n));
95: }
96: PetscCall(PetscFree(leaf_data));
97: PetscCall(PetscFree(gidxs));
98: PetscCall(ISRestoreIndices(subset, &idxs));
99: if (subset_mult) PetscCall(ISRestoreIndices(subset_mult, &idxs_mult));
100: PetscCall(PetscFree(ilocal));
101: PetscCall(PetscFree(ilocalneg));
102: PetscFunctionReturn(PETSC_SUCCESS);
103: }
105: /* split work */
106: N_n = gbounds[1] - gbounds[0] + 1;
107: PetscCall(PetscLayoutCreate(PetscObjectComm((PetscObject)subset), &map));
108: PetscCall(PetscLayoutSetBlockSize(map, 1));
109: PetscCall(PetscLayoutSetSize(map, N_n));
110: PetscCall(PetscLayoutSetUp(map));
111: PetscCall(PetscLayoutGetLocalSize(map, &Nl));
113: /* global indexes in layout */
114: for (i = 0; i < npos; i++) gidxs[i] = (ilocal ? idxs[ilocal[i]] : idxs[i]) - gbounds[0];
115: PetscCall(ISRestoreIndices(subset, &idxs));
116: PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)subset), &sf));
117: PetscCall(PetscSFSetGraphLayout(sf, map, npos, ilocal, PETSC_USE_POINTER, gidxs));
118: PetscCall(PetscLayoutDestroy(&map));
120: /* reduce from leaves to roots */
121: PetscCall(PetscCalloc1(Nl, &root_data));
122: PetscCall(PetscSFReduceBegin(sf, MPIU_INT, leaf_data, root_data, MPI_MAX));
123: PetscCall(PetscSFReduceEnd(sf, MPIU_INT, leaf_data, root_data, MPI_MAX));
125: /* count indexes in local part of layout */
126: nlocals = 0;
127: first_index = -1;
128: first_found = PETSC_FALSE;
129: for (i = 0; i < Nl; i++) {
130: if (!first_found && root_data[i]) {
131: first_found = PETSC_TRUE;
132: first_index = i;
133: }
134: nlocals += root_data[i];
135: }
137: /* cumulative of number of indexes and size of subset without holes */
138: #if defined(PETSC_HAVE_MPI_EXSCAN)
139: start = 0;
140: PetscCallMPI(MPI_Exscan(&nlocals, &start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)subset)));
141: #else
142: PetscCallMPI(MPI_Scan(&nlocals, &start, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)subset)));
143: start = start - nlocals;
144: #endif
146: if (N) { /* compute total size of new subset if requested */
147: *N = start + nlocals;
148: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)subset), &commsize));
149: PetscCallMPI(MPI_Bcast(N, 1, MPIU_INT, commsize - 1, PetscObjectComm((PetscObject)subset)));
150: }
152: if (!subset_n) {
153: PetscCall(PetscFree(gidxs));
154: PetscCall(PetscSFDestroy(&sf));
155: PetscCall(PetscFree(leaf_data));
156: PetscCall(PetscFree(root_data));
157: PetscCall(PetscFree(ilocal));
158: if (subset_mult) PetscCall(ISRestoreIndices(subset_mult, &idxs_mult));
159: PetscFunctionReturn(PETSC_SUCCESS);
160: }
162: /* adapt root data with cumulative */
163: if (first_found) {
164: PetscInt old_index;
166: root_data[first_index] += start;
167: old_index = first_index;
168: for (i = first_index + 1; i < Nl; i++) {
169: if (root_data[i]) {
170: root_data[i] += root_data[old_index];
171: old_index = i;
172: }
173: }
174: }
176: /* from roots to leaves */
177: PetscCall(PetscSFBcastBegin(sf, MPIU_INT, root_data, leaf_data, MPI_REPLACE));
178: PetscCall(PetscSFBcastEnd(sf, MPIU_INT, root_data, leaf_data, MPI_REPLACE));
179: PetscCall(PetscSFDestroy(&sf));
181: /* create new IS with global indexes without holes */
182: for (i = 0; i < n_n; i++) gidxs[i] = -1;
183: if (subset_mult) {
184: PetscInt cum;
186: isblock = PETSC_FALSE;
187: for (i = 0, cum = 0; i < n; i++)
188: for (PetscInt j = 0; j < idxs_mult[i]; j++) gidxs[cum++] = leaf_data[i] - idxs_mult[i] + j;
189: } else
190: for (i = 0; i < n; i++) gidxs[i] = leaf_data[i] - 1;
192: if (isblock) {
193: if (ibs > 1)
194: for (i = 0; i < n_n / ibs; i++) gidxs[i] = gidxs[i * ibs] / ibs;
195: PetscCall(ISCreateBlock(PetscObjectComm((PetscObject)subset), ibs, n_n / ibs, gidxs, PETSC_COPY_VALUES, subset_n));
196: } else {
197: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)subset), n_n, gidxs, PETSC_COPY_VALUES, subset_n));
198: }
199: if (subset_mult) PetscCall(ISRestoreIndices(subset_mult, &idxs_mult));
200: PetscCall(PetscFree(gidxs));
201: PetscCall(PetscFree(leaf_data));
202: PetscCall(PetscFree(root_data));
203: PetscCall(PetscFree(ilocal));
204: PetscFunctionReturn(PETSC_SUCCESS);
205: }
207: /*@
208: ISCreateSubIS - Create a sub index set from a global index set selecting some components.
210: Collective
212: Input Parameters:
213: + is - the index set
214: - comps - which components we will extract from is
216: Output Parameters:
217: . subis - the new sub index set
219: Level: intermediate
221: Example usage:
222: We have an index set (is) living on 3 processes with the following values:
223: | 4 9 0 | 2 6 7 | 10 11 1|
224: and another index set (comps) used to indicate which components of is we want to take,
225: | 7 5 | 1 2 | 0 4|
226: The output index set (subis) should look like:
227: | 11 7 | 9 0 | 4 6|
229: .seealso: `VecGetSubVector()`, `MatCreateSubMatrix()`
230: @*/
231: PetscErrorCode ISCreateSubIS(IS is, IS comps, IS *subis)
232: {
233: PetscSF sf;
234: const PetscInt *is_indices, *comps_indices;
235: PetscInt *subis_indices, nroots, nleaves, *mine, i, lidx;
236: PetscMPIInt owner;
237: PetscSFNode *remote;
238: MPI_Comm comm;
240: PetscFunctionBegin;
245: PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
246: PetscCall(ISGetLocalSize(comps, &nleaves));
247: PetscCall(ISGetLocalSize(is, &nroots));
248: PetscCall(PetscMalloc1(nleaves, &remote));
249: PetscCall(PetscMalloc1(nleaves, &mine));
250: PetscCall(ISGetIndices(comps, &comps_indices));
251: /*
252: * Construct a PetscSF in which "is" data serves as roots and "subis" is leaves.
253: * Root data are sent to leaves using PetscSFBcast().
254: * */
255: for (i = 0; i < nleaves; i++) {
256: mine[i] = i;
257: /* Connect a remote root with the current leaf. The value on the remote root
258: * will be received by the current local leaf.
259: * */
260: owner = -1;
261: lidx = -1;
262: PetscCall(PetscLayoutFindOwnerIndex(is->map, comps_indices[i], &owner, &lidx));
263: remote[i].rank = owner;
264: remote[i].index = lidx;
265: }
266: PetscCall(ISRestoreIndices(comps, &comps_indices));
267: PetscCall(PetscSFCreate(comm, &sf));
268: PetscCall(PetscSFSetFromOptions(sf));
269: PetscCall(PetscSFSetGraph(sf, nroots, nleaves, mine, PETSC_OWN_POINTER, remote, PETSC_OWN_POINTER));
271: PetscCall(PetscMalloc1(nleaves, &subis_indices));
272: PetscCall(ISGetIndices(is, &is_indices));
273: PetscCall(PetscSFBcastBegin(sf, MPIU_INT, is_indices, subis_indices, MPI_REPLACE));
274: PetscCall(PetscSFBcastEnd(sf, MPIU_INT, is_indices, subis_indices, MPI_REPLACE));
275: PetscCall(ISRestoreIndices(is, &is_indices));
276: PetscCall(PetscSFDestroy(&sf));
277: PetscCall(ISCreateGeneral(comm, nleaves, subis_indices, PETSC_OWN_POINTER, subis));
278: PetscFunctionReturn(PETSC_SUCCESS);
279: }
281: /*@
282: ISClearInfoCache - clear the cache of computed index set properties
284: Not collective
286: Input Parameters:
287: + is - the index set
288: - clear_permanent_local - whether to remove the permanent status of local properties
290: NOTE: because all processes must agree on the global permanent status of a property,
291: the permanent status can only be changed with ISSetInfo(), because this routine is not collective
293: Level: developer
295: .seealso: `ISInfo`, `ISInfoType`, `ISSetInfo()`, `ISClearInfoCache()`
297: @*/
298: PetscErrorCode ISClearInfoCache(IS is, PetscBool clear_permanent_local)
299: {
300: PetscInt i, j;
302: PetscFunctionBegin;
305: for (i = 0; i < IS_INFO_MAX; i++) {
306: if (clear_permanent_local) is->info_permanent[IS_LOCAL][i] = PETSC_FALSE;
307: for (j = 0; j < 2; j++) {
308: if (!is->info_permanent[j][i]) is->info[j][i] = IS_INFO_UNKNOWN;
309: }
310: }
311: PetscFunctionReturn(PETSC_SUCCESS);
312: }
314: static PetscErrorCode ISSetInfo_Internal(IS is, ISInfo info, ISInfoType type, ISInfoBool ipermanent, PetscBool flg)
315: {
316: ISInfoBool iflg = flg ? IS_INFO_TRUE : IS_INFO_FALSE;
317: PetscInt itype = (type == IS_LOCAL) ? 0 : 1;
318: PetscBool permanent_set = (ipermanent == IS_INFO_UNKNOWN) ? PETSC_FALSE : PETSC_TRUE;
319: PetscBool permanent = (ipermanent == IS_INFO_TRUE) ? PETSC_TRUE : PETSC_FALSE;
321: PetscFunctionBegin;
322: /* set this property */
323: is->info[itype][(int)info] = iflg;
324: if (permanent_set) is->info_permanent[itype][(int)info] = permanent;
325: /* set implications */
326: switch (info) {
327: case IS_SORTED:
328: if (flg && type == IS_GLOBAL) { /* an array that is globally sorted is also locally sorted */
329: is->info[IS_LOCAL][(int)info] = IS_INFO_TRUE;
330: /* global permanence implies local permanence */
331: if (permanent_set && permanent) is->info_permanent[IS_LOCAL][(int)info] = PETSC_TRUE;
332: }
333: if (!flg) { /* if an array is not sorted, it cannot be an interval or the identity */
334: is->info[itype][IS_INTERVAL] = IS_INFO_FALSE;
335: is->info[itype][IS_IDENTITY] = IS_INFO_FALSE;
336: if (permanent_set) {
337: is->info_permanent[itype][IS_INTERVAL] = permanent;
338: is->info_permanent[itype][IS_IDENTITY] = permanent;
339: }
340: }
341: break;
342: case IS_UNIQUE:
343: if (flg && type == IS_GLOBAL) { /* an array that is globally unique is also locally unique */
344: is->info[IS_LOCAL][(int)info] = IS_INFO_TRUE;
345: /* global permanence implies local permanence */
346: if (permanent_set && permanent) is->info_permanent[IS_LOCAL][(int)info] = PETSC_TRUE;
347: }
348: if (!flg) { /* if an array is not unique, it cannot be a permutation, and interval, or the identity */
349: is->info[itype][IS_PERMUTATION] = IS_INFO_FALSE;
350: is->info[itype][IS_INTERVAL] = IS_INFO_FALSE;
351: is->info[itype][IS_IDENTITY] = IS_INFO_FALSE;
352: if (permanent_set) {
353: is->info_permanent[itype][IS_PERMUTATION] = permanent;
354: is->info_permanent[itype][IS_INTERVAL] = permanent;
355: is->info_permanent[itype][IS_IDENTITY] = permanent;
356: }
357: }
358: break;
359: case IS_PERMUTATION:
360: if (flg) { /* an array that is a permutation is unique and is unique locally */
361: is->info[itype][IS_UNIQUE] = IS_INFO_TRUE;
362: is->info[IS_LOCAL][IS_UNIQUE] = IS_INFO_TRUE;
363: if (permanent_set && permanent) {
364: is->info_permanent[itype][IS_UNIQUE] = PETSC_TRUE;
365: is->info_permanent[IS_LOCAL][IS_UNIQUE] = PETSC_TRUE;
366: }
367: } else { /* an array that is not a permutation cannot be the identity */
368: is->info[itype][IS_IDENTITY] = IS_INFO_FALSE;
369: if (permanent_set) is->info_permanent[itype][IS_IDENTITY] = permanent;
370: }
371: break;
372: case IS_INTERVAL:
373: if (flg) { /* an array that is an interval is sorted and unique */
374: is->info[itype][IS_SORTED] = IS_INFO_TRUE;
375: is->info[IS_LOCAL][IS_SORTED] = IS_INFO_TRUE;
376: is->info[itype][IS_UNIQUE] = IS_INFO_TRUE;
377: is->info[IS_LOCAL][IS_UNIQUE] = IS_INFO_TRUE;
378: if (permanent_set && permanent) {
379: is->info_permanent[itype][IS_SORTED] = PETSC_TRUE;
380: is->info_permanent[IS_LOCAL][IS_SORTED] = PETSC_TRUE;
381: is->info_permanent[itype][IS_UNIQUE] = PETSC_TRUE;
382: is->info_permanent[IS_LOCAL][IS_UNIQUE] = PETSC_TRUE;
383: }
384: } else { /* an array that is not an interval cannot be the identity */
385: is->info[itype][IS_IDENTITY] = IS_INFO_FALSE;
386: if (permanent_set) is->info_permanent[itype][IS_IDENTITY] = permanent;
387: }
388: break;
389: case IS_IDENTITY:
390: if (flg) { /* an array that is the identity is sorted, unique, an interval, and a permutation */
391: is->info[itype][IS_SORTED] = IS_INFO_TRUE;
392: is->info[IS_LOCAL][IS_SORTED] = IS_INFO_TRUE;
393: is->info[itype][IS_UNIQUE] = IS_INFO_TRUE;
394: is->info[IS_LOCAL][IS_UNIQUE] = IS_INFO_TRUE;
395: is->info[itype][IS_PERMUTATION] = IS_INFO_TRUE;
396: is->info[itype][IS_INTERVAL] = IS_INFO_TRUE;
397: is->info[IS_LOCAL][IS_INTERVAL] = IS_INFO_TRUE;
398: if (permanent_set && permanent) {
399: is->info_permanent[itype][IS_SORTED] = PETSC_TRUE;
400: is->info_permanent[IS_LOCAL][IS_SORTED] = PETSC_TRUE;
401: is->info_permanent[itype][IS_UNIQUE] = PETSC_TRUE;
402: is->info_permanent[IS_LOCAL][IS_UNIQUE] = PETSC_TRUE;
403: is->info_permanent[itype][IS_PERMUTATION] = PETSC_TRUE;
404: is->info_permanent[itype][IS_INTERVAL] = PETSC_TRUE;
405: is->info_permanent[IS_LOCAL][IS_INTERVAL] = PETSC_TRUE;
406: }
407: }
408: break;
409: default:
410: PetscCheck(type != IS_LOCAL, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unknown IS property");
411: SETERRQ(PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_OUTOFRANGE, "Unknown IS property");
412: }
413: PetscFunctionReturn(PETSC_SUCCESS);
414: }
416: /*@
417: ISSetInfo - Set known information about an index set.
419: Logically Collective on IS if type is IS_GLOBAL
421: Input Parameters:
422: + is - the index set
423: . info - describing a property of the index set, one of those listed below,
424: . type - IS_LOCAL if the information describes the local portion of the index set,
425: IS_GLOBAL if it describes the whole index set
426: . permanent - PETSC_TRUE if it is known that the property will persist through changes to the index set, PETSC_FALSE otherwise
427: If the user sets a property as permanently known, it will bypass computation of that property
428: - flg - set the described property as true (PETSC_TRUE) or false (PETSC_FALSE)
430: Info Describing IS Structure:
431: + IS_SORTED - the [local part of the] index set is sorted in ascending order
432: . IS_UNIQUE - each entry in the [local part of the] index set is unique
433: . IS_PERMUTATION - the [local part of the] index set is a permutation of the integers {0, 1, ..., N-1}, where N is the size of the [local part of the] index set
434: . IS_INTERVAL - the [local part of the] index set is equal to a contiguous range of integers {f, f + 1, ..., f + N-1}
435: - IS_IDENTITY - the [local part of the] index set is equal to the integers {0, 1, ..., N-1}
437: Notes:
438: If type is IS_GLOBAL, all processes that share the index set must pass the same value in flg
440: It is possible to set a property with ISSetInfo() that contradicts what would be previously computed with ISGetInfo()
442: Level: advanced
444: .seealso: `ISInfo`, `ISInfoType`, `IS`
446: @*/
447: PetscErrorCode ISSetInfo(IS is, ISInfo info, ISInfoType type, PetscBool permanent, PetscBool flg)
448: {
449: MPI_Comm comm, errcomm;
450: PetscMPIInt size;
452: PetscFunctionBegin;
455: comm = PetscObjectComm((PetscObject)is);
456: if (type == IS_GLOBAL) {
460: errcomm = comm;
461: } else {
462: errcomm = PETSC_COMM_SELF;
463: }
465: PetscCheck(((int)info) > IS_INFO_MIN && ((int)info) < IS_INFO_MAX, errcomm, PETSC_ERR_ARG_OUTOFRANGE, "Options %d is out of range", (int)info);
467: PetscCallMPI(MPI_Comm_size(comm, &size));
468: /* do not use global values if size == 1: it makes it easier to keep the implications straight */
469: if (size == 1) type = IS_LOCAL;
470: PetscCall(ISSetInfo_Internal(is, info, type, permanent ? IS_INFO_TRUE : IS_INFO_FALSE, flg));
471: PetscFunctionReturn(PETSC_SUCCESS);
472: }
474: static PetscErrorCode ISGetInfo_Sorted(IS is, ISInfoType type, PetscBool *flg)
475: {
476: MPI_Comm comm;
477: PetscMPIInt size, rank;
479: PetscFunctionBegin;
480: comm = PetscObjectComm((PetscObject)is);
481: PetscCallMPI(MPI_Comm_size(comm, &size));
482: PetscCallMPI(MPI_Comm_size(comm, &rank));
483: if (type == IS_GLOBAL && is->ops->sortedglobal) {
484: PetscUseTypeMethod(is, sortedglobal, flg);
485: } else {
486: PetscBool sortedLocal = PETSC_FALSE;
488: /* determine if the array is locally sorted */
489: if (type == IS_GLOBAL && size > 1) {
490: /* call ISGetInfo so that a cached value will be used if possible */
491: PetscCall(ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, &sortedLocal));
492: } else if (is->ops->sortedlocal) {
493: PetscUseTypeMethod(is, sortedlocal, &sortedLocal);
494: } else {
495: /* default: get the local indices and directly check */
496: const PetscInt *idx;
497: PetscInt n;
499: PetscCall(ISGetIndices(is, &idx));
500: PetscCall(ISGetLocalSize(is, &n));
501: PetscCall(PetscSortedInt(n, idx, &sortedLocal));
502: PetscCall(ISRestoreIndices(is, &idx));
503: }
505: if (type == IS_LOCAL || size == 1) {
506: *flg = sortedLocal;
507: } else {
508: PetscCallMPI(MPI_Allreduce(&sortedLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
509: if (*flg) {
510: PetscInt n, min = PETSC_MAX_INT, max = PETSC_MIN_INT;
511: PetscInt maxprev;
513: PetscCall(ISGetLocalSize(is, &n));
514: if (n) PetscCall(ISGetMinMax(is, &min, &max));
515: maxprev = PETSC_MIN_INT;
516: PetscCallMPI(MPI_Exscan(&max, &maxprev, 1, MPIU_INT, MPI_MAX, comm));
517: if (rank && (maxprev > min)) sortedLocal = PETSC_FALSE;
518: PetscCallMPI(MPI_Allreduce(&sortedLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
519: }
520: }
521: }
522: PetscFunctionReturn(PETSC_SUCCESS);
523: }
525: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[]);
527: static PetscErrorCode ISGetInfo_Unique(IS is, ISInfoType type, PetscBool *flg)
528: {
529: MPI_Comm comm;
530: PetscMPIInt size, rank;
531: PetscInt i;
533: PetscFunctionBegin;
534: comm = PetscObjectComm((PetscObject)is);
535: PetscCallMPI(MPI_Comm_size(comm, &size));
536: PetscCallMPI(MPI_Comm_size(comm, &rank));
537: if (type == IS_GLOBAL && is->ops->uniqueglobal) {
538: PetscUseTypeMethod(is, uniqueglobal, flg);
539: } else {
540: PetscBool uniqueLocal;
541: PetscInt n = -1;
542: PetscInt *idx = NULL;
544: /* determine if the array is locally unique */
545: if (type == IS_GLOBAL && size > 1) {
546: /* call ISGetInfo so that a cached value will be used if possible */
547: PetscCall(ISGetInfo(is, IS_UNIQUE, IS_LOCAL, PETSC_TRUE, &uniqueLocal));
548: } else if (is->ops->uniquelocal) {
549: PetscUseTypeMethod(is, uniquelocal, &uniqueLocal);
550: } else {
551: /* default: get the local indices and directly check */
552: uniqueLocal = PETSC_TRUE;
553: PetscCall(ISGetLocalSize(is, &n));
554: PetscCall(PetscMalloc1(n, &idx));
555: PetscCall(ISGetIndicesCopy(is, idx));
556: PetscCall(PetscIntSortSemiOrdered(n, idx));
557: for (i = 1; i < n; i++)
558: if (idx[i] == idx[i - 1]) break;
559: if (i < n) uniqueLocal = PETSC_FALSE;
560: }
562: PetscCall(PetscFree(idx));
563: if (type == IS_LOCAL || size == 1) {
564: *flg = uniqueLocal;
565: } else {
566: PetscCallMPI(MPI_Allreduce(&uniqueLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
567: if (*flg) {
568: PetscInt min = PETSC_MAX_INT, max = PETSC_MIN_INT, maxprev;
570: if (!idx) {
571: PetscCall(ISGetLocalSize(is, &n));
572: PetscCall(PetscMalloc1(n, &idx));
573: PetscCall(ISGetIndicesCopy(is, idx));
574: }
575: PetscCall(PetscParallelSortInt(is->map, is->map, idx, idx));
576: if (n) {
577: min = idx[0];
578: max = idx[n - 1];
579: }
580: for (i = 1; i < n; i++)
581: if (idx[i] == idx[i - 1]) break;
582: if (i < n) uniqueLocal = PETSC_FALSE;
583: maxprev = PETSC_MIN_INT;
584: PetscCallMPI(MPI_Exscan(&max, &maxprev, 1, MPIU_INT, MPI_MAX, comm));
585: if (rank && (maxprev == min)) uniqueLocal = PETSC_FALSE;
586: PetscCallMPI(MPI_Allreduce(&uniqueLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
587: }
588: }
589: PetscCall(PetscFree(idx));
590: }
591: PetscFunctionReturn(PETSC_SUCCESS);
592: }
594: static PetscErrorCode ISGetInfo_Permutation(IS is, ISInfoType type, PetscBool *flg)
595: {
596: MPI_Comm comm;
597: PetscMPIInt size, rank;
599: PetscFunctionBegin;
600: comm = PetscObjectComm((PetscObject)is);
601: PetscCallMPI(MPI_Comm_size(comm, &size));
602: PetscCallMPI(MPI_Comm_size(comm, &rank));
603: if (type == IS_GLOBAL && is->ops->permglobal) {
604: PetscUseTypeMethod(is, permglobal, flg);
605: } else if (type == IS_LOCAL && is->ops->permlocal) {
606: PetscUseTypeMethod(is, permlocal, flg);
607: } else {
608: PetscBool permLocal;
609: PetscInt n, i, rStart;
610: PetscInt *idx;
612: PetscCall(ISGetLocalSize(is, &n));
613: PetscCall(PetscMalloc1(n, &idx));
614: PetscCall(ISGetIndicesCopy(is, idx));
615: if (type == IS_GLOBAL) {
616: PetscCall(PetscParallelSortInt(is->map, is->map, idx, idx));
617: PetscCall(PetscLayoutGetRange(is->map, &rStart, NULL));
618: } else {
619: PetscCall(PetscIntSortSemiOrdered(n, idx));
620: rStart = 0;
621: }
622: permLocal = PETSC_TRUE;
623: for (i = 0; i < n; i++) {
624: if (idx[i] != rStart + i) break;
625: }
626: if (i < n) permLocal = PETSC_FALSE;
627: if (type == IS_LOCAL || size == 1) {
628: *flg = permLocal;
629: } else {
630: PetscCallMPI(MPI_Allreduce(&permLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
631: }
632: PetscCall(PetscFree(idx));
633: }
634: PetscFunctionReturn(PETSC_SUCCESS);
635: }
637: static PetscErrorCode ISGetInfo_Interval(IS is, ISInfoType type, PetscBool *flg)
638: {
639: MPI_Comm comm;
640: PetscMPIInt size, rank;
641: PetscInt i;
643: PetscFunctionBegin;
644: comm = PetscObjectComm((PetscObject)is);
645: PetscCallMPI(MPI_Comm_size(comm, &size));
646: PetscCallMPI(MPI_Comm_size(comm, &rank));
647: if (type == IS_GLOBAL && is->ops->intervalglobal) {
648: PetscUseTypeMethod(is, intervalglobal, flg);
649: } else {
650: PetscBool intervalLocal;
652: /* determine if the array is locally an interval */
653: if (type == IS_GLOBAL && size > 1) {
654: /* call ISGetInfo so that a cached value will be used if possible */
655: PetscCall(ISGetInfo(is, IS_INTERVAL, IS_LOCAL, PETSC_TRUE, &intervalLocal));
656: } else if (is->ops->intervallocal) {
657: PetscUseTypeMethod(is, intervallocal, &intervalLocal);
658: } else {
659: PetscInt n;
660: const PetscInt *idx;
661: /* default: get the local indices and directly check */
662: intervalLocal = PETSC_TRUE;
663: PetscCall(ISGetLocalSize(is, &n));
664: PetscCall(ISGetIndices(is, &idx));
665: for (i = 1; i < n; i++)
666: if (idx[i] != idx[i - 1] + 1) break;
667: if (i < n) intervalLocal = PETSC_FALSE;
668: PetscCall(ISRestoreIndices(is, &idx));
669: }
671: if (type == IS_LOCAL || size == 1) {
672: *flg = intervalLocal;
673: } else {
674: PetscCallMPI(MPI_Allreduce(&intervalLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
675: if (*flg) {
676: PetscInt n, min = PETSC_MAX_INT, max = PETSC_MIN_INT;
677: PetscInt maxprev;
679: PetscCall(ISGetLocalSize(is, &n));
680: if (n) PetscCall(ISGetMinMax(is, &min, &max));
681: maxprev = PETSC_MIN_INT;
682: PetscCallMPI(MPI_Exscan(&max, &maxprev, 1, MPIU_INT, MPI_MAX, comm));
683: if (rank && n && (maxprev != min - 1)) intervalLocal = PETSC_FALSE;
684: PetscCallMPI(MPI_Allreduce(&intervalLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
685: }
686: }
687: }
688: PetscFunctionReturn(PETSC_SUCCESS);
689: }
691: static PetscErrorCode ISGetInfo_Identity(IS is, ISInfoType type, PetscBool *flg)
692: {
693: MPI_Comm comm;
694: PetscMPIInt size, rank;
696: PetscFunctionBegin;
697: comm = PetscObjectComm((PetscObject)is);
698: PetscCallMPI(MPI_Comm_size(comm, &size));
699: PetscCallMPI(MPI_Comm_size(comm, &rank));
700: if (type == IS_GLOBAL && is->ops->intervalglobal) {
701: PetscBool isinterval;
703: PetscUseTypeMethod(is, intervalglobal, &isinterval);
704: *flg = PETSC_FALSE;
705: if (isinterval) {
706: PetscInt min;
708: PetscCall(ISGetMinMax(is, &min, NULL));
709: PetscCallMPI(MPI_Bcast(&min, 1, MPIU_INT, 0, comm));
710: if (min == 0) *flg = PETSC_TRUE;
711: }
712: } else if (type == IS_LOCAL && is->ops->intervallocal) {
713: PetscBool isinterval;
715: PetscUseTypeMethod(is, intervallocal, &isinterval);
716: *flg = PETSC_FALSE;
717: if (isinterval) {
718: PetscInt min;
720: PetscCall(ISGetMinMax(is, &min, NULL));
721: if (min == 0) *flg = PETSC_TRUE;
722: }
723: } else {
724: PetscBool identLocal;
725: PetscInt n, i, rStart;
726: const PetscInt *idx;
728: PetscCall(ISGetLocalSize(is, &n));
729: PetscCall(ISGetIndices(is, &idx));
730: PetscCall(PetscLayoutGetRange(is->map, &rStart, NULL));
731: identLocal = PETSC_TRUE;
732: for (i = 0; i < n; i++) {
733: if (idx[i] != rStart + i) break;
734: }
735: if (i < n) identLocal = PETSC_FALSE;
736: if (type == IS_LOCAL || size == 1) {
737: *flg = identLocal;
738: } else {
739: PetscCallMPI(MPI_Allreduce(&identLocal, flg, 1, MPIU_BOOL, MPI_LAND, comm));
740: }
741: PetscCall(ISRestoreIndices(is, &idx));
742: }
743: PetscFunctionReturn(PETSC_SUCCESS);
744: }
746: /*@
747: ISGetInfo - Determine whether an index set satisfies a given property
749: Collective or logically collective on IS if the type is IS_GLOBAL (logically collective if the value of the property has been permanently set with ISSetInfo())
751: Input Parameters:
752: + is - the index set
753: . info - describing a property of the index set, one of those listed in the documentation of ISSetInfo()
754: . compute - if PETSC_FALSE, the property will not be computed if it is not already known and the property will be assumed to be false
755: - type - whether the property is local (IS_LOCAL) or global (IS_GLOBAL)
757: Output Parameter:
758: . flg - whether the property is true (PETSC_TRUE) or false (PETSC_FALSE)
760: Note: ISGetInfo uses cached values when possible, which will be incorrect if ISSetInfo() has been called with incorrect information. To clear cached values, use ISClearInfoCache().
762: Level: advanced
764: .seealso: `ISInfo`, `ISInfoType`, `ISSetInfo()`, `ISClearInfoCache()`
766: @*/
767: PetscErrorCode ISGetInfo(IS is, ISInfo info, ISInfoType type, PetscBool compute, PetscBool *flg)
768: {
769: MPI_Comm comm, errcomm;
770: PetscMPIInt rank, size;
771: PetscInt itype;
772: PetscBool hasprop;
773: PetscBool infer;
775: PetscFunctionBegin;
778: comm = PetscObjectComm((PetscObject)is);
779: if (type == IS_GLOBAL) {
781: errcomm = comm;
782: } else {
783: errcomm = PETSC_COMM_SELF;
784: }
786: PetscCallMPI(MPI_Comm_size(comm, &size));
787: PetscCallMPI(MPI_Comm_rank(comm, &rank));
789: PetscCheck(((int)info) > IS_INFO_MIN && ((int)info) < IS_INFO_MAX, errcomm, PETSC_ERR_ARG_OUTOFRANGE, "Options %d is out of range", (int)info);
790: if (size == 1) type = IS_LOCAL;
791: itype = (type == IS_LOCAL) ? 0 : 1;
792: hasprop = PETSC_FALSE;
793: infer = PETSC_FALSE;
794: if (is->info_permanent[itype][(int)info]) {
795: hasprop = (is->info[itype][(int)info] == IS_INFO_TRUE) ? PETSC_TRUE : PETSC_FALSE;
796: infer = PETSC_TRUE;
797: } else if ((itype == IS_LOCAL) && (is->info[IS_LOCAL][info] != IS_INFO_UNKNOWN)) {
798: /* we can cache local properties as long as we clear them when the IS changes */
799: /* NOTE: we only cache local values because there is no ISAssemblyBegin()/ISAssemblyEnd(),
800: so we have no way of knowing when a cached value has been invalidated by changes on a different process */
801: hasprop = (is->info[itype][(int)info] == IS_INFO_TRUE) ? PETSC_TRUE : PETSC_FALSE;
802: infer = PETSC_TRUE;
803: } else if (compute) {
804: switch (info) {
805: case IS_SORTED:
806: PetscCall(ISGetInfo_Sorted(is, type, &hasprop));
807: break;
808: case IS_UNIQUE:
809: PetscCall(ISGetInfo_Unique(is, type, &hasprop));
810: break;
811: case IS_PERMUTATION:
812: PetscCall(ISGetInfo_Permutation(is, type, &hasprop));
813: break;
814: case IS_INTERVAL:
815: PetscCall(ISGetInfo_Interval(is, type, &hasprop));
816: break;
817: case IS_IDENTITY:
818: PetscCall(ISGetInfo_Identity(is, type, &hasprop));
819: break;
820: default:
821: SETERRQ(errcomm, PETSC_ERR_ARG_OUTOFRANGE, "Unknown IS property");
822: }
823: infer = PETSC_TRUE;
824: }
825: /* call ISSetInfo_Internal to keep all of the implications straight */
826: if (infer) PetscCall(ISSetInfo_Internal(is, info, type, IS_INFO_UNKNOWN, hasprop));
827: *flg = hasprop;
828: PetscFunctionReturn(PETSC_SUCCESS);
829: }
831: static PetscErrorCode ISCopyInfo(IS source, IS dest)
832: {
833: PetscFunctionBegin;
834: PetscCall(PetscArraycpy(&dest->info[0], &source->info[0], 2));
835: PetscCall(PetscArraycpy(&dest->info_permanent[0], &source->info_permanent[0], 2));
836: PetscFunctionReturn(PETSC_SUCCESS);
837: }
839: /*@
840: ISIdentity - Determines whether index set is the identity mapping.
842: Collective
844: Input Parameters:
845: . is - the index set
847: Output Parameters:
848: . ident - PETSC_TRUE if an identity, else PETSC_FALSE
850: Level: intermediate
852: Note: If ISSetIdentity() (or ISSetInfo() for a permanent property) has been called,
853: ISIdentity() will return its answer without communication between processes, but
854: otherwise the output ident will be computed from ISGetInfo(),
855: which may require synchronization on the communicator of IS. To avoid this computation,
856: call ISGetInfo() directly with the compute flag set to PETSC_FALSE, and ident will be assumed false.
858: .seealso: `ISSetIdentity()`, `ISGetInfo()`
859: @*/
860: PetscErrorCode ISIdentity(IS is, PetscBool *ident)
861: {
862: PetscFunctionBegin;
865: PetscCall(ISGetInfo(is, IS_IDENTITY, IS_GLOBAL, PETSC_TRUE, ident));
866: PetscFunctionReturn(PETSC_SUCCESS);
867: }
869: /*@
870: ISSetIdentity - Informs the index set that it is an identity.
872: Logically Collective
874: Input Parameter:
875: . is - the index set
877: Level: intermediate
879: Note: The IS will be considered the identity permanently, even if indices have been changes (for example, with
880: ISGeneralSetIndices()). It's a good idea to only set this property if the IS will not change in the future.
881: To clear this property, use ISClearInfoCache().
883: .seealso: `ISIdentity()`, `ISSetInfo()`, `ISClearInfoCache()`
884: @*/
885: PetscErrorCode ISSetIdentity(IS is)
886: {
887: PetscFunctionBegin;
889: PetscCall(ISSetInfo(is, IS_IDENTITY, IS_GLOBAL, PETSC_TRUE, PETSC_TRUE));
890: PetscFunctionReturn(PETSC_SUCCESS);
891: }
893: /*@
894: ISContiguousLocal - Locates an index set with contiguous range within a global range, if possible
896: Not Collective
898: Input Parameters:
899: + is - the index set
900: . gstart - global start
901: - gend - global end
903: Output Parameters:
904: + start - start of contiguous block, as an offset from gstart
905: - contig - PETSC_TRUE if the index set refers to contiguous entries on this process, else PETSC_FALSE
907: Level: developer
909: .seealso: `ISGetLocalSize()`, `VecGetOwnershipRange()`
910: @*/
911: PetscErrorCode ISContiguousLocal(IS is, PetscInt gstart, PetscInt gend, PetscInt *start, PetscBool *contig)
912: {
913: PetscFunctionBegin;
917: *start = -1;
918: *contig = PETSC_FALSE;
919: PetscTryTypeMethod(is, contiguous, gstart, gend, start, contig);
920: PetscFunctionReturn(PETSC_SUCCESS);
921: }
923: /*@
924: ISPermutation - PETSC_TRUE or PETSC_FALSE depending on whether the
925: index set has been declared to be a permutation.
927: Logically Collective
929: Input Parameter:
930: . is - the index set
932: Output Parameter:
933: . perm - PETSC_TRUE if a permutation, else PETSC_FALSE
935: Level: intermediate
937: Note: If it is not already known that the IS is a permutation (if ISSetPermutation()
938: or ISSetInfo() has not been called), this routine will not attempt to compute
939: whether the index set is a permutation and will assume perm is PETSC_FALSE.
940: To compute the value when it is not already known, use ISGetInfo() with
941: the compute flag set to PETSC_TRUE.
943: .seealso: `ISSetPermutation()`, `ISGetInfo()`
944: @*/
945: PetscErrorCode ISPermutation(IS is, PetscBool *perm)
946: {
947: PetscFunctionBegin;
950: PetscCall(ISGetInfo(is, IS_PERMUTATION, IS_GLOBAL, PETSC_FALSE, perm));
951: PetscFunctionReturn(PETSC_SUCCESS);
952: }
954: /*@
955: ISSetPermutation - Informs the index set that it is a permutation.
957: Logically Collective
959: Input Parameter:
960: . is - the index set
962: Level: intermediate
964: The debug version of the libraries (./configure --with-debugging=1) checks if the
965: index set is actually a permutation. The optimized version just believes you.
967: Note: The IS will be considered a permutation permanently, even if indices have been changes (for example, with
968: ISGeneralSetIndices()). It's a good idea to only set this property if the IS will not change in the future.
969: To clear this property, use ISClearInfoCache().
971: .seealso: `ISPermutation()`, `ISSetInfo()`, `ISClearInfoCache().`
972: @*/
973: PetscErrorCode ISSetPermutation(IS is)
974: {
975: PetscFunctionBegin;
977: if (PetscDefined(USE_DEBUG)) {
978: PetscMPIInt size;
980: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is), &size));
981: if (size == 1) {
982: PetscInt i, n, *idx;
983: const PetscInt *iidx;
985: PetscCall(ISGetSize(is, &n));
986: PetscCall(PetscMalloc1(n, &idx));
987: PetscCall(ISGetIndices(is, &iidx));
988: PetscCall(PetscArraycpy(idx, iidx, n));
989: PetscCall(PetscIntSortSemiOrdered(n, idx));
990: for (i = 0; i < n; i++) PetscCheck(idx[i] == i, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index set is not a permutation");
991: PetscCall(PetscFree(idx));
992: PetscCall(ISRestoreIndices(is, &iidx));
993: }
994: }
995: PetscCall(ISSetInfo(is, IS_PERMUTATION, IS_GLOBAL, PETSC_TRUE, PETSC_TRUE));
996: PetscFunctionReturn(PETSC_SUCCESS);
997: }
999: /*@C
1000: ISDestroy - Destroys an index set.
1002: Collective
1004: Input Parameters:
1005: . is - the index set
1007: Level: beginner
1009: .seealso: `ISCreateGeneral()`, `ISCreateStride()`, `ISCreateBlocked()`
1010: @*/
1011: PetscErrorCode ISDestroy(IS *is)
1012: {
1013: PetscFunctionBegin;
1014: if (!*is) PetscFunctionReturn(PETSC_SUCCESS);
1016: if (--((PetscObject)(*is))->refct > 0) {
1017: *is = NULL;
1018: PetscFunctionReturn(PETSC_SUCCESS);
1019: }
1020: if ((*is)->complement) {
1021: PetscInt refcnt;
1022: PetscCall(PetscObjectGetReference((PetscObject)((*is)->complement), &refcnt));
1023: PetscCheck(refcnt <= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Nonlocal IS has not been restored");
1024: PetscCall(ISDestroy(&(*is)->complement));
1025: }
1026: if ((*is)->ops->destroy) PetscCall((*(*is)->ops->destroy)(*is));
1027: PetscCall(PetscLayoutDestroy(&(*is)->map));
1028: /* Destroy local representations of offproc data. */
1029: PetscCall(PetscFree((*is)->total));
1030: PetscCall(PetscFree((*is)->nonlocal));
1031: PetscCall(PetscHeaderDestroy(is));
1032: PetscFunctionReturn(PETSC_SUCCESS);
1033: }
1035: /*@
1036: ISInvertPermutation - Creates a new permutation that is the inverse of
1037: a given permutation.
1039: Collective
1041: Input Parameters:
1042: + is - the index set
1043: - nlocal - number of indices on this processor in result (ignored for 1 processor) or
1044: use PETSC_DECIDE
1046: Output Parameter:
1047: . isout - the inverse permutation
1049: Level: intermediate
1051: Notes:
1052: For parallel index sets this does the complete parallel permutation, but the
1053: code is not efficient for huge index sets (10,000,000 indices).
1055: @*/
1056: PetscErrorCode ISInvertPermutation(IS is, PetscInt nlocal, IS *isout)
1057: {
1058: PetscBool isperm, isidentity, issame;
1060: PetscFunctionBegin;
1063: PetscCall(ISGetInfo(is, IS_PERMUTATION, IS_GLOBAL, PETSC_TRUE, &isperm));
1064: PetscCheck(isperm, PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_WRONG, "Not a permutation");
1065: PetscCall(ISGetInfo(is, IS_IDENTITY, IS_GLOBAL, PETSC_TRUE, &isidentity));
1066: issame = PETSC_FALSE;
1067: if (isidentity) {
1068: PetscInt n;
1069: PetscBool isallsame;
1071: PetscCall(ISGetLocalSize(is, &n));
1072: issame = (PetscBool)(n == nlocal);
1073: PetscCallMPI(MPI_Allreduce(&issame, &isallsame, 1, MPIU_BOOL, MPI_LAND, PetscObjectComm((PetscObject)is)));
1074: issame = isallsame;
1075: }
1076: if (issame) {
1077: PetscCall(ISDuplicate(is, isout));
1078: } else {
1079: PetscUseTypeMethod(is, invertpermutation, nlocal, isout);
1080: PetscCall(ISSetPermutation(*isout));
1081: }
1082: PetscFunctionReturn(PETSC_SUCCESS);
1083: }
1085: /*@
1086: ISGetSize - Returns the global length of an index set.
1088: Not Collective
1090: Input Parameter:
1091: . is - the index set
1093: Output Parameter:
1094: . size - the global size
1096: Level: beginner
1098: @*/
1099: PetscErrorCode ISGetSize(IS is, PetscInt *size)
1100: {
1101: PetscFunctionBegin;
1104: *size = is->map->N;
1105: PetscFunctionReturn(PETSC_SUCCESS);
1106: }
1108: /*@
1109: ISGetLocalSize - Returns the local (processor) length of an index set.
1111: Not Collective
1113: Input Parameter:
1114: . is - the index set
1116: Output Parameter:
1117: . size - the local size
1119: Level: beginner
1121: @*/
1122: PetscErrorCode ISGetLocalSize(IS is, PetscInt *size)
1123: {
1124: PetscFunctionBegin;
1127: *size = is->map->n;
1128: PetscFunctionReturn(PETSC_SUCCESS);
1129: }
1131: /*@
1132: ISGetLayout - get PetscLayout describing index set layout
1134: Not Collective
1136: Input Parameter:
1137: . is - the index set
1139: Output Parameter:
1140: . map - the layout
1142: Level: developer
1144: .seealso: `ISSetLayout()`, `ISGetSize()`, `ISGetLocalSize()`
1145: @*/
1146: PetscErrorCode ISGetLayout(IS is, PetscLayout *map)
1147: {
1148: PetscFunctionBegin;
1151: *map = is->map;
1152: PetscFunctionReturn(PETSC_SUCCESS);
1153: }
1155: /*@
1156: ISSetLayout - set PetscLayout describing index set layout
1158: Collective
1160: Input Arguments:
1161: + is - the index set
1162: - map - the layout
1164: Level: developer
1166: Notes:
1167: Users should typically use higher level functions such as ISCreateGeneral().
1169: This function can be useful in some special cases of constructing a new IS, e.g. after ISCreate() and before ISLoad().
1170: Otherwise, it is only valid to replace the layout with a layout known to be equivalent.
1172: .seealso: `ISCreate()`, `ISGetLayout()`, `ISGetSize()`, `ISGetLocalSize()`
1173: @*/
1174: PetscErrorCode ISSetLayout(IS is, PetscLayout map)
1175: {
1176: PetscFunctionBegin;
1179: PetscCall(PetscLayoutReference(map, &is->map));
1180: PetscFunctionReturn(PETSC_SUCCESS);
1181: }
1183: /*@C
1184: ISGetIndices - Returns a pointer to the indices. The user should call
1185: `ISRestoreIndices()` after having looked at the indices. The user should
1186: NOT change the indices.
1188: Not Collective
1190: Input Parameter:
1191: . is - the index set
1193: Output Parameter:
1194: . ptr - the location to put the pointer to the indices
1196: Level: intermediate
1198: Fortran Note:
1199: `ISGetIndices()` Fortran binding is deprecated (since PETSc 3.19), use `ISGetIndicesF90()`
1201: .seealso: `ISRestoreIndices()`, `ISGetIndicesF90()`
1202: @*/
1203: PetscErrorCode ISGetIndices(IS is, const PetscInt *ptr[])
1204: {
1205: PetscFunctionBegin;
1208: PetscUseTypeMethod(is, getindices, ptr);
1209: PetscFunctionReturn(PETSC_SUCCESS);
1210: }
1212: /*@C
1213: ISGetMinMax - Gets the minimum and maximum values in an IS
1215: Not Collective
1217: Input Parameter:
1218: . is - the index set
1220: Output Parameters:
1221: + min - the minimum value
1222: - max - the maximum value
1224: Level: intermediate
1226: Notes:
1227: Empty index sets return min=PETSC_MAX_INT and max=PETSC_MIN_INT.
1228: In parallel, it returns the min and max of the local portion of the IS
1230: .seealso: `ISGetIndices()`, `ISRestoreIndices()`, `ISGetIndicesF90()`
1231: @*/
1232: PetscErrorCode ISGetMinMax(IS is, PetscInt *min, PetscInt *max)
1233: {
1234: PetscFunctionBegin;
1236: if (min) *min = is->min;
1237: if (max) *max = is->max;
1238: PetscFunctionReturn(PETSC_SUCCESS);
1239: }
1241: /*@
1242: ISLocate - determine the location of an index within the local component of an index set
1244: Not Collective
1246: Input Parameters:
1247: + is - the index set
1248: - key - the search key
1250: Output Parameter:
1251: . location - if >= 0, a location within the index set that is equal to the key, otherwise the key is not in the index set
1253: Level: intermediate
1254: @*/
1255: PetscErrorCode ISLocate(IS is, PetscInt key, PetscInt *location)
1256: {
1257: PetscFunctionBegin;
1258: if (is->ops->locate) {
1259: PetscUseTypeMethod(is, locate, key, location);
1260: } else {
1261: PetscInt numIdx;
1262: PetscBool sorted;
1263: const PetscInt *idx;
1265: PetscCall(ISGetLocalSize(is, &numIdx));
1266: PetscCall(ISGetIndices(is, &idx));
1267: PetscCall(ISSorted(is, &sorted));
1268: if (sorted) {
1269: PetscCall(PetscFindInt(key, numIdx, idx, location));
1270: } else {
1271: PetscInt i;
1273: *location = -1;
1274: for (i = 0; i < numIdx; i++) {
1275: if (idx[i] == key) {
1276: *location = i;
1277: break;
1278: }
1279: }
1280: }
1281: PetscCall(ISRestoreIndices(is, &idx));
1282: }
1283: PetscFunctionReturn(PETSC_SUCCESS);
1284: }
1286: /*@C
1287: ISRestoreIndices - Restores an index set to a usable state after a call to `ISGetIndices()`.
1289: Not Collective
1291: Input Parameters:
1292: + is - the index set
1293: - ptr - the pointer obtained by ISGetIndices()
1295: Level: intermediate
1297: Note:
1298: This routine zeros out ptr. This is to prevent accidental use of the array after it has been restored.
1300: Fortran Note:
1301: `ISRestoreIndices()` Fortran binding is deprecated (since PETSc 3.19), use `ISRestoreIndicesF90()`
1303: .seealso: `ISGetIndices()`, `ISRestoreIndicesF90()`
1304: @*/
1305: PetscErrorCode ISRestoreIndices(IS is, const PetscInt *ptr[])
1306: {
1307: PetscFunctionBegin;
1310: PetscTryTypeMethod(is, restoreindices, ptr);
1311: PetscFunctionReturn(PETSC_SUCCESS);
1312: }
1314: static PetscErrorCode ISGatherTotal_Private(IS is)
1315: {
1316: PetscInt i, n, N;
1317: const PetscInt *lindices;
1318: MPI_Comm comm;
1319: PetscMPIInt rank, size, *sizes = NULL, *offsets = NULL, nn;
1321: PetscFunctionBegin;
1324: PetscCall(PetscObjectGetComm((PetscObject)is, &comm));
1325: PetscCallMPI(MPI_Comm_size(comm, &size));
1326: PetscCallMPI(MPI_Comm_rank(comm, &rank));
1327: PetscCall(ISGetLocalSize(is, &n));
1328: PetscCall(PetscMalloc2(size, &sizes, size, &offsets));
1330: PetscCall(PetscMPIIntCast(n, &nn));
1331: PetscCallMPI(MPI_Allgather(&nn, 1, MPI_INT, sizes, 1, MPI_INT, comm));
1332: offsets[0] = 0;
1333: for (i = 1; i < size; ++i) offsets[i] = offsets[i - 1] + sizes[i - 1];
1334: N = offsets[size - 1] + sizes[size - 1];
1336: PetscCall(PetscMalloc1(N, &(is->total)));
1337: PetscCall(ISGetIndices(is, &lindices));
1338: PetscCallMPI(MPI_Allgatherv((void *)lindices, nn, MPIU_INT, is->total, sizes, offsets, MPIU_INT, comm));
1339: PetscCall(ISRestoreIndices(is, &lindices));
1340: is->local_offset = offsets[rank];
1341: PetscCall(PetscFree2(sizes, offsets));
1342: PetscFunctionReturn(PETSC_SUCCESS);
1343: }
1345: /*@C
1346: ISGetTotalIndices - Retrieve an array containing all indices across the communicator.
1348: Collective
1350: Input Parameter:
1351: . is - the index set
1353: Output Parameter:
1354: . indices - total indices with rank 0 indices first, and so on; total array size is
1355: the same as returned with ISGetSize().
1357: Level: intermediate
1359: Notes:
1360: this is potentially nonscalable, but depends on the size of the total index set
1361: and the size of the communicator. This may be feasible for index sets defined on
1362: subcommunicators, such that the set size does not grow with PETSC_WORLD_COMM.
1363: Note also that there is no way to tell where the local part of the indices starts
1364: (use ISGetIndices() and ISGetNonlocalIndices() to retrieve just the local and just
1365: the nonlocal part (complement), respectively).
1367: .seealso: `ISRestoreTotalIndices()`, `ISGetNonlocalIndices()`, `ISGetSize()`
1368: @*/
1369: PetscErrorCode ISGetTotalIndices(IS is, const PetscInt *indices[])
1370: {
1371: PetscMPIInt size;
1373: PetscFunctionBegin;
1376: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is), &size));
1377: if (size == 1) {
1378: PetscUseTypeMethod(is, getindices, indices);
1379: } else {
1380: if (!is->total) PetscCall(ISGatherTotal_Private(is));
1381: *indices = is->total;
1382: }
1383: PetscFunctionReturn(PETSC_SUCCESS);
1384: }
1386: /*@C
1387: ISRestoreTotalIndices - Restore the index array obtained with ISGetTotalIndices().
1389: Not Collective.
1391: Input Parameters:
1392: + is - the index set
1393: - indices - index array; must be the array obtained with ISGetTotalIndices()
1395: Level: intermediate
1397: .seealso: `ISRestoreTotalIndices()`, `ISGetNonlocalIndices()`
1398: @*/
1399: PetscErrorCode ISRestoreTotalIndices(IS is, const PetscInt *indices[])
1400: {
1401: PetscMPIInt size;
1403: PetscFunctionBegin;
1406: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is), &size));
1407: if (size == 1) {
1408: PetscUseTypeMethod(is, restoreindices, indices);
1409: } else {
1410: PetscCheck(is->total == *indices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
1411: }
1412: PetscFunctionReturn(PETSC_SUCCESS);
1413: }
1415: /*@C
1416: ISGetNonlocalIndices - Retrieve an array of indices from remote processors
1417: in this communicator.
1419: Collective
1421: Input Parameter:
1422: . is - the index set
1424: Output Parameter:
1425: . indices - indices with rank 0 indices first, and so on, omitting
1426: the current rank. Total number of indices is the difference
1427: total and local, obtained with ISGetSize() and ISGetLocalSize(),
1428: respectively.
1430: Level: intermediate
1432: Notes:
1433: restore the indices using ISRestoreNonlocalIndices().
1434: The same scalability considerations as those for ISGetTotalIndices
1435: apply here.
1437: .seealso: `ISGetTotalIndices()`, `ISRestoreNonlocalIndices()`, `ISGetSize()`, `ISGetLocalSize().`
1438: @*/
1439: PetscErrorCode ISGetNonlocalIndices(IS is, const PetscInt *indices[])
1440: {
1441: PetscMPIInt size;
1442: PetscInt n, N;
1444: PetscFunctionBegin;
1447: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)is), &size));
1448: if (size == 1) *indices = NULL;
1449: else {
1450: if (!is->total) PetscCall(ISGatherTotal_Private(is));
1451: PetscCall(ISGetLocalSize(is, &n));
1452: PetscCall(ISGetSize(is, &N));
1453: PetscCall(PetscMalloc1(N - n, &(is->nonlocal)));
1454: PetscCall(PetscArraycpy(is->nonlocal, is->total, is->local_offset));
1455: PetscCall(PetscArraycpy(is->nonlocal + is->local_offset, is->total + is->local_offset + n, N - is->local_offset - n));
1456: *indices = is->nonlocal;
1457: }
1458: PetscFunctionReturn(PETSC_SUCCESS);
1459: }
1461: /*@C
1462: ISRestoreNonlocalIndices - Restore the index array obtained with ISGetNonlocalIndices().
1464: Not Collective.
1466: Input Parameters:
1467: + is - the index set
1468: - indices - index array; must be the array obtained with ISGetNonlocalIndices()
1470: Level: intermediate
1472: .seealso: `ISGetTotalIndices()`, `ISGetNonlocalIndices()`, `ISRestoreTotalIndices()`
1473: @*/
1474: PetscErrorCode ISRestoreNonlocalIndices(IS is, const PetscInt *indices[])
1475: {
1476: PetscFunctionBegin;
1479: PetscCheck(is->nonlocal == *indices, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Index array pointer being restored does not point to the array obtained from the IS.");
1480: PetscFunctionReturn(PETSC_SUCCESS);
1481: }
1483: /*@
1484: ISGetNonlocalIS - Gather all nonlocal indices for this IS and present
1485: them as another sequential index set.
1487: Collective
1489: Input Parameter:
1490: . is - the index set
1492: Output Parameter:
1493: . complement - sequential IS with indices identical to the result of
1494: ISGetNonlocalIndices()
1496: Level: intermediate
1498: Notes:
1499: complement represents the result of ISGetNonlocalIndices as an IS.
1500: Therefore scalability issues similar to ISGetNonlocalIndices apply.
1501: The resulting IS must be restored using ISRestoreNonlocalIS().
1503: .seealso: `ISGetNonlocalIndices()`, `ISRestoreNonlocalIndices()`, `ISAllGather()`, `ISGetSize()`
1504: @*/
1505: PetscErrorCode ISGetNonlocalIS(IS is, IS *complement)
1506: {
1507: PetscFunctionBegin;
1510: /* Check if the complement exists already. */
1511: if (is->complement) {
1512: *complement = is->complement;
1513: PetscCall(PetscObjectReference((PetscObject)(is->complement)));
1514: } else {
1515: PetscInt N, n;
1516: const PetscInt *idx;
1517: PetscCall(ISGetSize(is, &N));
1518: PetscCall(ISGetLocalSize(is, &n));
1519: PetscCall(ISGetNonlocalIndices(is, &idx));
1520: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, N - n, idx, PETSC_USE_POINTER, &(is->complement)));
1521: PetscCall(PetscObjectReference((PetscObject)is->complement));
1522: *complement = is->complement;
1523: }
1524: PetscFunctionReturn(PETSC_SUCCESS);
1525: }
1527: /*@
1528: ISRestoreNonlocalIS - Restore the IS obtained with ISGetNonlocalIS().
1530: Not collective.
1532: Input Parameters:
1533: + is - the index set
1534: - complement - index set of is's nonlocal indices
1536: Level: intermediate
1538: .seealso: `ISGetNonlocalIS()`, `ISGetNonlocalIndices()`, `ISRestoreNonlocalIndices()`
1539: @*/
1540: PetscErrorCode ISRestoreNonlocalIS(IS is, IS *complement)
1541: {
1542: PetscInt refcnt;
1544: PetscFunctionBegin;
1547: PetscCheck(*complement == is->complement, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Complement IS being restored was not obtained with ISGetNonlocalIS()");
1548: PetscCall(PetscObjectGetReference((PetscObject)(is->complement), &refcnt));
1549: PetscCheck(refcnt > 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Duplicate call to ISRestoreNonlocalIS() detected");
1550: PetscCall(PetscObjectDereference((PetscObject)(is->complement)));
1551: PetscFunctionReturn(PETSC_SUCCESS);
1552: }
1554: /*@C
1555: ISViewFromOptions - View from Options
1557: Collective
1559: Input Parameters:
1560: + A - the index set
1561: . obj - Optional object
1562: - name - command line option
1564: Level: intermediate
1565: .seealso: `IS`, `ISView`, `PetscObjectViewFromOptions()`, `ISCreate()`
1566: @*/
1567: PetscErrorCode ISViewFromOptions(IS A, PetscObject obj, const char name[])
1568: {
1569: PetscFunctionBegin;
1571: PetscCall(PetscObjectViewFromOptions((PetscObject)A, obj, name));
1572: PetscFunctionReturn(PETSC_SUCCESS);
1573: }
1575: /*@C
1576: ISView - Displays an index set.
1578: Collective
1580: Input Parameters:
1581: + is - the index set
1582: - viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF.
1584: Level: intermediate
1586: .seealso: `PetscViewerASCIIOpen()`
1587: @*/
1588: PetscErrorCode ISView(IS is, PetscViewer viewer)
1589: {
1590: PetscFunctionBegin;
1592: if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)is), &viewer));
1594: PetscCheckSameComm(is, 1, viewer, 2);
1596: PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)is, viewer));
1597: PetscCall(PetscLogEventBegin(IS_View, is, viewer, 0, 0));
1598: PetscUseTypeMethod(is, view, viewer);
1599: PetscCall(PetscLogEventEnd(IS_View, is, viewer, 0, 0));
1600: PetscFunctionReturn(PETSC_SUCCESS);
1601: }
1603: /*@
1604: ISLoad - Loads a vector that has been stored in binary or HDF5 format with ISView().
1606: Collective
1608: Input Parameters:
1609: + is - the newly loaded vector, this needs to have been created with ISCreate() or some related function before a call to ISLoad().
1610: - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or HDF5 file viewer, obtained from PetscViewerHDF5Open()
1612: Level: intermediate
1614: Notes:
1615: IF using HDF5, you must assign the IS the same name as was used in the IS
1616: that was stored in the file using PetscObjectSetName(). Otherwise you will
1617: get the error message: "Cannot H5DOpen2() with Vec name NAMEOFOBJECT"
1619: .seealso: `PetscViewerBinaryOpen()`, `ISView()`, `MatLoad()`, `VecLoad()`
1620: @*/
1621: PetscErrorCode ISLoad(IS is, PetscViewer viewer)
1622: {
1623: PetscBool isbinary, ishdf5;
1625: PetscFunctionBegin;
1628: PetscCheckSameComm(is, 1, viewer, 2);
1629: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &isbinary));
1630: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERHDF5, &ishdf5));
1631: PetscCheck(isbinary || ishdf5, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Invalid viewer; open viewer with PetscViewerBinaryOpen()");
1632: if (!((PetscObject)is)->type_name) PetscCall(ISSetType(is, ISGENERAL));
1633: PetscCall(PetscLogEventBegin(IS_Load, is, viewer, 0, 0));
1634: PetscUseTypeMethod(is, load, viewer);
1635: PetscCall(PetscLogEventEnd(IS_Load, is, viewer, 0, 0));
1636: PetscFunctionReturn(PETSC_SUCCESS);
1637: }
1639: /*@
1640: ISSort - Sorts the indices of an index set.
1642: Collective
1644: Input Parameters:
1645: . is - the index set
1647: Level: intermediate
1649: .seealso: `ISSortRemoveDups()`, `ISSorted()`
1650: @*/
1651: PetscErrorCode ISSort(IS is)
1652: {
1653: PetscFunctionBegin;
1655: PetscUseTypeMethod(is, sort);
1656: PetscCall(ISSetInfo(is, IS_SORTED, IS_LOCAL, is->info_permanent[IS_LOCAL][IS_SORTED], PETSC_TRUE));
1657: PetscFunctionReturn(PETSC_SUCCESS);
1658: }
1660: /*@
1661: ISSortRemoveDups - Sorts the indices of an index set, removing duplicates.
1663: Collective
1665: Input Parameters:
1666: . is - the index set
1668: Level: intermediate
1670: .seealso: `ISSort()`, `ISSorted()`
1671: @*/
1672: PetscErrorCode ISSortRemoveDups(IS is)
1673: {
1674: PetscFunctionBegin;
1676: PetscCall(ISClearInfoCache(is, PETSC_FALSE));
1677: PetscUseTypeMethod(is, sortremovedups);
1678: PetscCall(ISSetInfo(is, IS_SORTED, IS_LOCAL, is->info_permanent[IS_LOCAL][IS_SORTED], PETSC_TRUE));
1679: PetscCall(ISSetInfo(is, IS_UNIQUE, IS_LOCAL, is->info_permanent[IS_LOCAL][IS_UNIQUE], PETSC_TRUE));
1680: PetscFunctionReturn(PETSC_SUCCESS);
1681: }
1683: /*@
1684: ISToGeneral - Converts an IS object of any type to ISGENERAL type
1686: Collective
1688: Input Parameters:
1689: . is - the index set
1691: Level: intermediate
1693: .seealso: `ISSorted()`
1694: @*/
1695: PetscErrorCode ISToGeneral(IS is)
1696: {
1697: PetscFunctionBegin;
1699: PetscUseTypeMethod(is, togeneral);
1700: PetscFunctionReturn(PETSC_SUCCESS);
1701: }
1703: /*@
1704: ISSorted - Checks the indices to determine whether they have been sorted.
1706: Collective
1708: Input Parameter:
1709: . is - the index set
1711: Output Parameter:
1712: . flg - output flag, either PETSC_TRUE if the index set is sorted,
1713: or PETSC_FALSE otherwise.
1715: Notes:
1716: For parallel IS objects this only indicates if the local part of the IS
1717: is sorted. So some processors may return PETSC_TRUE while others may
1718: return PETSC_FALSE.
1720: Level: intermediate
1722: .seealso: `ISSort()`, `ISSortRemoveDups()`
1723: @*/
1724: PetscErrorCode ISSorted(IS is, PetscBool *flg)
1725: {
1726: PetscFunctionBegin;
1729: PetscCall(ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, flg));
1730: PetscFunctionReturn(PETSC_SUCCESS);
1731: }
1733: /*@
1734: ISDuplicate - Creates a duplicate copy of an index set.
1736: Collective
1738: Input Parameter:
1739: . is - the index set
1741: Output Parameter:
1742: . isnew - the copy of the index set
1744: Level: beginner
1746: .seealso: `ISCreateGeneral()`, `ISCopy()`
1747: @*/
1748: PetscErrorCode ISDuplicate(IS is, IS *newIS)
1749: {
1750: PetscFunctionBegin;
1753: PetscUseTypeMethod(is, duplicate, newIS);
1754: PetscCall(ISCopyInfo(is, *newIS));
1755: PetscFunctionReturn(PETSC_SUCCESS);
1756: }
1758: /*@
1759: ISCopy - Copies an index set.
1761: Collective
1763: Input Parameter:
1764: . is - the index set
1766: Output Parameter:
1767: . isy - the copy of the index set
1769: Level: beginner
1771: .seealso: `ISDuplicate()`, `ISShift()`
1772: @*/
1773: PetscErrorCode ISCopy(IS is, IS isy)
1774: {
1775: PetscInt bs, bsy;
1777: PetscFunctionBegin;
1780: PetscCheckSameComm(is, 1, isy, 2);
1781: if (is == isy) PetscFunctionReturn(PETSC_SUCCESS);
1782: PetscCall(PetscLayoutGetBlockSize(is->map, &bs));
1783: PetscCall(PetscLayoutGetBlockSize(isy->map, &bsy));
1784: PetscCheck(is->map->N == isy->map->N, PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_INCOMP, "Index sets have different global size %" PetscInt_FMT " != %" PetscInt_FMT, is->map->N, isy->map->N);
1785: PetscCheck(is->map->n == isy->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Index sets have different local size %" PetscInt_FMT " != %" PetscInt_FMT, is->map->n, isy->map->n);
1786: PetscCheck(bs == bsy, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Index sets have different block size %" PetscInt_FMT " != %" PetscInt_FMT, bs, bsy);
1787: PetscCall(ISCopyInfo(is, isy));
1788: isy->max = is->max;
1789: isy->min = is->min;
1790: PetscUseTypeMethod(is, copy, isy);
1791: PetscFunctionReturn(PETSC_SUCCESS);
1792: }
1794: /*@
1795: ISShift - Shift all indices by given offset
1797: Collective
1799: Input Parameters:
1800: + is - the index set
1801: - offset - the offset
1803: Output Parameter:
1804: . isy - the shifted copy of the input index set
1806: Notes:
1807: The offset can be different across processes.
1808: IS is and isy can be the same.
1810: Level: beginner
1812: .seealso: `ISDuplicate()`, `ISCopy()`
1813: @*/
1814: PetscErrorCode ISShift(IS is, PetscInt offset, IS isy)
1815: {
1816: PetscFunctionBegin;
1819: PetscCheckSameComm(is, 1, isy, 3);
1820: if (!offset) {
1821: PetscCall(ISCopy(is, isy));
1822: PetscFunctionReturn(PETSC_SUCCESS);
1823: }
1824: PetscCheck(is->map->N == isy->map->N, PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_INCOMP, "Index sets have different global size %" PetscInt_FMT " != %" PetscInt_FMT, is->map->N, isy->map->N);
1825: PetscCheck(is->map->n == isy->map->n, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Index sets have different local size %" PetscInt_FMT " != %" PetscInt_FMT, is->map->n, isy->map->n);
1826: PetscCheck(is->map->bs == isy->map->bs, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Index sets have different block size %" PetscInt_FMT " != %" PetscInt_FMT, is->map->bs, isy->map->bs);
1827: PetscCall(ISCopyInfo(is, isy));
1828: isy->max = is->max + offset;
1829: isy->min = is->min + offset;
1830: PetscUseMethod(is, "ISShift_C", (IS, PetscInt, IS), (is, offset, isy));
1831: PetscFunctionReturn(PETSC_SUCCESS);
1832: }
1834: /*@
1835: ISOnComm - Split a parallel IS on subcomms (usually self) or concatenate index sets on subcomms into a parallel index set
1837: Collective
1839: Input Parameters:
1840: + is - index set
1841: . comm - communicator for new index set
1842: - mode - copy semantics, PETSC_USE_POINTER for no-copy if possible, otherwise PETSC_COPY_VALUES
1844: Output Parameter:
1845: . newis - new IS on comm
1847: Level: advanced
1849: Notes:
1850: It is usually desirable to create a parallel IS and look at the local part when necessary.
1852: This function is useful if serial ISs must be created independently, or to view many
1853: logically independent serial ISs.
1855: The input IS must have the same type on every process.
1856: @*/
1857: PetscErrorCode ISOnComm(IS is, MPI_Comm comm, PetscCopyMode mode, IS *newis)
1858: {
1859: PetscMPIInt match;
1861: PetscFunctionBegin;
1864: PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)is), comm, &match));
1865: if (mode != PETSC_COPY_VALUES && (match == MPI_IDENT || match == MPI_CONGRUENT)) {
1866: PetscCall(PetscObjectReference((PetscObject)is));
1867: *newis = is;
1868: } else PetscUseTypeMethod(is, oncomm, comm, mode, newis);
1869: PetscFunctionReturn(PETSC_SUCCESS);
1870: }
1872: /*@
1873: ISSetBlockSize - informs an index set that it has a given block size
1875: Logicall Collective
1877: Input Parameters:
1878: + is - index set
1879: - bs - block size
1881: Level: intermediate
1883: Notes:
1884: This is much like the block size for Vecs. It indicates that one can think of the indices as
1885: being in a collection of equal size blocks. For ISBlock() these collections of blocks are all contiquous
1886: within a block but this is not the case for other IS.
1887: ISBlockGetIndices() only works for ISBlock IS, not others.
1889: .seealso: `ISGetBlockSize()`, `ISCreateBlock()`, `ISBlockGetIndices()`,
1890: @*/
1891: PetscErrorCode ISSetBlockSize(IS is, PetscInt bs)
1892: {
1893: PetscFunctionBegin;
1896: PetscCheck(bs >= 1, PetscObjectComm((PetscObject)is), PETSC_ERR_ARG_OUTOFRANGE, "Block size %" PetscInt_FMT ", must be positive", bs);
1897: if (PetscDefined(USE_DEBUG)) {
1898: const PetscInt *indices;
1899: PetscInt length, i, j;
1900: PetscCall(ISGetIndices(is, &indices));
1901: if (indices) {
1902: PetscCall(ISGetLocalSize(is, &length));
1903: PetscCheck(length % bs == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Local size %" PetscInt_FMT " not compatible with block size %" PetscInt_FMT, length, bs);
1904: for (i = 0; i < length / bs; i += bs) {
1905: for (j = 0; j < bs - 1; j++) {
1906: PetscCheck(indices[i * bs + j] == indices[i * bs + j + 1] - 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Block size %" PetscInt_FMT " is incompatible with the indices: non consecutive indices %" PetscInt_FMT " %" PetscInt_FMT, bs, indices[i * bs + j], indices[i * bs + j + 1]);
1907: }
1908: }
1909: }
1910: PetscCall(ISRestoreIndices(is, &indices));
1911: }
1912: PetscUseTypeMethod(is, setblocksize, bs);
1913: PetscFunctionReturn(PETSC_SUCCESS);
1914: }
1916: /*@
1917: ISGetBlockSize - Returns the number of elements in a block.
1919: Not Collective
1921: Input Parameter:
1922: . is - the index set
1924: Output Parameter:
1925: . size - the number of elements in a block
1927: Level: intermediate
1929: Notes:
1930: This is much like the block size for Vecs. It indicates that one can think of the indices as
1931: being in a collection of equal size blocks. For ISBlock() these collections of blocks are all contiquous
1932: within a block but this is not the case for other IS.
1933: ISBlockGetIndices() only works for ISBlock IS, not others.
1935: .seealso: `ISBlockGetSize()`, `ISGetSize()`, `ISCreateBlock()`, `ISSetBlockSize()`
1936: @*/
1937: PetscErrorCode ISGetBlockSize(IS is, PetscInt *size)
1938: {
1939: PetscFunctionBegin;
1940: PetscCall(PetscLayoutGetBlockSize(is->map, size));
1941: PetscFunctionReturn(PETSC_SUCCESS);
1942: }
1944: PetscErrorCode ISGetIndicesCopy(IS is, PetscInt idx[])
1945: {
1946: PetscInt len, i;
1947: const PetscInt *ptr;
1949: PetscFunctionBegin;
1950: PetscCall(ISGetLocalSize(is, &len));
1951: PetscCall(ISGetIndices(is, &ptr));
1952: for (i = 0; i < len; i++) idx[i] = ptr[i];
1953: PetscCall(ISRestoreIndices(is, &ptr));
1954: PetscFunctionReturn(PETSC_SUCCESS);
1955: }
1957: /*MC
1958: ISGetIndicesF90 - Accesses the elements of an index set from Fortran90.
1959: The users should call ISRestoreIndicesF90() after having looked at the
1960: indices. The user should NOT change the indices.
1962: Synopsis:
1963: ISGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
1965: Not collective
1967: Input Parameter:
1968: . x - index set
1970: Output Parameters:
1971: + xx_v - the Fortran90 pointer to the array
1972: - ierr - error code
1974: Example of Usage:
1975: .vb
1976: PetscInt, pointer xx_v(:)
1977: ....
1978: call ISGetIndicesF90(x,xx_v,ierr)
1979: a = xx_v(3)
1980: call ISRestoreIndicesF90(x,xx_v,ierr)
1981: .ve
1983: Level: intermediate
1985: .seealso: `ISRestoreIndicesF90()`, `ISGetIndices()`, `ISRestoreIndices()`
1987: M*/
1989: /*MC
1990: ISRestoreIndicesF90 - Restores an index set to a usable state after
1991: a call to ISGetIndicesF90().
1993: Synopsis:
1994: ISRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
1996: Not collective
1998: Input Parameters:
1999: + x - index set
2000: - xx_v - the Fortran90 pointer to the array
2002: Output Parameter:
2003: . ierr - error code
2005: Example of Usage:
2006: .vb
2007: PetscInt, pointer xx_v(:)
2008: ....
2009: call ISGetIndicesF90(x,xx_v,ierr)
2010: a = xx_v(3)
2011: call ISRestoreIndicesF90(x,xx_v,ierr)
2012: .ve
2014: Level: intermediate
2016: .seealso: `ISGetIndicesF90()`, `ISGetIndices()`, `ISRestoreIndices()`
2018: M*/
2020: /*MC
2021: ISBlockGetIndicesF90 - Accesses the elements of an index set from Fortran90.
2022: The users should call ISBlockRestoreIndicesF90() after having looked at the
2023: indices. The user should NOT change the indices.
2025: Synopsis:
2026: ISBlockGetIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
2028: Not collective
2030: Input Parameter:
2031: . x - index set
2033: Output Parameters:
2034: + xx_v - the Fortran90 pointer to the array
2035: - ierr - error code
2036: Example of Usage:
2037: .vb
2038: PetscInt, pointer xx_v(:)
2039: ....
2040: call ISBlockGetIndicesF90(x,xx_v,ierr)
2041: a = xx_v(3)
2042: call ISBlockRestoreIndicesF90(x,xx_v,ierr)
2043: .ve
2045: Level: intermediate
2047: .seealso: `ISBlockRestoreIndicesF90()`, `ISGetIndices()`, `ISRestoreIndices()`,
2048: `ISRestoreIndices()`
2050: M*/
2052: /*MC
2053: ISBlockRestoreIndicesF90 - Restores an index set to a usable state after
2054: a call to ISBlockGetIndicesF90().
2056: Synopsis:
2057: ISBlockRestoreIndicesF90(IS x,{integer, pointer :: xx_v(:)},integer ierr)
2059: Not Collective
2061: Input Parameters:
2062: + x - index set
2063: - xx_v - the Fortran90 pointer to the array
2065: Output Parameter:
2066: . ierr - error code
2068: Example of Usage:
2069: .vb
2070: PetscInt, pointer xx_v(:)
2071: ....
2072: call ISBlockGetIndicesF90(x,xx_v,ierr)
2073: a = xx_v(3)
2074: call ISBlockRestoreIndicesF90(x,xx_v,ierr)
2075: .ve
2077: Notes:
2078: Not yet supported for all F90 compilers
2080: Level: intermediate
2082: .seealso: `ISBlockGetIndicesF90()`, `ISGetIndices()`, `ISRestoreIndices()`, `ISRestoreIndicesF90()`
2084: M*/