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