Actual source code: block.c
2: /*
3: Provides the functions for index sets (IS) defined by a list of integers.
4: These are for blocks of data, each block is indicated with a single integer.
5: */
6: #include <petsc/private/isimpl.h>
7: #include <petscviewer.h>
9: typedef struct {
10: PetscBool sorted; /* are the blocks sorted? */
11: PetscBool allocated; /* did we allocate the index array ourselves? */
12: PetscInt *idx;
13: } IS_Block;
15: static PetscErrorCode ISDestroy_Block(IS is)
16: {
17: IS_Block *sub = (IS_Block *)is->data;
19: if (sub->allocated) PetscFree(sub->idx);
20: PetscObjectComposeFunction((PetscObject)is, "ISBlockSetIndices_C", NULL);
21: PetscObjectComposeFunction((PetscObject)is, "ISBlockGetIndices_C", NULL);
22: PetscObjectComposeFunction((PetscObject)is, "ISBlockRestoreIndices_C", NULL);
23: PetscObjectComposeFunction((PetscObject)is, "ISBlockGetSize_C", NULL);
24: PetscObjectComposeFunction((PetscObject)is, "ISBlockGetLocalSize_C", NULL);
25: PetscObjectComposeFunction((PetscObject)is, "ISShift_C", NULL);
26: PetscFree(is->data);
27: return 0;
28: }
30: static PetscErrorCode ISLocate_Block(IS is, PetscInt key, PetscInt *location)
31: {
32: IS_Block *sub = (IS_Block *)is->data;
33: PetscInt numIdx, i, bs, bkey, mkey;
34: PetscBool sorted;
36: PetscLayoutGetBlockSize(is->map, &bs);
37: PetscLayoutGetSize(is->map, &numIdx);
38: numIdx /= bs;
39: bkey = key / bs;
40: mkey = key % bs;
41: if (mkey < 0) {
42: bkey--;
43: mkey += bs;
44: }
45: ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, &sorted);
46: if (sorted) {
47: PetscFindInt(bkey, numIdx, sub->idx, location);
48: } else {
49: const PetscInt *idx = sub->idx;
51: *location = -1;
52: for (i = 0; i < numIdx; i++) {
53: if (idx[i] == bkey) {
54: *location = i;
55: break;
56: }
57: }
58: }
59: if (*location >= 0) *location = *location * bs + mkey;
60: return 0;
61: }
63: static PetscErrorCode ISGetIndices_Block(IS in, const PetscInt *idx[])
64: {
65: IS_Block *sub = (IS_Block *)in->data;
66: PetscInt i, j, k, bs, n, *ii, *jj;
68: PetscLayoutGetBlockSize(in->map, &bs);
69: PetscLayoutGetLocalSize(in->map, &n);
70: n /= bs;
71: if (bs == 1) *idx = sub->idx;
72: else {
73: if (n) {
74: PetscMalloc1(bs * n, &jj);
75: *idx = jj;
76: k = 0;
77: ii = sub->idx;
78: for (i = 0; i < n; i++)
79: for (j = 0; j < bs; j++) jj[k++] = bs * ii[i] + j;
80: } else {
81: /* do not malloc for zero size because F90Array1dCreate() inside ISRestoreArrayF90() does not keep array when zero length array */
82: *idx = NULL;
83: }
84: }
85: return 0;
86: }
88: static PetscErrorCode ISRestoreIndices_Block(IS is, const PetscInt *idx[])
89: {
90: IS_Block *sub = (IS_Block *)is->data;
91: PetscInt bs;
93: PetscLayoutGetBlockSize(is->map, &bs);
94: if (bs != 1) {
95: PetscFree(*(void **)idx);
96: } else {
97: /* F90Array1dCreate() inside ISRestoreArrayF90() does not keep array when zero length array */
99: }
100: return 0;
101: }
103: static PetscErrorCode ISInvertPermutation_Block(IS is, PetscInt nlocal, IS *isout)
104: {
105: IS_Block *sub = (IS_Block *)is->data;
106: PetscInt i, *ii, bs, n, *idx = sub->idx;
107: PetscMPIInt size;
109: MPI_Comm_size(PetscObjectComm((PetscObject)is), &size);
110: PetscLayoutGetBlockSize(is->map, &bs);
111: PetscLayoutGetLocalSize(is->map, &n);
112: n /= bs;
113: if (size == 1) {
114: PetscMalloc1(n, &ii);
115: for (i = 0; i < n; i++) ii[idx[i]] = i;
116: ISCreateBlock(PETSC_COMM_SELF, bs, n, ii, PETSC_OWN_POINTER, isout);
117: ISSetPermutation(*isout);
118: } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No inversion written yet for block IS");
119: return 0;
120: }
122: static PetscErrorCode ISView_Block(IS is, PetscViewer viewer)
123: {
124: IS_Block *sub = (IS_Block *)is->data;
125: PetscInt i, bs, n, *idx = sub->idx;
126: PetscBool iascii, ibinary;
128: PetscLayoutGetBlockSize(is->map, &bs);
129: PetscLayoutGetLocalSize(is->map, &n);
130: n /= bs;
131: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii);
132: PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERBINARY, &ibinary);
133: if (iascii) {
134: PetscViewerFormat fmt;
136: PetscViewerGetFormat(viewer, &fmt);
137: if (fmt == PETSC_VIEWER_ASCII_MATLAB) {
138: IS ist;
139: const char *name;
140: const PetscInt *idx;
141: PetscInt n;
143: PetscObjectGetName((PetscObject)is, &name);
144: ISGetLocalSize(is, &n);
145: ISGetIndices(is, &idx);
146: ISCreateGeneral(PetscObjectComm((PetscObject)is), n, idx, PETSC_USE_POINTER, &ist);
147: PetscObjectSetName((PetscObject)ist, name);
148: ISView(ist, viewer);
149: ISDestroy(&ist);
150: ISRestoreIndices(is, &idx);
151: } else {
152: PetscBool isperm;
154: ISGetInfo(is, IS_PERMUTATION, IS_GLOBAL, PETSC_FALSE, &isperm);
155: if (isperm) PetscViewerASCIIPrintf(viewer, "Block Index set is permutation\n");
156: PetscViewerASCIIPushSynchronized(viewer);
157: PetscViewerASCIISynchronizedPrintf(viewer, "Block size %" PetscInt_FMT "\n", bs);
158: PetscViewerASCIISynchronizedPrintf(viewer, "Number of block indices in set %" PetscInt_FMT "\n", n);
159: PetscViewerASCIISynchronizedPrintf(viewer, "The first indices of each block are\n");
160: for (i = 0; i < n; i++) PetscViewerASCIISynchronizedPrintf(viewer, "Block %" PetscInt_FMT " Index %" PetscInt_FMT "\n", i, idx[i]);
161: PetscViewerFlush(viewer);
162: PetscViewerASCIIPopSynchronized(viewer);
163: }
164: } else if (ibinary) ISView_Binary(is, viewer);
165: return 0;
166: }
168: static PetscErrorCode ISSort_Block(IS is)
169: {
170: IS_Block *sub = (IS_Block *)is->data;
171: PetscInt bs, n;
173: PetscLayoutGetBlockSize(is->map, &bs);
174: PetscLayoutGetLocalSize(is->map, &n);
175: PetscIntSortSemiOrdered(n / bs, sub->idx);
176: return 0;
177: }
179: static PetscErrorCode ISSortRemoveDups_Block(IS is)
180: {
181: IS_Block *sub = (IS_Block *)is->data;
182: PetscInt bs, n, nb;
183: PetscBool sorted;
185: PetscLayoutGetBlockSize(is->map, &bs);
186: PetscLayoutGetLocalSize(is->map, &n);
187: nb = n / bs;
188: ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, &sorted);
189: if (sorted) {
190: PetscSortedRemoveDupsInt(&nb, sub->idx);
191: } else {
192: PetscSortRemoveDupsInt(&nb, sub->idx);
193: }
194: PetscLayoutDestroy(&is->map);
195: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), nb * bs, PETSC_DECIDE, bs, &is->map);
196: return 0;
197: }
199: static PetscErrorCode ISSorted_Block(IS is, PetscBool *flg)
200: {
201: ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, flg);
202: return 0;
203: }
205: static PetscErrorCode ISSortedLocal_Block(IS is, PetscBool *flg)
206: {
207: IS_Block *sub = (IS_Block *)is->data;
208: PetscInt n, bs, i, *idx;
210: PetscLayoutGetLocalSize(is->map, &n);
211: PetscLayoutGetBlockSize(is->map, &bs);
212: n /= bs;
213: idx = sub->idx;
214: for (i = 1; i < n; i++)
215: if (idx[i] < idx[i - 1]) break;
216: if (i < n) *flg = PETSC_FALSE;
217: else *flg = PETSC_TRUE;
218: return 0;
219: }
221: static PetscErrorCode ISUniqueLocal_Block(IS is, PetscBool *flg)
222: {
223: IS_Block *sub = (IS_Block *)is->data;
224: PetscInt n, bs, i, *idx, *idxcopy = NULL;
225: PetscBool sortedLocal;
227: PetscLayoutGetLocalSize(is->map, &n);
228: PetscLayoutGetBlockSize(is->map, &bs);
229: n /= bs;
230: idx = sub->idx;
231: ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, &sortedLocal);
232: if (!sortedLocal) {
233: PetscMalloc1(n, &idxcopy);
234: PetscArraycpy(idxcopy, idx, n);
235: PetscIntSortSemiOrdered(n, idxcopy);
236: idx = idxcopy;
237: }
238: for (i = 1; i < n; i++)
239: if (idx[i] == idx[i - 1]) break;
240: if (i < n) *flg = PETSC_FALSE;
241: else *flg = PETSC_TRUE;
242: PetscFree(idxcopy);
243: return 0;
244: }
246: static PetscErrorCode ISPermutationLocal_Block(IS is, PetscBool *flg)
247: {
248: IS_Block *sub = (IS_Block *)is->data;
249: PetscInt n, bs, i, *idx, *idxcopy = NULL;
250: PetscBool sortedLocal;
252: PetscLayoutGetLocalSize(is->map, &n);
253: PetscLayoutGetBlockSize(is->map, &bs);
254: n /= bs;
255: idx = sub->idx;
256: ISGetInfo(is, IS_SORTED, IS_LOCAL, PETSC_TRUE, &sortedLocal);
257: if (!sortedLocal) {
258: PetscMalloc1(n, &idxcopy);
259: PetscArraycpy(idxcopy, idx, n);
260: PetscIntSortSemiOrdered(n, idxcopy);
261: idx = idxcopy;
262: }
263: for (i = 0; i < n; i++)
264: if (idx[i] != i) break;
265: if (i < n) *flg = PETSC_FALSE;
266: else *flg = PETSC_TRUE;
267: PetscFree(idxcopy);
268: return 0;
269: }
271: static PetscErrorCode ISIntervalLocal_Block(IS is, PetscBool *flg)
272: {
273: IS_Block *sub = (IS_Block *)is->data;
274: PetscInt n, bs, i, *idx;
276: PetscLayoutGetLocalSize(is->map, &n);
277: PetscLayoutGetBlockSize(is->map, &bs);
278: n /= bs;
279: idx = sub->idx;
280: for (i = 1; i < n; i++)
281: if (idx[i] != idx[i - 1] + 1) break;
282: if (i < n) *flg = PETSC_FALSE;
283: else *flg = PETSC_TRUE;
284: return 0;
285: }
287: static PetscErrorCode ISDuplicate_Block(IS is, IS *newIS)
288: {
289: IS_Block *sub = (IS_Block *)is->data;
290: PetscInt bs, n;
292: PetscLayoutGetBlockSize(is->map, &bs);
293: PetscLayoutGetLocalSize(is->map, &n);
294: n /= bs;
295: ISCreateBlock(PetscObjectComm((PetscObject)is), bs, n, sub->idx, PETSC_COPY_VALUES, newIS);
296: return 0;
297: }
299: static PetscErrorCode ISCopy_Block(IS is, IS isy)
300: {
301: IS_Block *is_block = (IS_Block *)is->data, *isy_block = (IS_Block *)isy->data;
302: PetscInt bs, n;
304: PetscLayoutGetBlockSize(is->map, &bs);
305: PetscLayoutGetLocalSize(is->map, &n);
306: PetscArraycpy(isy_block->idx, is_block->idx, n / bs);
307: return 0;
308: }
310: static PetscErrorCode ISOnComm_Block(IS is, MPI_Comm comm, PetscCopyMode mode, IS *newis)
311: {
312: IS_Block *sub = (IS_Block *)is->data;
313: PetscInt bs, n;
316: PetscLayoutGetBlockSize(is->map, &bs);
317: PetscLayoutGetLocalSize(is->map, &n);
318: ISCreateBlock(comm, bs, n / bs, sub->idx, mode, newis);
319: return 0;
320: }
322: static PetscErrorCode ISShift_Block(IS is, PetscInt shift, IS isy)
323: {
324: IS_Block *isb = (IS_Block *)is->data;
325: IS_Block *isby = (IS_Block *)isy->data;
326: PetscInt i, n, bs;
328: PetscLayoutGetLocalSize(is->map, &n);
329: PetscLayoutGetBlockSize(is->map, &bs);
330: shift /= bs;
331: for (i = 0; i < n / bs; i++) isby->idx[i] = isb->idx[i] + shift;
332: return 0;
333: }
335: static PetscErrorCode ISSetBlockSize_Block(IS is, PetscInt bs)
336: {
338: PetscLayoutSetBlockSize(is->map, bs);
339: return 0;
340: }
342: static PetscErrorCode ISToGeneral_Block(IS inis)
343: {
344: IS_Block *sub = (IS_Block *)inis->data;
345: PetscInt bs, n;
346: const PetscInt *idx;
348: ISGetBlockSize(inis, &bs);
349: ISGetLocalSize(inis, &n);
350: ISGetIndices(inis, &idx);
351: if (bs == 1) {
352: PetscCopyMode mode = sub->allocated ? PETSC_OWN_POINTER : PETSC_USE_POINTER;
353: sub->allocated = PETSC_FALSE; /* prevent deallocation when changing the subtype*/
354: ISSetType(inis, ISGENERAL);
355: ISGeneralSetIndices(inis, n, idx, mode);
356: } else {
357: ISSetType(inis, ISGENERAL);
358: ISGeneralSetIndices(inis, n, idx, PETSC_OWN_POINTER);
359: }
360: return 0;
361: }
363: static struct _ISOps myops = {ISGetIndices_Block, ISRestoreIndices_Block, ISInvertPermutation_Block, ISSort_Block, ISSortRemoveDups_Block, ISSorted_Block, ISDuplicate_Block, ISDestroy_Block, ISView_Block, ISLoad_Default, ISCopy_Block, ISToGeneral_Block, ISOnComm_Block, ISSetBlockSize_Block, NULL, ISLocate_Block,
364: /* we can have specialized local routines for determining properties,
365: * but unless the block size is the same on each process (which is not guaranteed at
366: * the moment), then trying to do something specialized for global properties is too
367: * complicated */
368: ISSortedLocal_Block, NULL, ISUniqueLocal_Block, NULL, ISPermutationLocal_Block, NULL, ISIntervalLocal_Block, NULL};
370: /*@
371: ISBlockSetIndices - Set integers representing blocks of indices in an index set.
373: Collective on IS
375: Input Parameters:
376: + is - the index set
377: . bs - number of elements in each block
378: . n - the length of the index set (the number of blocks)
379: . idx - the list of integers, one for each block, the integers contain the index of the first index of each block divided by the block size
380: - mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported
382: Notes:
383: When the communicator is not MPI_COMM_SELF, the operations on the
384: index sets, IS, are NOT conceptually the same as MPI_Group operations.
385: The index sets are then distributed sets of indices and thus certain operations
386: on them are collective.
388: Example:
389: If you wish to index the values {0,1,4,5}, then use
390: a block size of 2 and idx of {0,2}.
392: Level: beginner
394: .seealso: `ISCreateStride()`, `ISCreateGeneral()`, `ISAllGather()`, `ISCreateBlock()`, `ISBLOCK`, `ISGeneralSetIndices()`
395: @*/
396: PetscErrorCode ISBlockSetIndices(IS is, PetscInt bs, PetscInt n, const PetscInt idx[], PetscCopyMode mode)
397: {
398: ISClearInfoCache(is, PETSC_FALSE);
399: PetscUseMethod(is, "ISBlockSetIndices_C", (IS, PetscInt, PetscInt, const PetscInt[], PetscCopyMode), (is, bs, n, idx, mode));
400: return 0;
401: }
403: static PetscErrorCode ISBlockSetIndices_Block(IS is, PetscInt bs, PetscInt n, const PetscInt idx[], PetscCopyMode mode)
404: {
405: PetscInt i, min, max;
406: IS_Block *sub = (IS_Block *)is->data;
407: PetscLayout map;
413: PetscLayoutCreateFromSizes(PetscObjectComm((PetscObject)is), n * bs, is->map->N, bs, &map);
414: PetscLayoutDestroy(&is->map);
415: is->map = map;
417: if (sub->allocated) PetscFree(sub->idx);
418: if (mode == PETSC_COPY_VALUES) {
419: PetscMalloc1(n, &sub->idx);
420: PetscArraycpy(sub->idx, idx, n);
421: sub->allocated = PETSC_TRUE;
422: } else if (mode == PETSC_OWN_POINTER) {
423: sub->idx = (PetscInt *)idx;
424: sub->allocated = PETSC_TRUE;
425: } else if (mode == PETSC_USE_POINTER) {
426: sub->idx = (PetscInt *)idx;
427: sub->allocated = PETSC_FALSE;
428: }
430: if (n) {
431: min = max = idx[0];
432: for (i = 1; i < n; i++) {
433: if (idx[i] < min) min = idx[i];
434: if (idx[i] > max) max = idx[i];
435: }
436: is->min = bs * min;
437: is->max = bs * max + bs - 1;
438: } else {
439: is->min = PETSC_MAX_INT;
440: is->max = PETSC_MIN_INT;
441: }
442: return 0;
443: }
445: /*@
446: ISCreateBlock - Creates a data structure for an index set containing
447: a list of integers. Each integer represents a fixed block size set of indices.
449: Collective
451: Input Parameters:
452: + comm - the MPI communicator
453: . bs - number of elements in each block
454: . n - the length of the index set (the number of blocks)
455: . idx - the list of integers, one for each block, the integers contain the index of the first entry of each block divided by the block size
456: - mode - see PetscCopyMode, only PETSC_COPY_VALUES and PETSC_OWN_POINTER are supported in this routine
458: Output Parameter:
459: . is - the new index set
461: Notes:
462: When the communicator is not MPI_COMM_SELF, the operations on the
463: index sets, IS, are NOT conceptually the same as MPI_Group operations.
464: The index sets are then distributed sets of indices and thus certain operations
465: on them are collective.
467: Example:
468: If you wish to index the values {0,1,6,7}, then use
469: a block size of 2 and idx of {0,3}.
471: Level: beginner
473: .seealso: `ISCreateStride()`, `ISCreateGeneral()`, `ISAllGather()`, `ISBlockSetIndices()`, `ISBLOCK`, `ISGENERAL`
474: @*/
475: PetscErrorCode ISCreateBlock(MPI_Comm comm, PetscInt bs, PetscInt n, const PetscInt idx[], PetscCopyMode mode, IS *is)
476: {
482: ISCreate(comm, is);
483: ISSetType(*is, ISBLOCK);
484: ISBlockSetIndices(*is, bs, n, idx, mode);
485: return 0;
486: }
488: static PetscErrorCode ISBlockGetIndices_Block(IS is, const PetscInt *idx[])
489: {
490: IS_Block *sub = (IS_Block *)is->data;
492: *idx = sub->idx;
493: return 0;
494: }
496: static PetscErrorCode ISBlockRestoreIndices_Block(IS is, const PetscInt *idx[])
497: {
498: return 0;
499: }
501: /*@C
502: ISBlockGetIndices - Gets the indices associated with each block.
504: Not Collective
506: Input Parameter:
507: . is - the index set
509: Output Parameter:
510: . idx - the integer indices, one for each block and count of block not indices
512: Level: intermediate
514: .seealso: `ISGetIndices()`, `ISBlockRestoreIndices()`, `ISBLOCK`, `ISBlockSetIndices()`, `ISCreateBlock()`
515: @*/
516: PetscErrorCode ISBlockGetIndices(IS is, const PetscInt *idx[])
517: {
518: PetscUseMethod(is, "ISBlockGetIndices_C", (IS, const PetscInt *[]), (is, idx));
519: return 0;
520: }
522: /*@C
523: ISBlockRestoreIndices - Restores the indices associated with each block.
525: Not Collective
527: Input Parameter:
528: . is - the index set
530: Output Parameter:
531: . idx - the integer indices
533: Level: intermediate
535: .seealso: `ISRestoreIndices()`, `ISBlockGetIndices()`
536: @*/
537: PetscErrorCode ISBlockRestoreIndices(IS is, const PetscInt *idx[])
538: {
539: PetscUseMethod(is, "ISBlockRestoreIndices_C", (IS, const PetscInt *[]), (is, idx));
540: return 0;
541: }
543: /*@
544: ISBlockGetLocalSize - Returns the local number of blocks in the index set.
546: Not Collective
548: Input Parameter:
549: . is - the index set
551: Output Parameter:
552: . size - the local number of blocks
554: Level: intermediate
556: .seealso: `ISGetBlockSize()`, `ISBlockGetSize()`, `ISGetSize()`, `ISCreateBlock()`, `ISBLOCK`
557: @*/
558: PetscErrorCode ISBlockGetLocalSize(IS is, PetscInt *size)
559: {
560: PetscUseMethod(is, "ISBlockGetLocalSize_C", (IS, PetscInt *), (is, size));
561: return 0;
562: }
564: static PetscErrorCode ISBlockGetLocalSize_Block(IS is, PetscInt *size)
565: {
566: PetscInt bs, n;
568: PetscLayoutGetBlockSize(is->map, &bs);
569: PetscLayoutGetLocalSize(is->map, &n);
570: *size = n / bs;
571: return 0;
572: }
574: /*@
575: ISBlockGetSize - Returns the global number of blocks in the index set.
577: Not Collective
579: Input Parameter:
580: . is - the index set
582: Output Parameter:
583: . size - the global number of blocks
585: Level: intermediate
587: .seealso: `ISGetBlockSize()`, `ISBlockGetLocalSize()`, `ISGetSize()`, `ISCreateBlock()`, `ISBLOCK`
588: @*/
589: PetscErrorCode ISBlockGetSize(IS is, PetscInt *size)
590: {
591: PetscUseMethod(is, "ISBlockGetSize_C", (IS, PetscInt *), (is, size));
592: return 0;
593: }
595: static PetscErrorCode ISBlockGetSize_Block(IS is, PetscInt *size)
596: {
597: PetscInt bs, N;
599: PetscLayoutGetBlockSize(is->map, &bs);
600: PetscLayoutGetSize(is->map, &N);
601: *size = N / bs;
602: return 0;
603: }
605: PETSC_EXTERN PetscErrorCode ISCreate_Block(IS is)
606: {
607: IS_Block *sub;
609: PetscNew(&sub);
610: is->data = (void *)sub;
611: PetscMemcpy(is->ops, &myops, sizeof(myops));
612: PetscObjectComposeFunction((PetscObject)is, "ISBlockSetIndices_C", ISBlockSetIndices_Block);
613: PetscObjectComposeFunction((PetscObject)is, "ISBlockGetIndices_C", ISBlockGetIndices_Block);
614: PetscObjectComposeFunction((PetscObject)is, "ISBlockRestoreIndices_C", ISBlockRestoreIndices_Block);
615: PetscObjectComposeFunction((PetscObject)is, "ISBlockGetSize_C", ISBlockGetSize_Block);
616: PetscObjectComposeFunction((PetscObject)is, "ISBlockGetLocalSize_C", ISBlockGetLocalSize_Block);
617: PetscObjectComposeFunction((PetscObject)is, "ISShift_C", ISShift_Block);
618: return 0;
619: }